gsl_linalg__choleskyc.c 5.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204
  1. /* linalg/choleskyc.c
  2. *
  3. * Copyright (C) 2007 Patrick Alken
  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_matrix.h"
  21. #include "gsl_vector.h"
  22. #include "gsl_linalg.h"
  23. #include "gsl_math.h"
  24. #include "gsl_complex.h"
  25. #include "gsl_complex_math.h"
  26. #include "gsl_blas.h"
  27. #include "gsl_errno.h"
  28. /*
  29. * This module contains routines related to the Cholesky decomposition
  30. * of a complex Hermitian positive definite matrix.
  31. */
  32. static void cholesky_complex_conj_vector(gsl_vector_complex *v);
  33. /*
  34. gsl_linalg_complex_cholesky_decomp()
  35. Perform the Cholesky decomposition on a Hermitian positive definite
  36. matrix. See Golub & Van Loan, "Matrix Computations" (3rd ed),
  37. algorithm 4.2.2.
  38. Inputs: A - (input/output) complex postive definite matrix
  39. Return: success or error
  40. The lower triangle of A is overwritten with the Cholesky decomposition
  41. */
  42. int
  43. gsl_linalg_complex_cholesky_decomp(gsl_matrix_complex *A)
  44. {
  45. const size_t N = A->size1;
  46. if (N != A->size2)
  47. {
  48. GSL_ERROR("cholesky decomposition requires square matrix", GSL_ENOTSQR);
  49. }
  50. else
  51. {
  52. size_t j;
  53. gsl_complex z;
  54. double ajj;
  55. for (j = 0; j < N; ++j)
  56. {
  57. z = gsl_matrix_complex_get(A, j, j);
  58. ajj = GSL_REAL(z);
  59. if (j > 0)
  60. {
  61. gsl_vector_complex_const_view aj =
  62. gsl_matrix_complex_const_subrow(A, j, 0, j);
  63. gsl_blas_zdotc(&aj.vector, &aj.vector, &z);
  64. ajj -= GSL_REAL(z);
  65. }
  66. if (ajj <= 0.0)
  67. {
  68. GSL_ERROR("matrix is not positive definite", GSL_EDOM);
  69. }
  70. ajj = sqrt(ajj);
  71. GSL_SET_COMPLEX(&z, ajj, 0.0);
  72. gsl_matrix_complex_set(A, j, j, z);
  73. if (j < N - 1)
  74. {
  75. gsl_vector_complex_view av =
  76. gsl_matrix_complex_subcolumn(A, j, j + 1, N - j - 1);
  77. if (j > 0)
  78. {
  79. gsl_vector_complex_view aj =
  80. gsl_matrix_complex_subrow(A, j, 0, j);
  81. gsl_matrix_complex_view am =
  82. gsl_matrix_complex_submatrix(A, j + 1, 0, N - j - 1, j);
  83. cholesky_complex_conj_vector(&aj.vector);
  84. gsl_blas_zgemv(CblasNoTrans,
  85. GSL_COMPLEX_NEGONE,
  86. &am.matrix,
  87. &aj.vector,
  88. GSL_COMPLEX_ONE,
  89. &av.vector);
  90. cholesky_complex_conj_vector(&aj.vector);
  91. }
  92. gsl_blas_zdscal(1.0 / ajj, &av.vector);
  93. }
  94. }
  95. return GSL_SUCCESS;
  96. }
  97. } /* gsl_linalg_complex_cholesky_decomp() */
  98. /*
  99. gsl_linalg_complex_cholesky_solve()
  100. Solve A x = b where A is in cholesky form
  101. */
  102. int
  103. gsl_linalg_complex_cholesky_solve (const gsl_matrix_complex * cholesky,
  104. const gsl_vector_complex * b,
  105. gsl_vector_complex * x)
  106. {
  107. if (cholesky->size1 != cholesky->size2)
  108. {
  109. GSL_ERROR ("cholesky matrix must be square", GSL_ENOTSQR);
  110. }
  111. else if (cholesky->size1 != b->size)
  112. {
  113. GSL_ERROR ("matrix size must match b size", GSL_EBADLEN);
  114. }
  115. else if (cholesky->size2 != x->size)
  116. {
  117. GSL_ERROR ("matrix size must match solution size", GSL_EBADLEN);
  118. }
  119. else
  120. {
  121. gsl_vector_complex_memcpy (x, b);
  122. /* solve for y using forward-substitution, L y = b */
  123. gsl_blas_ztrsv (CblasLower, CblasNoTrans, CblasNonUnit, cholesky, x);
  124. /* perform back-substitution, L^H x = y */
  125. gsl_blas_ztrsv (CblasLower, CblasConjTrans, CblasNonUnit, cholesky, x);
  126. return GSL_SUCCESS;
  127. }
  128. } /* gsl_linalg_complex_cholesky_solve() */
  129. /*
  130. gsl_linalg_complex_cholesky_svx()
  131. Solve A x = b in place where A is in cholesky form
  132. */
  133. int
  134. gsl_linalg_complex_cholesky_svx (const gsl_matrix_complex * cholesky,
  135. gsl_vector_complex * x)
  136. {
  137. if (cholesky->size1 != cholesky->size2)
  138. {
  139. GSL_ERROR ("cholesky matrix must be square", GSL_ENOTSQR);
  140. }
  141. else if (cholesky->size2 != x->size)
  142. {
  143. GSL_ERROR ("matrix size must match solution size", GSL_EBADLEN);
  144. }
  145. else
  146. {
  147. /* solve for y using forward-substitution, L y = b */
  148. gsl_blas_ztrsv (CblasLower, CblasNoTrans, CblasNonUnit, cholesky, x);
  149. /* perform back-substitution, L^H x = y */
  150. gsl_blas_ztrsv (CblasLower, CblasConjTrans, CblasNonUnit, cholesky, x);
  151. return GSL_SUCCESS;
  152. }
  153. } /* gsl_linalg_complex_cholesky_svx() */
  154. /********************************************
  155. * INTERNAL ROUTINES *
  156. ********************************************/
  157. static void
  158. cholesky_complex_conj_vector(gsl_vector_complex *v)
  159. {
  160. size_t i;
  161. for (i = 0; i < v->size; ++i)
  162. {
  163. gsl_complex z = gsl_vector_complex_get(v, i);
  164. gsl_vector_complex_set(v, i, gsl_complex_conjugate(z));
  165. }
  166. } /* cholesky_complex_conj_vector() */