gsl_linalg__svd.c 17 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630
  1. /* linalg/svd.c
  2. *
  3. * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2004, 2007 Gerard Jungman, 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 <stdlib.h>
  21. #include <string.h>
  22. #include "gsl_math.h"
  23. #include "gsl_vector.h"
  24. #include "gsl_matrix.h"
  25. #include "gsl_blas.h"
  26. #include "gsl_linalg.h"
  27. #include "gsl_linalg__givens.c"
  28. #include "gsl_linalg__svdstep.c"
  29. /* Factorise a general M x N matrix A into,
  30. *
  31. * A = U D V^T
  32. *
  33. * where U is a column-orthogonal M x N matrix (U^T U = I),
  34. * D is a diagonal N x N matrix,
  35. * and V is an N x N orthogonal matrix (V^T V = V V^T = I)
  36. *
  37. * U is stored in the original matrix A, which has the same size
  38. *
  39. * V is stored as a separate matrix (not V^T). You must take the
  40. * transpose to form the product above.
  41. *
  42. * The diagonal matrix D is stored in the vector S, D_ii = S_i
  43. */
  44. int
  45. gsl_linalg_SV_decomp (gsl_matrix * A, gsl_matrix * V, gsl_vector * S,
  46. gsl_vector * work)
  47. {
  48. size_t a, b, i, j, iter;
  49. const size_t M = A->size1;
  50. const size_t N = A->size2;
  51. const size_t K = GSL_MIN (M, N);
  52. if (M < N)
  53. {
  54. GSL_ERROR ("svd of MxN matrix, M<N, is not implemented", GSL_EUNIMPL);
  55. }
  56. else if (V->size1 != N)
  57. {
  58. GSL_ERROR ("square matrix V must match second dimension of matrix A",
  59. GSL_EBADLEN);
  60. }
  61. else if (V->size1 != V->size2)
  62. {
  63. GSL_ERROR ("matrix V must be square", GSL_ENOTSQR);
  64. }
  65. else if (S->size != N)
  66. {
  67. GSL_ERROR ("length of vector S must match second dimension of matrix A",
  68. GSL_EBADLEN);
  69. }
  70. else if (work->size != N)
  71. {
  72. GSL_ERROR ("length of workspace must match second dimension of matrix A",
  73. GSL_EBADLEN);
  74. }
  75. /* Handle the case of N = 1 (SVD of a column vector) */
  76. if (N == 1)
  77. {
  78. gsl_vector_view column = gsl_matrix_column (A, 0);
  79. double norm = gsl_blas_dnrm2 (&column.vector);
  80. gsl_vector_set (S, 0, norm);
  81. gsl_matrix_set (V, 0, 0, 1.0);
  82. if (norm != 0.0)
  83. {
  84. gsl_blas_dscal (1.0/norm, &column.vector);
  85. }
  86. return GSL_SUCCESS;
  87. }
  88. {
  89. gsl_vector_view f = gsl_vector_subvector (work, 0, K - 1);
  90. /* bidiagonalize matrix A, unpack A into U S V */
  91. gsl_linalg_bidiag_decomp (A, S, &f.vector);
  92. gsl_linalg_bidiag_unpack2 (A, S, &f.vector, V);
  93. /* apply reduction steps to B=(S,Sd) */
  94. chop_small_elements (S, &f.vector);
  95. /* Progressively reduce the matrix until it is diagonal */
  96. b = N - 1;
  97. iter = 0;
  98. while (b > 0)
  99. {
  100. double fbm1 = gsl_vector_get (&f.vector, b - 1);
  101. if (fbm1 == 0.0 || gsl_isnan (fbm1))
  102. {
  103. b--;
  104. continue;
  105. }
  106. /* Find the largest unreduced block (a,b) starting from b
  107. and working backwards */
  108. a = b - 1;
  109. while (a > 0)
  110. {
  111. double fam1 = gsl_vector_get (&f.vector, a - 1);
  112. if (fam1 == 0.0 || gsl_isnan (fam1))
  113. {
  114. break;
  115. }
  116. a--;
  117. }
  118. iter++;
  119. if (iter > 100 * N)
  120. {
  121. GSL_ERROR("SVD decomposition failed to converge", GSL_EMAXITER);
  122. }
  123. {
  124. const size_t n_block = b - a + 1;
  125. gsl_vector_view S_block = gsl_vector_subvector (S, a, n_block);
  126. gsl_vector_view f_block = gsl_vector_subvector (&f.vector, a, n_block - 1);
  127. gsl_matrix_view U_block =
  128. gsl_matrix_submatrix (A, 0, a, A->size1, n_block);
  129. gsl_matrix_view V_block =
  130. gsl_matrix_submatrix (V, 0, a, V->size1, n_block);
  131. qrstep (&S_block.vector, &f_block.vector, &U_block.matrix, &V_block.matrix);
  132. /* remove any small off-diagonal elements */
  133. chop_small_elements (&S_block.vector, &f_block.vector);
  134. }
  135. }
  136. }
  137. /* Make singular values positive by reflections if necessary */
  138. for (j = 0; j < K; j++)
  139. {
  140. double Sj = gsl_vector_get (S, j);
  141. if (Sj < 0.0)
  142. {
  143. for (i = 0; i < N; i++)
  144. {
  145. double Vij = gsl_matrix_get (V, i, j);
  146. gsl_matrix_set (V, i, j, -Vij);
  147. }
  148. gsl_vector_set (S, j, -Sj);
  149. }
  150. }
  151. /* Sort singular values into decreasing order */
  152. for (i = 0; i < K; i++)
  153. {
  154. double S_max = gsl_vector_get (S, i);
  155. size_t i_max = i;
  156. for (j = i + 1; j < K; j++)
  157. {
  158. double Sj = gsl_vector_get (S, j);
  159. if (Sj > S_max)
  160. {
  161. S_max = Sj;
  162. i_max = j;
  163. }
  164. }
  165. if (i_max != i)
  166. {
  167. /* swap eigenvalues */
  168. gsl_vector_swap_elements (S, i, i_max);
  169. /* swap eigenvectors */
  170. gsl_matrix_swap_columns (A, i, i_max);
  171. gsl_matrix_swap_columns (V, i, i_max);
  172. }
  173. }
  174. return GSL_SUCCESS;
  175. }
  176. /* Modified algorithm which is better for M>>N */
  177. int
  178. gsl_linalg_SV_decomp_mod (gsl_matrix * A,
  179. gsl_matrix * X,
  180. gsl_matrix * V, gsl_vector * S, gsl_vector * work)
  181. {
  182. size_t i, j;
  183. const size_t M = A->size1;
  184. const size_t N = A->size2;
  185. if (M < N)
  186. {
  187. GSL_ERROR ("svd of MxN matrix, M<N, is not implemented", GSL_EUNIMPL);
  188. }
  189. else if (V->size1 != N)
  190. {
  191. GSL_ERROR ("square matrix V must match second dimension of matrix A",
  192. GSL_EBADLEN);
  193. }
  194. else if (V->size1 != V->size2)
  195. {
  196. GSL_ERROR ("matrix V must be square", GSL_ENOTSQR);
  197. }
  198. else if (X->size1 != N)
  199. {
  200. GSL_ERROR ("square matrix X must match second dimension of matrix A",
  201. GSL_EBADLEN);
  202. }
  203. else if (X->size1 != X->size2)
  204. {
  205. GSL_ERROR ("matrix X must be square", GSL_ENOTSQR);
  206. }
  207. else if (S->size != N)
  208. {
  209. GSL_ERROR ("length of vector S must match second dimension of matrix A",
  210. GSL_EBADLEN);
  211. }
  212. else if (work->size != N)
  213. {
  214. GSL_ERROR ("length of workspace must match second dimension of matrix A",
  215. GSL_EBADLEN);
  216. }
  217. if (N == 1)
  218. {
  219. gsl_vector_view column = gsl_matrix_column (A, 0);
  220. double norm = gsl_blas_dnrm2 (&column.vector);
  221. gsl_vector_set (S, 0, norm);
  222. gsl_matrix_set (V, 0, 0, 1.0);
  223. if (norm != 0.0)
  224. {
  225. gsl_blas_dscal (1.0/norm, &column.vector);
  226. }
  227. return GSL_SUCCESS;
  228. }
  229. /* Convert A into an upper triangular matrix R */
  230. for (i = 0; i < N; i++)
  231. {
  232. gsl_vector_view c = gsl_matrix_column (A, i);
  233. gsl_vector_view v = gsl_vector_subvector (&c.vector, i, M - i);
  234. double tau_i = gsl_linalg_householder_transform (&v.vector);
  235. /* Apply the transformation to the remaining columns */
  236. if (i + 1 < N)
  237. {
  238. gsl_matrix_view m =
  239. gsl_matrix_submatrix (A, i, i + 1, M - i, N - (i + 1));
  240. gsl_linalg_householder_hm (tau_i, &v.vector, &m.matrix);
  241. }
  242. gsl_vector_set (S, i, tau_i);
  243. }
  244. /* Copy the upper triangular part of A into X */
  245. for (i = 0; i < N; i++)
  246. {
  247. for (j = 0; j < i; j++)
  248. {
  249. gsl_matrix_set (X, i, j, 0.0);
  250. }
  251. {
  252. double Aii = gsl_matrix_get (A, i, i);
  253. gsl_matrix_set (X, i, i, Aii);
  254. }
  255. for (j = i + 1; j < N; j++)
  256. {
  257. double Aij = gsl_matrix_get (A, i, j);
  258. gsl_matrix_set (X, i, j, Aij);
  259. }
  260. }
  261. /* Convert A into an orthogonal matrix L */
  262. for (j = N; j > 0 && j--;)
  263. {
  264. /* Householder column transformation to accumulate L */
  265. double tj = gsl_vector_get (S, j);
  266. gsl_matrix_view m = gsl_matrix_submatrix (A, j, j, M - j, N - j);
  267. gsl_linalg_householder_hm1 (tj, &m.matrix);
  268. }
  269. /* unpack R into X V S */
  270. gsl_linalg_SV_decomp (X, V, S, work);
  271. /* Multiply L by X, to obtain U = L X, stored in U */
  272. {
  273. gsl_vector_view sum = gsl_vector_subvector (work, 0, N);
  274. for (i = 0; i < M; i++)
  275. {
  276. gsl_vector_view L_i = gsl_matrix_row (A, i);
  277. gsl_vector_set_zero (&sum.vector);
  278. for (j = 0; j < N; j++)
  279. {
  280. double Lij = gsl_vector_get (&L_i.vector, j);
  281. gsl_vector_view X_j = gsl_matrix_row (X, j);
  282. gsl_blas_daxpy (Lij, &X_j.vector, &sum.vector);
  283. }
  284. gsl_vector_memcpy (&L_i.vector, &sum.vector);
  285. }
  286. }
  287. return GSL_SUCCESS;
  288. }
  289. /* Solves the system A x = b using the SVD factorization
  290. *
  291. * A = U S V^T
  292. *
  293. * to obtain x. For M x N systems it finds the solution in the least
  294. * squares sense.
  295. */
  296. int
  297. gsl_linalg_SV_solve (const gsl_matrix * U,
  298. const gsl_matrix * V,
  299. const gsl_vector * S,
  300. const gsl_vector * b, gsl_vector * x)
  301. {
  302. if (U->size1 != b->size)
  303. {
  304. GSL_ERROR ("first dimension of matrix U must size of vector b",
  305. GSL_EBADLEN);
  306. }
  307. else if (U->size2 != S->size)
  308. {
  309. GSL_ERROR ("length of vector S must match second dimension of matrix U",
  310. GSL_EBADLEN);
  311. }
  312. else if (V->size1 != V->size2)
  313. {
  314. GSL_ERROR ("matrix V must be square", GSL_ENOTSQR);
  315. }
  316. else if (S->size != V->size1)
  317. {
  318. GSL_ERROR ("length of vector S must match size of matrix V",
  319. GSL_EBADLEN);
  320. }
  321. else if (V->size2 != x->size)
  322. {
  323. GSL_ERROR ("size of matrix V must match size of vector x", GSL_EBADLEN);
  324. }
  325. else
  326. {
  327. const size_t N = U->size2;
  328. size_t i;
  329. gsl_vector *w = gsl_vector_calloc (N);
  330. gsl_blas_dgemv (CblasTrans, 1.0, U, b, 0.0, w);
  331. for (i = 0; i < N; i++)
  332. {
  333. double wi = gsl_vector_get (w, i);
  334. double alpha = gsl_vector_get (S, i);
  335. if (alpha != 0)
  336. alpha = 1.0 / alpha;
  337. gsl_vector_set (w, i, alpha * wi);
  338. }
  339. gsl_blas_dgemv (CblasNoTrans, 1.0, V, w, 0.0, x);
  340. gsl_vector_free (w);
  341. return GSL_SUCCESS;
  342. }
  343. }
  344. /* This is a the jacobi version */
  345. /* Author: G. Jungman */
  346. /*
  347. * Algorithm due to J.C. Nash, Compact Numerical Methods for
  348. * Computers (New York: Wiley and Sons, 1979), chapter 3.
  349. * See also Algorithm 4.1 in
  350. * James Demmel, Kresimir Veselic, "Jacobi's Method is more
  351. * accurate than QR", Lapack Working Note 15 (LAWN15), October 1989.
  352. * Available from netlib.
  353. *
  354. * Based on code by Arthur Kosowsky, Rutgers University
  355. * kosowsky@physics.rutgers.edu
  356. *
  357. * Another relevant paper is, P.P.M. De Rijk, "A One-Sided Jacobi
  358. * Algorithm for computing the singular value decomposition on a
  359. * vector computer", SIAM Journal of Scientific and Statistical
  360. * Computing, Vol 10, No 2, pp 359-371, March 1989.
  361. *
  362. */
  363. int
  364. gsl_linalg_SV_decomp_jacobi (gsl_matrix * A, gsl_matrix * Q, gsl_vector * S)
  365. {
  366. if (A->size1 < A->size2)
  367. {
  368. /* FIXME: only implemented M>=N case so far */
  369. GSL_ERROR ("svd of MxN matrix, M<N, is not implemented", GSL_EUNIMPL);
  370. }
  371. else if (Q->size1 != A->size2)
  372. {
  373. GSL_ERROR ("square matrix Q must match second dimension of matrix A",
  374. GSL_EBADLEN);
  375. }
  376. else if (Q->size1 != Q->size2)
  377. {
  378. GSL_ERROR ("matrix Q must be square", GSL_ENOTSQR);
  379. }
  380. else if (S->size != A->size2)
  381. {
  382. GSL_ERROR ("length of vector S must match second dimension of matrix A",
  383. GSL_EBADLEN);
  384. }
  385. else
  386. {
  387. const size_t M = A->size1;
  388. const size_t N = A->size2;
  389. size_t i, j, k;
  390. /* Initialize the rotation counter and the sweep counter. */
  391. int count = 1;
  392. int sweep = 0;
  393. int sweepmax = 5*N;
  394. double tolerance = 10 * M * GSL_DBL_EPSILON;
  395. /* Always do at least 12 sweeps. */
  396. sweepmax = GSL_MAX (sweepmax, 12);
  397. /* Set Q to the identity matrix. */
  398. gsl_matrix_set_identity (Q);
  399. /* Store the column error estimates in S, for use during the
  400. orthogonalization */
  401. for (j = 0; j < N; j++)
  402. {
  403. gsl_vector_view cj = gsl_matrix_column (A, j);
  404. double sj = gsl_blas_dnrm2 (&cj.vector);
  405. gsl_vector_set(S, j, GSL_DBL_EPSILON * sj);
  406. }
  407. /* Orthogonalize A by plane rotations. */
  408. while (count > 0 && sweep <= sweepmax)
  409. {
  410. /* Initialize rotation counter. */
  411. count = N * (N - 1) / 2;
  412. for (j = 0; j < N - 1; j++)
  413. {
  414. for (k = j + 1; k < N; k++)
  415. {
  416. double a = 0.0;
  417. double b = 0.0;
  418. double p = 0.0;
  419. double q = 0.0;
  420. double cosine, sine;
  421. double v;
  422. double abserr_a, abserr_b;
  423. int sorted, orthog, noisya, noisyb;
  424. gsl_vector_view cj = gsl_matrix_column (A, j);
  425. gsl_vector_view ck = gsl_matrix_column (A, k);
  426. gsl_blas_ddot (&cj.vector, &ck.vector, &p);
  427. p *= 2.0 ; /* equation 9a: p = 2 x.y */
  428. a = gsl_blas_dnrm2 (&cj.vector);
  429. b = gsl_blas_dnrm2 (&ck.vector);
  430. q = a * a - b * b;
  431. v = hypot(p, q);
  432. /* test for columns j,k orthogonal, or dominant errors */
  433. abserr_a = gsl_vector_get(S,j);
  434. abserr_b = gsl_vector_get(S,k);
  435. sorted = (gsl_coerce_double(a) >= gsl_coerce_double(b));
  436. orthog = (fabs (p) <= tolerance * gsl_coerce_double(a * b));
  437. noisya = (a < abserr_a);
  438. noisyb = (b < abserr_b);
  439. if (sorted && (orthog || noisya || noisyb))
  440. {
  441. count--;
  442. continue;
  443. }
  444. /* calculate rotation angles */
  445. if (v == 0 || !sorted)
  446. {
  447. cosine = 0.0;
  448. sine = 1.0;
  449. }
  450. else
  451. {
  452. cosine = sqrt((v + q) / (2.0 * v));
  453. sine = p / (2.0 * v * cosine);
  454. }
  455. /* apply rotation to A */
  456. for (i = 0; i < M; i++)
  457. {
  458. const double Aik = gsl_matrix_get (A, i, k);
  459. const double Aij = gsl_matrix_get (A, i, j);
  460. gsl_matrix_set (A, i, j, Aij * cosine + Aik * sine);
  461. gsl_matrix_set (A, i, k, -Aij * sine + Aik * cosine);
  462. }
  463. gsl_vector_set(S, j, fabs(cosine) * abserr_a + fabs(sine) * abserr_b);
  464. gsl_vector_set(S, k, fabs(sine) * abserr_a + fabs(cosine) * abserr_b);
  465. /* apply rotation to Q */
  466. for (i = 0; i < N; i++)
  467. {
  468. const double Qij = gsl_matrix_get (Q, i, j);
  469. const double Qik = gsl_matrix_get (Q, i, k);
  470. gsl_matrix_set (Q, i, j, Qij * cosine + Qik * sine);
  471. gsl_matrix_set (Q, i, k, -Qij * sine + Qik * cosine);
  472. }
  473. }
  474. }
  475. /* Sweep completed. */
  476. sweep++;
  477. }
  478. /*
  479. * Orthogonalization complete. Compute singular values.
  480. */
  481. {
  482. double prev_norm = -1.0;
  483. for (j = 0; j < N; j++)
  484. {
  485. gsl_vector_view column = gsl_matrix_column (A, j);
  486. double norm = gsl_blas_dnrm2 (&column.vector);
  487. /* Determine if singular value is zero, according to the
  488. criteria used in the main loop above (i.e. comparison
  489. with norm of previous column). */
  490. if (norm == 0.0 || prev_norm == 0.0
  491. || (j > 0 && norm <= tolerance * prev_norm))
  492. {
  493. gsl_vector_set (S, j, 0.0); /* singular */
  494. gsl_vector_set_zero (&column.vector); /* annihilate column */
  495. prev_norm = 0.0;
  496. }
  497. else
  498. {
  499. gsl_vector_set (S, j, norm); /* non-singular */
  500. gsl_vector_scale (&column.vector, 1.0 / norm); /* normalize column */
  501. prev_norm = norm;
  502. }
  503. }
  504. }
  505. if (count > 0)
  506. {
  507. /* reached sweep limit */
  508. GSL_ERROR ("Jacobi iterations did not reach desired tolerance",
  509. GSL_ETOL);
  510. }
  511. return GSL_SUCCESS;
  512. }
  513. }