multiple_regression.c 12 KB

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