multiple_regression.c 10.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323
  1. /* StarPU --- Runtime system for heterogeneous multicore architectures.
  2. *
  3. * Copyright (C) 2009, 2010, 2011, 2015 Université de Bordeaux
  4. * Copyright (C) 2010, 2011 CNRS
  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. int dgels_(char *trans, integer *m, integer *n, integer *nrhs, doublereal *a, integer *lda, doublereal *b, integer *ldb, doublereal *work, integer *lwork, integer *info);
  22. static long count_file_lines(FILE *f)
  23. {
  24. int ch, lines=0;
  25. while(!feof(f))
  26. {
  27. ch = fgetc(f);
  28. if(ch == '\n')
  29. {
  30. lines++;
  31. }
  32. }
  33. rewind(f);
  34. return lines;
  35. }
  36. static void dump_multiple_regression_list(double *mpar, double *my, int start, unsigned nparameters, struct starpu_perfmodel_history_list *list_history)
  37. {
  38. struct starpu_perfmodel_history_list *ptr = list_history;
  39. int i = start;
  40. while (ptr)
  41. {
  42. my[i] = ptr->entry->duration;
  43. for(int j=0; j<nparameters; j++)
  44. mpar[i*nparameters+j] = ptr->entry->parameters[j];
  45. ptr = ptr->next;
  46. i++;
  47. }
  48. }
  49. static void load_old_calibration(double *mx, double *my, unsigned nparameters, FILE *f)
  50. {
  51. char buffer[1024];
  52. char *record,*line;
  53. int i=0,j=0;
  54. line=fgets(buffer,sizeof(buffer),f);//skipping first line
  55. while((line=fgets(buffer,sizeof(buffer),f))!=NULL)
  56. {
  57. record = strtok(line,",");
  58. my[i] = atof(record);
  59. record = strtok(NULL,",");
  60. j=0;
  61. while(record != NULL)
  62. {
  63. mx[i*nparameters+j] = atof(record) ;
  64. ++j;
  65. record = strtok(NULL,",");
  66. }
  67. ++i ;
  68. }
  69. }
  70. static long find_long_list_size(struct starpu_perfmodel_history_list *list_history)
  71. {
  72. long cnt = 0;
  73. struct starpu_perfmodel_history_list *ptr = list_history;
  74. while (ptr)
  75. {
  76. cnt++;
  77. ptr = ptr->next;
  78. }
  79. return cnt;
  80. }
  81. int dgels_multiple_reg_coeff(double *mpar, double *my, long nn, unsigned ncoeff, unsigned nparameters, double *coeff, unsigned **combinations)
  82. {
  83. /* Arguments */
  84. /* ========= */
  85. /* TRANS (input) CHARACTER*1 */
  86. /* = 'N': the linear system involves A; */
  87. /* = 'T': the linear system involves A**T. */
  88. /* M (input) INTEGER */
  89. /* The number of rows of the matrix A. M >= 0. */
  90. /* N (input) INTEGER */
  91. /* The number of columns of the matrix A. N >= 0. */
  92. /* NRHS (input) INTEGER */
  93. /* The number of right hand sides, i.e., the number of */
  94. /* columns of the matrices B and X. NRHS >=0. */
  95. /* A (input/output) DOUBLE PRECISION array, dimension (LDA,N) */
  96. /* On entry, the M-by-N matrix A. */
  97. /* On exit, */
  98. /* if M >= N, A is overwritten by details of its QR */
  99. /* factorization as returned by DGEQRF; */
  100. /* if M < N, A is overwritten by details of its LQ */
  101. /* factorization as returned by DGELQF. */
  102. /* LDA (input) INTEGER */
  103. /* The leading dimension of the array A. LDA >= max(1,M). */
  104. /* B (input/output) DOUBLE PRECISION array, dimension (LDB,NRHS) */
  105. /* On entry, the matrix B of right hand side vectors, stored */
  106. /* columnwise; B is M-by-NRHS if TRANS = 'N', or N-by-NRHS */
  107. /* if TRANS = 'T'. */
  108. /* On exit, if INFO = 0, B is overwritten by the solution */
  109. /* vectors, stored columnwise: */
  110. /* if TRANS = 'N' and m >= n, rows 1 to n of B contain the least */
  111. /* squares solution vectors; the residual sum of squares for the */
  112. /* solution in each column is given by the sum of squares of */
  113. /* elements N+1 to M in that column; */
  114. /* if TRANS = 'N' and m < n, rows 1 to N of B contain the */
  115. /* minimum norm solution vectors; */
  116. /* if TRANS = 'T' and m >= n, rows 1 to M of B contain the */
  117. /* minimum norm solution vectors; */
  118. /* if TRANS = 'T' and m < n, rows 1 to M of B contain the */
  119. /* least squares solution vectors; the residual sum of squares */
  120. /* for the solution in each column is given by the sum of */
  121. /* squares of elements M+1 to N in that column. */
  122. /* LDB (input) INTEGER */
  123. /* The leading dimension of the array B. LDB >= MAX(1,M,N). */
  124. /* WORK (workspace/output) DOUBLE PRECISION array, dimension (MAX(1,LWORK)) */
  125. /* On exit, if INFO = 0, WORK(1) returns the optimal LWORK. */
  126. /* LWORK (input) INTEGER */
  127. /* The dimension of the array WORK. */
  128. /* LWORK >= max( 1, MN + max( MN, NRHS ) ). */
  129. /* For optimal performance, */
  130. /* LWORK >= max( 1, MN + max( MN, NRHS )*NB ). */
  131. /* where MN = min(M,N) and NB is the optimum block size. */
  132. /* If LWORK = -1, then a workspace query is assumed; the routine */
  133. /* only calculates the optimal size of the WORK array, returns */
  134. /* this value as the first entry of the WORK array, and no error */
  135. /* message related to LWORK is issued by XERBLA. */
  136. /* INFO (output) INTEGER */
  137. /* = 0: successful exit */
  138. /* < 0: if INFO = -i, the i-th argument had an illegal value */
  139. /* > 0: if INFO = i, the i-th diagonal element of the */
  140. /* triangular factor of A is zero, so that A does not have */
  141. /* full rank; the least squares solution could not be */
  142. /* computed. */
  143. /* ===================================================================== */
  144. if(nn <= ncoeff)
  145. {
  146. _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");
  147. return 1;
  148. }
  149. char trans = 'N';
  150. integer m = nn;
  151. integer n = ncoeff;
  152. integer nrhs = 1; // number of columns of B and X (wich are vectors therefore nrhs=1)
  153. doublereal *X = malloc(sizeof(double)*n*m); // (/!\ modified at the output) contain the model and the different values of pararmters
  154. doublereal *Y = malloc(sizeof(double)*m);
  155. double coefficient;
  156. for (int i=0; i < m; i++)
  157. {
  158. Y[i] = my[i];
  159. X[i*n] = 1.;
  160. for (int j=1; j < n; j++)
  161. {
  162. coefficient = 1.;
  163. for(int k=0; k < nparameters; k++)
  164. {
  165. coefficient *= pow(mpar[i*nparameters+k],combinations[j-1][k]);
  166. }
  167. X[i*n+j] = coefficient;
  168. }
  169. }
  170. integer lda = m;
  171. integer ldb = m; //
  172. integer info = 0;
  173. integer lwork = n*2;
  174. doublereal *work = malloc(sizeof(double)*lwork); // (output)
  175. /* // Running LAPACK dgels_ */
  176. dgels_(&trans, &m, &n, &nrhs, X, &lda, Y, &ldb, work, &lwork, &info);
  177. /* Check for the full rank */
  178. if( info != 0 )
  179. {
  180. _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);
  181. return 1;
  182. }
  183. /* Copy computed coefficients */
  184. for(int i=0; i<ncoeff; i++)
  185. coeff[i] = Y[i];
  186. free(X);
  187. free(Y);
  188. free(work);
  189. return 0;
  190. }
  191. /*
  192. Validating the accuracy of the coefficients.
  193. For the the validation is extremely basic, but it should be improved.
  194. */
  195. void validate(double *coeff, unsigned ncoeff)
  196. {
  197. if (coeff[0] < 0)
  198. _STARPU_DISP("Warning: Constant computed by least square method is negative (%f). The model is likely to be inaccurate.\n", coeff[0]);
  199. for(int i=1; i<ncoeff; i++)
  200. if(coeff[i] < 1E-10)
  201. _STARPU_DISP("Warning: Coefficient computed by least square method is extremelly small (%f). The model is likely to be inaccurate.\n", coeff[i]);
  202. }
  203. int _starpu_multiple_regression(struct starpu_perfmodel_history_list *ptr, double *coeff, unsigned ncoeff, unsigned nparameters, unsigned **combinations, const char *codelet_name)
  204. {
  205. /* Computing number of rows */
  206. long n=find_long_list_size(ptr);
  207. STARPU_ASSERT(n);
  208. /* Reading old calibrations if necessary */
  209. FILE *f=NULL;
  210. char directory[100];
  211. snprintf(directory, 100, "%s/.starpu/sampling/codelets/tmp", _starpu_get_home_path());
  212. _starpu_mkpath_and_check(directory, S_IRWXU);
  213. char filepath[100];
  214. snprintf(filepath, 100, "%s/%s.out", directory,codelet_name);
  215. long old_lines=0;
  216. int calibrate = starpu_get_env_number("STARPU_CALIBRATE");
  217. if (calibrate==1)
  218. {
  219. f = fopen(filepath, "a+");
  220. STARPU_ASSERT_MSG(f, "Could not save performance model into the file %s\n", filepath);
  221. old_lines=count_file_lines(f);
  222. STARPU_ASSERT(old_lines);
  223. n+=old_lines;
  224. }
  225. /* Allocating X and Y matrices */
  226. double *mpar = (double *) malloc(nparameters*n*sizeof(double));
  227. STARPU_ASSERT(mpar);
  228. double *my = (double *) malloc(n*sizeof(double));
  229. STARPU_ASSERT(my);
  230. /* Loading old calibration */
  231. if (calibrate==1)
  232. load_old_calibration(mpar, my, nparameters, f);
  233. /* Filling X and Y matrices with measured values */
  234. dump_multiple_regression_list(mpar, my, old_lines, nparameters, ptr);
  235. /* Computing coefficients using multiple linear regression */
  236. if(dgels_multiple_reg_coeff(mpar, my, n, ncoeff, nparameters, coeff, combinations))
  237. return 1;
  238. /* Basic validation of the model accuracy */
  239. validate(coeff, ncoeff);
  240. /* Preparing new output calibration file */
  241. if (calibrate==2)
  242. {
  243. f = fopen(filepath, "w+");
  244. STARPU_ASSERT_MSG(f, "Could not save performance model into the file %s\n", filepath);
  245. fprintf(f, "Duration");
  246. for(int k=0; k < nparameters; k++)
  247. {
  248. fprintf(f, ", P%d", k);
  249. }
  250. }
  251. /* Writing parameters to calibration file */
  252. if (calibrate==1 || calibrate==2)
  253. {
  254. for(int i=old_lines; i<n; i++)
  255. {
  256. fprintf(f, "\n%f", my[i]);
  257. for(int j=0; j<nparameters;j++)
  258. fprintf(f, ", %f", mpar[i*nparameters+j]);
  259. }
  260. fclose(f);
  261. }
  262. /* Cleanup */
  263. free(mpar);
  264. free(my);
  265. return 0;
  266. }