multiple_regression.c 10 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348
  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. #ifdef TESTGSL
  20. #include <gsl/gsl_matrix.h>
  21. #include <gsl/gsl_math.h>
  22. #include <gsl/gsl_multifit.h>
  23. #endif //TESTGSL
  24. // From additional-lapack-func
  25. #ifdef MLR_MODEL
  26. #include "mindgels.h"
  27. #endif //MLR_MODEL
  28. typedef struct { int h, w; double *x;} matrix_t, *matrix;
  29. static long count_file_lines(FILE *f)
  30. {
  31. int ch, lines=0;
  32. while(!feof(f))
  33. {
  34. ch = fgetc(f);
  35. if(ch == '\n')
  36. {
  37. lines++;
  38. }
  39. }
  40. rewind(f);
  41. return lines;
  42. }
  43. static void dump_multiple_regression_list(double *mpar, double *my, int start, unsigned nparameters, struct starpu_perfmodel_history_list *list_history)
  44. {
  45. struct starpu_perfmodel_history_list *ptr = list_history;
  46. int i = start;
  47. while (ptr)
  48. {
  49. my[i] = ptr->entry->duration;
  50. for(int 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, FILE *f)
  57. {
  58. char buffer[1024];
  59. char *record,*line;
  60. int i=0,j=0;
  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. }
  77. static long find_long_list_size(struct starpu_perfmodel_history_list *list_history)
  78. {
  79. long cnt = 0;
  80. struct starpu_perfmodel_history_list *ptr = list_history;
  81. while (ptr)
  82. {
  83. cnt++;
  84. ptr = ptr->next;
  85. }
  86. return cnt;
  87. }
  88. #ifdef TESTGSL
  89. void gsl_multiple_reg_coeff(double *mpar, double *my, long n, unsigned ncoeff, unsigned nparameters, double *coeff, unsigned **combinations)
  90. {
  91. double coefficient;
  92. gsl_matrix *X = gsl_matrix_calloc(n, ncoeff);
  93. gsl_vector *Y = gsl_vector_alloc(n);
  94. gsl_vector *beta = gsl_vector_alloc(ncoeff);
  95. for (int i = 0; i < n; i++) {
  96. gsl_vector_set(Y, i, my[i]);
  97. gsl_matrix_set(X, i, 0, 1.);
  98. for (int j = 1; j < ncoeff; j++)
  99. {
  100. coefficient = 1.;
  101. for(int k=0; k < nparameters; k++)
  102. {
  103. coefficient *= pow(mpar[i*nparameters+k],combinations[j-1][k]);
  104. }
  105. gsl_matrix_set(X, i, j, coefficient);
  106. }
  107. }
  108. double chisq;
  109. gsl_matrix *cov = gsl_matrix_alloc(ncoeff, ncoeff);
  110. gsl_multifit_linear_workspace * wspc = gsl_multifit_linear_alloc(n, ncoeff);
  111. gsl_multifit_linear(X, Y, beta, cov, &chisq, wspc);
  112. for(int i=0; i<ncoeff; i++)
  113. coeff[i] = gsl_vector_get(beta, i);
  114. gsl_matrix_free(X);
  115. gsl_matrix_free(cov);
  116. gsl_vector_free(Y);
  117. gsl_vector_free(beta);
  118. gsl_multifit_linear_free(wspc);
  119. }
  120. #endif //TESTGSL
  121. #ifdef MLR_MODEL
  122. void dgels_multiple_reg_coeff(double *mpar, double *my, long nn, unsigned ncoeff, unsigned nparameters, double *coeff, unsigned **combinations)
  123. {
  124. /* Arguments */
  125. /* ========= */
  126. /* TRANS (input) CHARACTER*1 */
  127. /* = 'N': the linear system involves A; */
  128. /* = 'T': the linear system involves A**T. */
  129. /* M (input) INTEGER */
  130. /* The number of rows of the matrix A. M >= 0. */
  131. /* N (input) INTEGER */
  132. /* The number of columns of the matrix A. N >= 0. */
  133. /* NRHS (input) INTEGER */
  134. /* The number of right hand sides, i.e., the number of */
  135. /* columns of the matrices B and X. NRHS >=0. */
  136. /* A (input/output) DOUBLE PRECISION array, dimension (LDA,N) */
  137. /* On entry, the M-by-N matrix A. */
  138. /* On exit, */
  139. /* if M >= N, A is overwritten by details of its QR */
  140. /* factorization as returned by DGEQRF; */
  141. /* if M < N, A is overwritten by details of its LQ */
  142. /* factorization as returned by DGELQF. */
  143. /* LDA (input) INTEGER */
  144. /* The leading dimension of the array A. LDA >= max(1,M). */
  145. /* B (input/output) DOUBLE PRECISION array, dimension (LDB,NRHS) */
  146. /* On entry, the matrix B of right hand side vectors, stored */
  147. /* columnwise; B is M-by-NRHS if TRANS = 'N', or N-by-NRHS */
  148. /* if TRANS = 'T'. */
  149. /* On exit, if INFO = 0, B is overwritten by the solution */
  150. /* vectors, stored columnwise: */
  151. /* if TRANS = 'N' and m >= n, rows 1 to n of B contain the least */
  152. /* squares solution vectors; the residual sum of squares for the */
  153. /* solution in each column is given by the sum of squares of */
  154. /* elements N+1 to M in that column; */
  155. /* if TRANS = 'N' and m < n, rows 1 to N of B contain the */
  156. /* minimum norm solution vectors; */
  157. /* if TRANS = 'T' and m >= n, rows 1 to M of B contain the */
  158. /* minimum norm solution vectors; */
  159. /* if TRANS = 'T' and m < n, rows 1 to M of B contain the */
  160. /* least squares solution vectors; the residual sum of squares */
  161. /* for the solution in each column is given by the sum of */
  162. /* squares of elements M+1 to N in that column. */
  163. /* LDB (input) INTEGER */
  164. /* The leading dimension of the array B. LDB >= MAX(1,M,N). */
  165. /* WORK (workspace/output) DOUBLE PRECISION array, dimension (MAX(1,LWORK)) */
  166. /* On exit, if INFO = 0, WORK(1) returns the optimal LWORK. */
  167. /* LWORK (input) INTEGER */
  168. /* The dimension of the array WORK. */
  169. /* LWORK >= max( 1, MN + max( MN, NRHS ) ). */
  170. /* For optimal performance, */
  171. /* LWORK >= max( 1, MN + max( MN, NRHS )*NB ). */
  172. /* where MN = min(M,N) and NB is the optimum block size. */
  173. /* If LWORK = -1, then a workspace query is assumed; the routine */
  174. /* only calculates the optimal size of the WORK array, returns */
  175. /* this value as the first entry of the WORK array, and no error */
  176. /* message related to LWORK is issued by XERBLA. */
  177. /* INFO (output) INTEGER */
  178. /* = 0: successful exit */
  179. /* < 0: if INFO = -i, the i-th argument had an illegal value */
  180. /* > 0: if INFO = i, the i-th diagonal element of the */
  181. /* triangular factor of A is zero, so that A does not have */
  182. /* full rank; the least squares solution could not be */
  183. /* computed. */
  184. /* ===================================================================== */
  185. if(nn <= ncoeff)
  186. {
  187. printf("\nERROR: 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");
  188. exit(1);
  189. }
  190. char trans = 'N';
  191. integer m = nn;
  192. integer n = ncoeff;
  193. integer nrhs = 1; // number of columns of B and X (wich are vectors therefore nrhs=1)
  194. doublereal *X = malloc(sizeof(double)*n*m); // (/!\ modified at the output) contain the model and the different values of pararmters
  195. doublereal *Y = malloc(sizeof(double)*m);
  196. double coefficient;
  197. for (int i=0; i < m; i++)
  198. {
  199. Y[i] = my[i];
  200. X[i*n] = 1.;
  201. for (int j=1; j < n; j++)
  202. {
  203. coefficient = 1.;
  204. for(int k=0; k < nparameters; k++)
  205. {
  206. coefficient *= pow(mpar[i*nparameters+k],combinations[j-1][k]);
  207. }
  208. X[i*n+j] = coefficient;
  209. }
  210. }
  211. integer lda = m;
  212. integer ldb = m; //
  213. integer info;
  214. integer lwork = n*2;
  215. doublereal *work = malloc(sizeof(double)*lwork); // (output)
  216. /* // Running CLAPACK */
  217. dgels_(&trans, &m, &n, &nrhs, X, &lda, Y, &ldb, work, &lwork, &info);
  218. /* Check for the full rank */
  219. if( info != 0 )
  220. {
  221. printf( "Problems with DGELS; info=%ld\n", info);
  222. exit(1);
  223. }
  224. for(int i=0; i<ncoeff; i++)
  225. coeff[i] = Y[i];
  226. free(X);
  227. free(Y);
  228. free(work);
  229. }
  230. #endif //MLR_MODEL
  231. int _starpu_multiple_regression(struct starpu_perfmodel_history_list *ptr, double *coeff, unsigned ncoeff, unsigned nparameters, unsigned **combinations, char *codelet_name)
  232. {
  233. // Computing number of rows
  234. long n=find_long_list_size(ptr);
  235. STARPU_ASSERT(n);
  236. // Reading old calibrations if necessary
  237. FILE *f;
  238. char filepath[50];
  239. snprintf(filepath, 50, "/tmp/%s.out", codelet_name);
  240. long old_lines=0;
  241. int calibrate = starpu_get_env_number("STARPU_CALIBRATE");
  242. if (calibrate==1)
  243. {
  244. f = fopen(filepath, "a+");
  245. STARPU_ASSERT_MSG(f, "Could not save performance model %s\n", filepath);
  246. old_lines=count_file_lines(f);
  247. STARPU_ASSERT(old_lines);
  248. n+=old_lines;
  249. }
  250. // Allocating X and Y matrices
  251. double *mpar = (double *) malloc(nparameters*n*sizeof(double));
  252. STARPU_ASSERT(mpar);
  253. double *my = (double *) malloc(n*sizeof(double));
  254. STARPU_ASSERT(my);
  255. // Loading old calibration
  256. if (calibrate==1)
  257. load_old_calibration(mpar, my, nparameters, f);
  258. // Filling X and Y matrices with measured values
  259. dump_multiple_regression_list(mpar, my, old_lines, nparameters, ptr);
  260. // Computing coefficients using multiple linear regression
  261. #ifdef TESTGSL
  262. gsl_multiple_reg_coeff(mpar, my, n, ncoeff, nparameters, coeff, combinations);
  263. #endif
  264. #ifdef MLR_MODEL
  265. dgels_multiple_reg_coeff(mpar, my, n, ncoeff, nparameters, coeff, combinations);
  266. #endif //MLR_MODEL
  267. // Preparing new output calibration file
  268. if (calibrate==2)
  269. {
  270. f = fopen(filepath, "w+");
  271. STARPU_ASSERT_MSG(f, "Could not save performance model %s\n", filepath);
  272. fprintf(f, "Duration");
  273. for(int k=0; k < nparameters; k++)
  274. {
  275. fprintf(f, ", P%d", k);
  276. }
  277. }
  278. // Writing parameters to calibration file
  279. if (calibrate>0)
  280. {
  281. for(int i=old_lines; i<n; i++)
  282. {
  283. fprintf(f, "\n%f", my[i]);
  284. for(int j=0; j<nparameters;j++)
  285. fprintf(f, ", %f", mpar[i*nparameters+j]);
  286. }
  287. fclose(f);
  288. }
  289. // Cleanup
  290. free(mpar);
  291. free(my);
  292. return 0;
  293. }