gsl_multifit__multilinear.c 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433
  1. /* multifit/multilinear.c
  2. *
  3. * Copyright (C) 2000, 2007 Brian Gough
  4. *
  5. * This program is free software; you can redistribute it and/or modify
  6. * it under the terms of the GNU General Public License as published by
  7. * the Free Software Foundation; either version 3 of the License, or (at
  8. * your option) any later version.
  9. *
  10. * This program is distributed in the hope that it will be useful, but
  11. * WITHOUT ANY WARRANTY; without even the implied warranty of
  12. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
  13. * General Public License for more details.
  14. *
  15. * You should have received a copy of the GNU General Public License
  16. * along with this program; if not, write to the Free Software
  17. * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
  18. */
  19. #include "gsl__config.h"
  20. #include "gsl_errno.h"
  21. #include "gsl_multifit.h"
  22. #include "gsl_blas.h"
  23. #include "gsl_linalg.h"
  24. /* Fit
  25. *
  26. * y = X c
  27. *
  28. * where X is an M x N matrix of M observations for N variables.
  29. *
  30. */
  31. int
  32. gsl_multifit_linear (const gsl_matrix * X,
  33. const gsl_vector * y,
  34. gsl_vector * c,
  35. gsl_matrix * cov,
  36. double *chisq, gsl_multifit_linear_workspace * work)
  37. {
  38. size_t rank;
  39. int status = gsl_multifit_linear_svd (X, y, GSL_DBL_EPSILON, &rank, c,
  40. cov, chisq, work);
  41. return status;
  42. }
  43. /* Handle the general case of the SVD with tolerance and rank */
  44. int
  45. gsl_multifit_linear_svd (const gsl_matrix * X,
  46. const gsl_vector * y,
  47. double tol,
  48. size_t * rank,
  49. gsl_vector * c,
  50. gsl_matrix * cov,
  51. double *chisq, gsl_multifit_linear_workspace * work)
  52. {
  53. if (X->size1 != y->size)
  54. {
  55. GSL_ERROR
  56. ("number of observations in y does not match rows of matrix X",
  57. GSL_EBADLEN);
  58. }
  59. else if (X->size2 != c->size)
  60. {
  61. GSL_ERROR ("number of parameters c does not match columns of matrix X",
  62. GSL_EBADLEN);
  63. }
  64. else if (cov->size1 != cov->size2)
  65. {
  66. GSL_ERROR ("covariance matrix is not square", GSL_ENOTSQR);
  67. }
  68. else if (c->size != cov->size1)
  69. {
  70. GSL_ERROR
  71. ("number of parameters does not match size of covariance matrix",
  72. GSL_EBADLEN);
  73. }
  74. else if (X->size1 != work->n || X->size2 != work->p)
  75. {
  76. GSL_ERROR
  77. ("size of workspace does not match size of observation matrix",
  78. GSL_EBADLEN);
  79. }
  80. else if (tol <= 0)
  81. {
  82. GSL_ERROR ("tolerance must be positive", GSL_EINVAL);
  83. }
  84. else
  85. {
  86. const size_t n = X->size1;
  87. const size_t p = X->size2;
  88. size_t i, j, p_eff;
  89. gsl_matrix *A = work->A;
  90. gsl_matrix *Q = work->Q;
  91. gsl_matrix *QSI = work->QSI;
  92. gsl_vector *S = work->S;
  93. gsl_vector *xt = work->xt;
  94. gsl_vector *D = work->D;
  95. /* Copy X to workspace, A <= X */
  96. gsl_matrix_memcpy (A, X);
  97. /* Balance the columns of the matrix A */
  98. gsl_linalg_balance_columns (A, D);
  99. /* Decompose A into U S Q^T */
  100. gsl_linalg_SV_decomp_mod (A, QSI, Q, S, xt);
  101. /* Solve y = A c for c */
  102. gsl_blas_dgemv (CblasTrans, 1.0, A, y, 0.0, xt);
  103. /* Scale the matrix Q, Q' = Q S^-1 */
  104. gsl_matrix_memcpy (QSI, Q);
  105. {
  106. double alpha0 = gsl_vector_get (S, 0);
  107. p_eff = 0;
  108. for (j = 0; j < p; j++)
  109. {
  110. gsl_vector_view column = gsl_matrix_column (QSI, j);
  111. double alpha = gsl_vector_get (S, j);
  112. if (alpha <= tol * alpha0) {
  113. alpha = 0.0;
  114. } else {
  115. alpha = 1.0 / alpha;
  116. p_eff++;
  117. }
  118. gsl_vector_scale (&column.vector, alpha);
  119. }
  120. *rank = p_eff;
  121. }
  122. gsl_vector_set_zero (c);
  123. gsl_blas_dgemv (CblasNoTrans, 1.0, QSI, xt, 0.0, c);
  124. /* Unscale the balancing factors */
  125. gsl_vector_div (c, D);
  126. /* Compute chisq, from residual r = y - X c */
  127. {
  128. double s2 = 0, r2 = 0;
  129. for (i = 0; i < n; i++)
  130. {
  131. double yi = gsl_vector_get (y, i);
  132. gsl_vector_const_view row = gsl_matrix_const_row (X, i);
  133. double y_est, ri;
  134. gsl_blas_ddot (&row.vector, c, &y_est);
  135. ri = yi - y_est;
  136. r2 += ri * ri;
  137. }
  138. s2 = r2 / (n - p_eff); /* p_eff == rank */
  139. *chisq = r2;
  140. /* Form variance-covariance matrix cov = s2 * (Q S^-1) (Q S^-1)^T */
  141. for (i = 0; i < p; i++)
  142. {
  143. gsl_vector_view row_i = gsl_matrix_row (QSI, i);
  144. double d_i = gsl_vector_get (D, i);
  145. for (j = i; j < p; j++)
  146. {
  147. gsl_vector_view row_j = gsl_matrix_row (QSI, j);
  148. double d_j = gsl_vector_get (D, j);
  149. double s;
  150. gsl_blas_ddot (&row_i.vector, &row_j.vector, &s);
  151. gsl_matrix_set (cov, i, j, s * s2 / (d_i * d_j));
  152. gsl_matrix_set (cov, j, i, s * s2 / (d_i * d_j));
  153. }
  154. }
  155. }
  156. return GSL_SUCCESS;
  157. }
  158. }
  159. int
  160. gsl_multifit_wlinear (const gsl_matrix * X,
  161. const gsl_vector * w,
  162. const gsl_vector * y,
  163. gsl_vector * c,
  164. gsl_matrix * cov,
  165. double *chisq, gsl_multifit_linear_workspace * work)
  166. {
  167. size_t rank;
  168. int status = gsl_multifit_wlinear_svd (X, w, y, GSL_DBL_EPSILON, &rank, c,
  169. cov, chisq, work);
  170. return status;
  171. }
  172. int
  173. gsl_multifit_wlinear_svd (const gsl_matrix * X,
  174. const gsl_vector * w,
  175. const gsl_vector * y,
  176. double tol,
  177. size_t * rank,
  178. gsl_vector * c,
  179. gsl_matrix * cov,
  180. double *chisq, gsl_multifit_linear_workspace * work)
  181. {
  182. if (X->size1 != y->size)
  183. {
  184. GSL_ERROR
  185. ("number of observations in y does not match rows of matrix X",
  186. GSL_EBADLEN);
  187. }
  188. else if (X->size2 != c->size)
  189. {
  190. GSL_ERROR ("number of parameters c does not match columns of matrix X",
  191. GSL_EBADLEN);
  192. }
  193. else if (w->size != y->size)
  194. {
  195. GSL_ERROR ("number of weights does not match number of observations",
  196. GSL_EBADLEN);
  197. }
  198. else if (cov->size1 != cov->size2)
  199. {
  200. GSL_ERROR ("covariance matrix is not square", GSL_ENOTSQR);
  201. }
  202. else if (c->size != cov->size1)
  203. {
  204. GSL_ERROR
  205. ("number of parameters does not match size of covariance matrix",
  206. GSL_EBADLEN);
  207. }
  208. else if (X->size1 != work->n || X->size2 != work->p)
  209. {
  210. GSL_ERROR
  211. ("size of workspace does not match size of observation matrix",
  212. GSL_EBADLEN);
  213. }
  214. else
  215. {
  216. const size_t n = X->size1;
  217. const size_t p = X->size2;
  218. size_t i, j, p_eff;
  219. gsl_matrix *A = work->A;
  220. gsl_matrix *Q = work->Q;
  221. gsl_matrix *QSI = work->QSI;
  222. gsl_vector *S = work->S;
  223. gsl_vector *t = work->t;
  224. gsl_vector *xt = work->xt;
  225. gsl_vector *D = work->D;
  226. /* Scale X, A = sqrt(w) X */
  227. gsl_matrix_memcpy (A, X);
  228. for (i = 0; i < n; i++)
  229. {
  230. double wi = gsl_vector_get (w, i);
  231. if (wi < 0)
  232. wi = 0;
  233. {
  234. gsl_vector_view row = gsl_matrix_row (A, i);
  235. gsl_vector_scale (&row.vector, sqrt (wi));
  236. }
  237. }
  238. /* Balance the columns of the matrix A */
  239. gsl_linalg_balance_columns (A, D);
  240. /* Decompose A into U S Q^T */
  241. gsl_linalg_SV_decomp_mod (A, QSI, Q, S, xt);
  242. /* Solve sqrt(w) y = A c for c, by first computing t = sqrt(w) y */
  243. for (i = 0; i < n; i++)
  244. {
  245. double wi = gsl_vector_get (w, i);
  246. double yi = gsl_vector_get (y, i);
  247. if (wi < 0)
  248. wi = 0;
  249. gsl_vector_set (t, i, sqrt (wi) * yi);
  250. }
  251. gsl_blas_dgemv (CblasTrans, 1.0, A, t, 0.0, xt);
  252. /* Scale the matrix Q, Q' = Q S^-1 */
  253. gsl_matrix_memcpy (QSI, Q);
  254. {
  255. double alpha0 = gsl_vector_get (S, 0);
  256. p_eff = 0;
  257. for (j = 0; j < p; j++)
  258. {
  259. gsl_vector_view column = gsl_matrix_column (QSI, j);
  260. double alpha = gsl_vector_get (S, j);
  261. if (alpha <= tol * alpha0) {
  262. alpha = 0.0;
  263. } else {
  264. alpha = 1.0 / alpha;
  265. p_eff++;
  266. }
  267. gsl_vector_scale (&column.vector, alpha);
  268. }
  269. *rank = p_eff;
  270. }
  271. gsl_vector_set_zero (c);
  272. /* Solution */
  273. gsl_blas_dgemv (CblasNoTrans, 1.0, QSI, xt, 0.0, c);
  274. /* Unscale the balancing factors */
  275. gsl_vector_div (c, D);
  276. /* Form covariance matrix cov = (Q S^-1) (Q S^-1)^T */
  277. for (i = 0; i < p; i++)
  278. {
  279. gsl_vector_view row_i = gsl_matrix_row (QSI, i);
  280. double d_i = gsl_vector_get (D, i);
  281. for (j = i; j < p; j++)
  282. {
  283. gsl_vector_view row_j = gsl_matrix_row (QSI, j);
  284. double d_j = gsl_vector_get (D, j);
  285. double s;
  286. gsl_blas_ddot (&row_i.vector, &row_j.vector, &s);
  287. gsl_matrix_set (cov, i, j, s / (d_i * d_j));
  288. gsl_matrix_set (cov, j, i, s / (d_i * d_j));
  289. }
  290. }
  291. /* Compute chisq, from residual r = y - X c */
  292. {
  293. double r2 = 0;
  294. for (i = 0; i < n; i++)
  295. {
  296. double yi = gsl_vector_get (y, i);
  297. double wi = gsl_vector_get (w, i);
  298. gsl_vector_const_view row = gsl_matrix_const_row (X, i);
  299. double y_est, ri;
  300. gsl_blas_ddot (&row.vector, c, &y_est);
  301. ri = yi - y_est;
  302. r2 += wi * ri * ri;
  303. }
  304. *chisq = r2;
  305. }
  306. return GSL_SUCCESS;
  307. }
  308. }
  309. int
  310. gsl_multifit_linear_est (const gsl_vector * x,
  311. const gsl_vector * c,
  312. const gsl_matrix * cov, double *y, double *y_err)
  313. {
  314. if (x->size != c->size)
  315. {
  316. GSL_ERROR ("number of parameters c does not match number of observations x",
  317. GSL_EBADLEN);
  318. }
  319. else if (cov->size1 != cov->size2)
  320. {
  321. GSL_ERROR ("covariance matrix is not square", GSL_ENOTSQR);
  322. }
  323. else if (c->size != cov->size1)
  324. {
  325. GSL_ERROR ("number of parameters c does not match size of covariance matrix cov",
  326. GSL_EBADLEN);
  327. }
  328. else
  329. {
  330. size_t i, j;
  331. double var = 0;
  332. gsl_blas_ddot(x, c, y); /* y = x.c */
  333. /* var = x' cov x */
  334. for (i = 0; i < x->size; i++)
  335. {
  336. const double xi = gsl_vector_get (x, i);
  337. var += xi * xi * gsl_matrix_get (cov, i, i);
  338. for (j = 0; j < i; j++)
  339. {
  340. const double xj = gsl_vector_get (x, j);
  341. var += 2 * xi * xj * gsl_matrix_get (cov, i, j);
  342. }
  343. }
  344. *y_err = sqrt (var);
  345. return GSL_SUCCESS;
  346. }
  347. }