multiple_regression.c 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2009, 2010, 2011, 2015-2016 Université de Bordeaux
  4. * Copyright (C) 2010, 2011 CNRS
  5. * Copyright (C) 2016 Inria
  6. *
  7. * StarPU is free software; you can redistribute it and/or modify
  8. * it under the terms of the GNU Lesser General Public License as published by
  9. * the Free Software Foundation; either version 2.1 of the License, or (at
  10. * your option) any later version.
  11. *
  12. * StarPU is distributed in the hope that it will be useful, but
  13. * WITHOUT ANY WARRANTY; without even the implied warranty of
  14. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
  15. *
  16. * See the GNU Lesser General Public License in COPYING.LGPL for more details.
  17. */
  18. /* Code for computing multiple linear regression */
  19. #include <core/perfmodel/multiple_regression.h>
  20. typedef long int integer;
  21. typedef double doublereal;
  22. #ifdef STARPU_MLR_MODEL
  23. int dgels_(char *trans, integer *m, integer *n, integer *nrhs, doublereal *a, integer *lda, doublereal *b, integer *ldb, doublereal *work, integer *lwork, integer *info);
  24. #endif //STARPU_MLR_MODEL
  25. static long count_file_lines(FILE *f)
  26. {
  27. int ch, lines=0;
  28. while(!feof(f))
  29. {
  30. ch = fgetc(f);
  31. if(ch == '\n')
  32. {
  33. lines++;
  34. }
  35. }
  36. rewind(f);
  37. return lines;
  38. }
  39. static void dump_multiple_regression_list(double *mpar, double *my, int start, unsigned nparameters, struct starpu_perfmodel_history_list *list_history)
  40. {
  41. struct starpu_perfmodel_history_list *ptr = list_history;
  42. int i = start;
  43. unsigned j;
  44. while (ptr)
  45. {
  46. my[i] = ptr->entry->duration;
  47. for(j=0; j<nparameters; j++)
  48. mpar[i*nparameters+j] = ptr->entry->parameters[j];
  49. ptr = ptr->next;
  50. i++;
  51. }
  52. }
  53. static void load_old_calibration(double *mx, double *my, unsigned nparameters, char *filepath)
  54. {
  55. char buffer[1024];
  56. char *record,*line;
  57. int i=0,j=0;
  58. FILE *f=NULL;
  59. f = fopen(filepath, "a+");
  60. STARPU_ASSERT_MSG(f, "Could not save performance model into the file %s\n", filepath);
  61. line=fgets(buffer,sizeof(buffer),f);//skipping first line
  62. while((line=fgets(buffer,sizeof(buffer),f))!=NULL)
  63. {
  64. record = strtok(line,",");
  65. my[i] = atof(record);
  66. record = strtok(NULL,",");
  67. j=0;
  68. while(record != NULL)
  69. {
  70. mx[i*nparameters+j] = atof(record) ;
  71. ++j;
  72. record = strtok(NULL,",");
  73. }
  74. ++i ;
  75. }
  76. fclose(f);
  77. }
  78. static long find_long_list_size(struct starpu_perfmodel_history_list *list_history)
  79. {
  80. long cnt = 0;
  81. struct starpu_perfmodel_history_list *ptr = list_history;
  82. while (ptr)
  83. {
  84. cnt++;
  85. ptr = ptr->next;
  86. }
  87. return cnt;
  88. }
  89. #ifdef STARPU_MLR_MODEL
  90. int dgels_multiple_reg_coeff(double *mpar, double *my, long nn, unsigned ncoeff, unsigned nparameters, double *coeff, unsigned **combinations)
  91. {
  92. /* Arguments */
  93. /* ========= */
  94. /* TRANS (input) CHARACTER*1 */
  95. /* = 'N': the linear system involves A; */
  96. /* = 'T': the linear system involves A**T. */
  97. /* M (input) INTEGER */
  98. /* The number of rows of the matrix A. M >= 0. */
  99. /* N (input) INTEGER */
  100. /* The number of columns of the matrix A. N >= 0. */
  101. /* NRHS (input) INTEGER */
  102. /* The number of right hand sides, i.e., the number of */
  103. /* columns of the matrices B and X. NRHS >=0. */
  104. /* A (input/output) DOUBLE PRECISION array, dimension (LDA,N) */
  105. /* On entry, the M-by-N matrix A. */
  106. /* On exit, */
  107. /* if M >= N, A is overwritten by details of its QR */
  108. /* factorization as returned by DGEQRF; */
  109. /* if M < N, A is overwritten by details of its LQ */
  110. /* factorization as returned by DGELQF. */
  111. /* LDA (input) INTEGER */
  112. /* The leading dimension of the array A. LDA >= max(1,M). */
  113. /* B (input/output) DOUBLE PRECISION array, dimension (LDB,NRHS) */
  114. /* On entry, the matrix B of right hand side vectors, stored */
  115. /* columnwise; B is M-by-NRHS if TRANS = 'N', or N-by-NRHS */
  116. /* if TRANS = 'T'. */
  117. /* On exit, if INFO = 0, B is overwritten by the solution */
  118. /* vectors, stored columnwise: */
  119. /* if TRANS = 'N' and m >= n, rows 1 to n of B contain the least */
  120. /* squares solution vectors; the residual sum of squares for the */
  121. /* solution in each column is given by the sum of squares of */
  122. /* elements N+1 to M in that column; */
  123. /* if TRANS = 'N' and m < n, rows 1 to N of B contain the */
  124. /* minimum norm solution vectors; */
  125. /* if TRANS = 'T' and m >= n, rows 1 to M of B contain the */
  126. /* minimum norm solution vectors; */
  127. /* if TRANS = 'T' and m < n, rows 1 to M of B contain the */
  128. /* least squares solution vectors; the residual sum of squares */
  129. /* for the solution in each column is given by the sum of */
  130. /* squares of elements M+1 to N in that column. */
  131. /* LDB (input) INTEGER */
  132. /* The leading dimension of the array B. LDB >= MAX(1,M,N). */
  133. /* WORK (workspace/output) DOUBLE PRECISION array, dimension (MAX(1,LWORK)) */
  134. /* On exit, if INFO = 0, WORK(1) returns the optimal LWORK. */
  135. /* LWORK (input) INTEGER */
  136. /* The dimension of the array WORK. */
  137. /* LWORK >= max( 1, MN + max( MN, NRHS ) ). */
  138. /* For optimal performance, */
  139. /* LWORK >= max( 1, MN + max( MN, NRHS )*NB ). */
  140. /* where MN = min(M,N) and NB is the optimum block size. */
  141. /* If LWORK = -1, then a workspace query is assumed; the routine */
  142. /* only calculates the optimal size of the WORK array, returns */
  143. /* this value as the first entry of the WORK array, and no error */
  144. /* message related to LWORK is issued by XERBLA. */
  145. /* INFO (output) INTEGER */
  146. /* = 0: successful exit */
  147. /* < 0: if INFO = -i, the i-th argument had an illegal value */
  148. /* > 0: if INFO = i, the i-th diagonal element of the */
  149. /* triangular factor of A is zero, so that A does not have */
  150. /* full rank; the least squares solution could not be */
  151. /* computed. */
  152. /* ===================================================================== */
  153. if(nn <= ncoeff)
  154. {
  155. _STARPU_DISP("Warning: This function is not intended for the use when number of parameters is larger than the number of observations. Check how your matrices A and B were allocated or simply add more benchmarks.\n Multiple linear regression model will not be written into perfmodel file.\n");
  156. return 1;
  157. }
  158. char trans = 'N';
  159. integer m = nn;
  160. integer n = ncoeff;
  161. integer nrhs = 1; // number of columns of B and X (wich are vectors therefore nrhs=1)
  162. doublereal *X = malloc(sizeof(double)*n*m); // (/!\ modified at the output) contain the model and the different values of pararmters
  163. doublereal *Y = malloc(sizeof(double)*m);
  164. double coefficient;
  165. int i;
  166. unsigned j, k;
  167. for (i=0; i < m; i++)
  168. {
  169. Y[i] = my[i];
  170. X[i*n] = 1.;
  171. for (j=1; j < n; j++)
  172. {
  173. coefficient = 1.;
  174. for(k=0; k < nparameters; k++)
  175. {
  176. coefficient *= pow(mpar[i*nparameters+k],combinations[j-1][k]);
  177. }
  178. X[i*n+j] = coefficient;
  179. }
  180. }
  181. integer lda = m;
  182. integer ldb = m; //
  183. integer info = 0;
  184. integer lwork = n*2;
  185. doublereal *work = malloc(sizeof(double)*lwork); // (output)
  186. /* // Running LAPACK dgels_ */
  187. dgels_(&trans, &m, &n, &nrhs, X, &lda, Y, &ldb, work, &lwork, &info);
  188. /* Check for the full rank */
  189. if( info != 0 )
  190. {
  191. _STARPU_DISP("Warning: Problems when executing dgels_ function. It seems like the diagonal element %ld is zero.\n Multiple linear regression model will not be written into perfmodel file.\n", info);
  192. return 1;
  193. }
  194. /* Copy computed coefficients */
  195. for(i=0; i<ncoeff; i++)
  196. coeff[i] = Y[i];
  197. free(X);
  198. free(Y);
  199. free(work);
  200. return 0;
  201. }
  202. #endif //STARPU_MLR_MODEL
  203. /*
  204. Validating the accuracy of the coefficients.
  205. For the the validation is extremely basic, but it should be improved.
  206. */
  207. void validate(double *coeff, unsigned ncoeff, const char *codelet_name)
  208. {
  209. unsigned i;
  210. if (coeff[0] < 0)
  211. _STARPU_DISP("Warning: Constant computed by least square method is negative (%f). The model %s is likely to be inaccurate.\n", coeff[0], codelet_name);
  212. for(i=1; i<ncoeff; i++)
  213. if(fabs(coeff[i]) < 1E-10)
  214. _STARPU_DISP("Warning: Coefficient computed by least square method is extremelly small (%f). The model %s is likely to be inaccurate.\n", coeff[i], codelet_name);
  215. }
  216. int _starpu_multiple_regression(struct starpu_perfmodel_history_list *ptr, double *coeff, unsigned ncoeff, unsigned nparameters, const char **parameters_names, unsigned **combinations, const char *codelet_name)
  217. {
  218. long i;
  219. unsigned j;
  220. /* Computing number of rows */
  221. long n=find_long_list_size(ptr);
  222. STARPU_ASSERT(n);
  223. /* Reading old calibrations if necessary */
  224. FILE *f=NULL;
  225. char directory[300];
  226. snprintf(directory, 300, "%s/.starpu/sampling/codelets/tmp", _starpu_get_home_path());
  227. _starpu_mkpath_and_check(directory, S_IRWXU);
  228. char filepath[300];
  229. snprintf(filepath, 300, "%s/%s.out", directory,codelet_name);
  230. long old_lines=0;
  231. int calibrate = _starpu_get_calibrate_flag();
  232. if (calibrate==1)
  233. {
  234. f = fopen(filepath, "a+");
  235. STARPU_ASSERT_MSG(f, "Could not save performance model into the file %s\n", filepath);
  236. old_lines=count_file_lines(f);
  237. /* If the program is run for the first time the old_lines will be 0 */
  238. //STARPU_ASSERT(old_lines);
  239. n+=old_lines;
  240. fclose(f);
  241. }
  242. /* Allocating X and Y matrices */
  243. double *mpar = (double *) malloc(nparameters*n*sizeof(double));
  244. STARPU_ASSERT(mpar);
  245. double *my = (double *) malloc(n*sizeof(double));
  246. STARPU_ASSERT(my);
  247. /* Loading old calibration */
  248. if (calibrate==1)
  249. load_old_calibration(mpar, my, nparameters, filepath);
  250. /* Filling X and Y matrices with measured values */
  251. dump_multiple_regression_list(mpar, my, old_lines, nparameters, ptr);
  252. if (ncoeff!=0 && combinations!=NULL)
  253. {
  254. #ifdef STARPU_MLR_MODEL
  255. /* Computing coefficients using multiple linear regression */
  256. if(dgels_multiple_reg_coeff(mpar, my, n, ncoeff, nparameters, coeff, combinations))
  257. return 1;
  258. /* Basic validation of the model accuracy */
  259. validate(coeff, ncoeff, codelet_name);
  260. #else
  261. _STARPU_DISP("Warning: StarPU was compiled with '--disable-mlr' option or on Windows machine, thus multiple linear regression model will not be computed.\n");
  262. for(i=0; i<ncoeff; i++)
  263. coeff[i] = 0.;
  264. #endif //STARPU_MLR_MODEL
  265. }
  266. /* Preparing new output calibration file */
  267. if (calibrate==1 || calibrate==2)
  268. {
  269. if (old_lines > 0)
  270. {
  271. f = fopen(filepath, "a+");
  272. STARPU_ASSERT_MSG(f, "Could not save performance model into the file %s\n", filepath);
  273. }
  274. else
  275. {
  276. f = fopen(filepath, "w+");
  277. STARPU_ASSERT_MSG(f, "Could not save performance model into the file %s\n", filepath);
  278. fprintf(f, "Duration");
  279. for(j=0; j<nparameters; j++)
  280. {
  281. if(parameters_names != NULL && parameters_names[j]!= NULL)
  282. fprintf(f, ", %s", parameters_names[j]);
  283. else
  284. fprintf(f, ", P%d", j);
  285. }
  286. }
  287. }
  288. /* Writing parameters to calibration file */
  289. if (calibrate==1 || calibrate==2)
  290. {
  291. for(i=old_lines; i<n; i++)
  292. {
  293. fprintf(f, "\n%f", my[i]);
  294. for(j=0; j<nparameters; j++)
  295. fprintf(f, ", %f", mpar[i*nparameters+j]);
  296. }
  297. fclose(f);
  298. }
  299. /* Cleanup */
  300. free(mpar);
  301. free(my);
  302. return 0;
  303. }