gsl_linalg__cholesky.c 6.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267
  1. /* Cholesky Decomposition
  2. *
  3. * Copyright (C) 2000 Thomas Walter
  4. *
  5. * 03 May 2000: Modified for GSL by Brian Gough
  6. * 29 Jul 2005: Additions by Gerard Jungman
  7. * Copyright (C) 2000,2001, 2002, 2003, 2005, 2007 Brian Gough, Gerard Jungman
  8. *
  9. * This is free software; you can redistribute it and/or modify it
  10. * under the terms of the GNU General Public License as published by the
  11. * Free Software Foundation; either version 3, or (at your option) any
  12. * later version.
  13. *
  14. * This source is distributed in the hope that it will be useful, but WITHOUT
  15. * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
  16. * FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
  17. * for more details.
  18. */
  19. /*
  20. * Cholesky decomposition of a symmetrix positive definite matrix.
  21. * This is useful to solve the matrix arising in
  22. * periodic cubic splines
  23. * approximating splines
  24. *
  25. * This algorithm does:
  26. * A = L * L'
  27. * with
  28. * L := lower left triangle matrix
  29. * L' := the transposed form of L.
  30. *
  31. */
  32. #include "gsl__config.h"
  33. #include "gsl_math.h"
  34. #include "gsl_errno.h"
  35. #include "gsl_vector.h"
  36. #include "gsl_matrix.h"
  37. #include "gsl_blas.h"
  38. #include "gsl_linalg.h"
  39. static inline
  40. double
  41. quiet_sqrt (double x)
  42. /* avoids runtime error, for checking matrix for positive definiteness */
  43. {
  44. return (x >= 0) ? sqrt(x) : GSL_NAN;
  45. }
  46. int
  47. gsl_linalg_cholesky_decomp (gsl_matrix * A)
  48. {
  49. const size_t M = A->size1;
  50. const size_t N = A->size2;
  51. if (M != N)
  52. {
  53. GSL_ERROR("cholesky decomposition requires square matrix", GSL_ENOTSQR);
  54. }
  55. else
  56. {
  57. size_t i,j,k;
  58. int status = 0;
  59. /* Do the first 2 rows explicitly. It is simple, and faster. And
  60. * one can return if the matrix has only 1 or 2 rows.
  61. */
  62. double A_00 = gsl_matrix_get (A, 0, 0);
  63. double L_00 = quiet_sqrt(A_00);
  64. if (A_00 <= 0)
  65. {
  66. status = GSL_EDOM ;
  67. }
  68. gsl_matrix_set (A, 0, 0, L_00);
  69. if (M > 1)
  70. {
  71. double A_10 = gsl_matrix_get (A, 1, 0);
  72. double A_11 = gsl_matrix_get (A, 1, 1);
  73. double L_10 = A_10 / L_00;
  74. double diag = A_11 - L_10 * L_10;
  75. double L_11 = quiet_sqrt(diag);
  76. if (diag <= 0)
  77. {
  78. status = GSL_EDOM;
  79. }
  80. gsl_matrix_set (A, 1, 0, L_10);
  81. gsl_matrix_set (A, 1, 1, L_11);
  82. }
  83. for (k = 2; k < M; k++)
  84. {
  85. double A_kk = gsl_matrix_get (A, k, k);
  86. for (i = 0; i < k; i++)
  87. {
  88. double sum = 0;
  89. double A_ki = gsl_matrix_get (A, k, i);
  90. double A_ii = gsl_matrix_get (A, i, i);
  91. gsl_vector_view ci = gsl_matrix_row (A, i);
  92. gsl_vector_view ck = gsl_matrix_row (A, k);
  93. if (i > 0) {
  94. gsl_vector_view di = gsl_vector_subvector(&ci.vector, 0, i);
  95. gsl_vector_view dk = gsl_vector_subvector(&ck.vector, 0, i);
  96. gsl_blas_ddot (&di.vector, &dk.vector, &sum);
  97. }
  98. A_ki = (A_ki - sum) / A_ii;
  99. gsl_matrix_set (A, k, i, A_ki);
  100. }
  101. {
  102. gsl_vector_view ck = gsl_matrix_row (A, k);
  103. gsl_vector_view dk = gsl_vector_subvector (&ck.vector, 0, k);
  104. double sum = gsl_blas_dnrm2 (&dk.vector);
  105. double diag = A_kk - sum * sum;
  106. double L_kk = quiet_sqrt(diag);
  107. if (diag <= 0)
  108. {
  109. status = GSL_EDOM;
  110. }
  111. gsl_matrix_set (A, k, k, L_kk);
  112. }
  113. }
  114. /* Now copy the transposed lower triangle to the upper triangle,
  115. * the diagonal is common.
  116. */
  117. for (i = 1; i < M; i++)
  118. {
  119. for (j = 0; j < i; j++)
  120. {
  121. double A_ij = gsl_matrix_get (A, i, j);
  122. gsl_matrix_set (A, j, i, A_ij);
  123. }
  124. }
  125. if (status == GSL_EDOM)
  126. {
  127. GSL_ERROR ("matrix must be positive definite", GSL_EDOM);
  128. }
  129. return GSL_SUCCESS;
  130. }
  131. }
  132. int
  133. gsl_linalg_cholesky_solve (const gsl_matrix * LLT,
  134. const gsl_vector * b,
  135. gsl_vector * x)
  136. {
  137. if (LLT->size1 != LLT->size2)
  138. {
  139. GSL_ERROR ("cholesky matrix must be square", GSL_ENOTSQR);
  140. }
  141. else if (LLT->size1 != b->size)
  142. {
  143. GSL_ERROR ("matrix size must match b size", GSL_EBADLEN);
  144. }
  145. else if (LLT->size2 != x->size)
  146. {
  147. GSL_ERROR ("matrix size must match solution size", GSL_EBADLEN);
  148. }
  149. else
  150. {
  151. /* Copy x <- b */
  152. gsl_vector_memcpy (x, b);
  153. /* Solve for c using forward-substitution, L c = b */
  154. gsl_blas_dtrsv (CblasLower, CblasNoTrans, CblasNonUnit, LLT, x);
  155. /* Perform back-substitution, U x = c */
  156. gsl_blas_dtrsv (CblasUpper, CblasNoTrans, CblasNonUnit, LLT, x);
  157. return GSL_SUCCESS;
  158. }
  159. }
  160. int
  161. gsl_linalg_cholesky_svx (const gsl_matrix * LLT,
  162. gsl_vector * x)
  163. {
  164. if (LLT->size1 != LLT->size2)
  165. {
  166. GSL_ERROR ("cholesky matrix must be square", GSL_ENOTSQR);
  167. }
  168. else if (LLT->size2 != x->size)
  169. {
  170. GSL_ERROR ("matrix size must match solution size", GSL_EBADLEN);
  171. }
  172. else
  173. {
  174. /* Solve for c using forward-substitution, L c = b */
  175. gsl_blas_dtrsv (CblasLower, CblasNoTrans, CblasNonUnit, LLT, x);
  176. /* Perform back-substitution, U x = c */
  177. gsl_blas_dtrsv (CblasUpper, CblasNoTrans, CblasNonUnit, LLT, x);
  178. return GSL_SUCCESS;
  179. }
  180. }
  181. int
  182. gsl_linalg_cholesky_decomp_unit(gsl_matrix * A, gsl_vector * D)
  183. {
  184. const size_t N = A->size1;
  185. size_t i, j;
  186. /* initial Cholesky */
  187. int stat_chol = gsl_linalg_cholesky_decomp(A);
  188. if(stat_chol == GSL_SUCCESS)
  189. {
  190. /* calculate D from diagonal part of initial Cholesky */
  191. for(i = 0; i < N; ++i)
  192. {
  193. const double C_ii = gsl_matrix_get(A, i, i);
  194. gsl_vector_set(D, i, C_ii*C_ii);
  195. }
  196. /* multiply initial Cholesky by 1/sqrt(D) on the right */
  197. for(i = 0; i < N; ++i)
  198. {
  199. for(j = 0; j < N; ++j)
  200. {
  201. gsl_matrix_set(A, i, j, gsl_matrix_get(A, i, j) / sqrt(gsl_vector_get(D, j)));
  202. }
  203. }
  204. /* Because the initial Cholesky contained both L and transpose(L),
  205. the result of the multiplication is not symmetric anymore;
  206. but the lower triangle _is_ correct. Therefore we reflect
  207. it to the upper triangle and declare victory.
  208. */
  209. for(i = 0; i < N; ++i)
  210. for(j = i + 1; j < N; ++j)
  211. gsl_matrix_set(A, i, j, gsl_matrix_get(A, j, i));
  212. }
  213. return stat_chol;
  214. }