multiple_regression.c 11 KB

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