gsl_linalg__luc.c 8.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335
  1. /* linalg/luc.c
  2. *
  3. * Copyright (C) 2001, 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 <stdlib.h>
  21. #include <string.h>
  22. #include "gsl_math.h"
  23. #include "gsl_vector.h"
  24. #include "gsl_matrix.h"
  25. #include "gsl_complex.h"
  26. #include "gsl_complex_math.h"
  27. #include "gsl_permute_vector.h"
  28. #include "gsl_blas.h"
  29. #include "gsl_complex_math.h"
  30. #include "gsl_linalg.h"
  31. /* Factorise a general N x N complex matrix A into,
  32. *
  33. * P A = L U
  34. *
  35. * where P is a permutation matrix, L is unit lower triangular and U
  36. * is upper triangular.
  37. *
  38. * L is stored in the strict lower triangular part of the input
  39. * matrix. The diagonal elements of L are unity and are not stored.
  40. *
  41. * U is stored in the diagonal and upper triangular part of the
  42. * input matrix.
  43. *
  44. * P is stored in the permutation p. Column j of P is column k of the
  45. * identity matrix, where k = permutation->data[j]
  46. *
  47. * signum gives the sign of the permutation, (-1)^n, where n is the
  48. * number of interchanges in the permutation.
  49. *
  50. * See Golub & Van Loan, Matrix Computations, Algorithm 3.4.1 (Gauss
  51. * Elimination with Partial Pivoting).
  52. */
  53. int
  54. gsl_linalg_complex_LU_decomp (gsl_matrix_complex * A, gsl_permutation * p, int *signum)
  55. {
  56. if (A->size1 != A->size2)
  57. {
  58. GSL_ERROR ("LU decomposition requires square matrix", GSL_ENOTSQR);
  59. }
  60. else if (p->size != A->size1)
  61. {
  62. GSL_ERROR ("permutation length must match matrix size", GSL_EBADLEN);
  63. }
  64. else
  65. {
  66. const size_t N = A->size1;
  67. size_t i, j, k;
  68. *signum = 1;
  69. gsl_permutation_init (p);
  70. for (j = 0; j < N - 1; j++)
  71. {
  72. /* Find maximum in the j-th column */
  73. gsl_complex ajj = gsl_matrix_complex_get (A, j, j);
  74. double max = gsl_complex_abs (ajj);
  75. size_t i_pivot = j;
  76. for (i = j + 1; i < N; i++)
  77. {
  78. gsl_complex aij = gsl_matrix_complex_get (A, i, j);
  79. double ai = gsl_complex_abs (aij);
  80. if (ai > max)
  81. {
  82. max = ai;
  83. i_pivot = i;
  84. }
  85. }
  86. if (i_pivot != j)
  87. {
  88. gsl_matrix_complex_swap_rows (A, j, i_pivot);
  89. gsl_permutation_swap (p, j, i_pivot);
  90. *signum = -(*signum);
  91. }
  92. ajj = gsl_matrix_complex_get (A, j, j);
  93. if (!(GSL_REAL(ajj) == 0.0 && GSL_IMAG(ajj) == 0.0))
  94. {
  95. for (i = j + 1; i < N; i++)
  96. {
  97. gsl_complex aij_orig = gsl_matrix_complex_get (A, i, j);
  98. gsl_complex aij = gsl_complex_div (aij_orig, ajj);
  99. gsl_matrix_complex_set (A, i, j, aij);
  100. for (k = j + 1; k < N; k++)
  101. {
  102. gsl_complex aik = gsl_matrix_complex_get (A, i, k);
  103. gsl_complex ajk = gsl_matrix_complex_get (A, j, k);
  104. /* aik = aik - aij * ajk */
  105. gsl_complex aijajk = gsl_complex_mul (aij, ajk);
  106. gsl_complex aik_new = gsl_complex_sub (aik, aijajk);
  107. gsl_matrix_complex_set (A, i, k, aik_new);
  108. }
  109. }
  110. }
  111. }
  112. return GSL_SUCCESS;
  113. }
  114. }
  115. int
  116. gsl_linalg_complex_LU_solve (const gsl_matrix_complex * LU, const gsl_permutation * p, const gsl_vector_complex * b, gsl_vector_complex * x)
  117. {
  118. if (LU->size1 != LU->size2)
  119. {
  120. GSL_ERROR ("LU matrix must be square", GSL_ENOTSQR);
  121. }
  122. else if (LU->size1 != p->size)
  123. {
  124. GSL_ERROR ("permutation length must match matrix size", GSL_EBADLEN);
  125. }
  126. else if (LU->size1 != b->size)
  127. {
  128. GSL_ERROR ("matrix size must match b size", GSL_EBADLEN);
  129. }
  130. else if (LU->size2 != x->size)
  131. {
  132. GSL_ERROR ("matrix size must match solution size", GSL_EBADLEN);
  133. }
  134. else
  135. {
  136. /* Copy x <- b */
  137. gsl_vector_complex_memcpy (x, b);
  138. /* Solve for x */
  139. gsl_linalg_complex_LU_svx (LU, p, x);
  140. return GSL_SUCCESS;
  141. }
  142. }
  143. int
  144. gsl_linalg_complex_LU_svx (const gsl_matrix_complex * LU, const gsl_permutation * p, gsl_vector_complex * x)
  145. {
  146. if (LU->size1 != LU->size2)
  147. {
  148. GSL_ERROR ("LU matrix must be square", GSL_ENOTSQR);
  149. }
  150. else if (LU->size1 != p->size)
  151. {
  152. GSL_ERROR ("permutation length must match matrix size", GSL_EBADLEN);
  153. }
  154. else if (LU->size1 != x->size)
  155. {
  156. GSL_ERROR ("matrix size must match solution/rhs size", GSL_EBADLEN);
  157. }
  158. else
  159. {
  160. /* Apply permutation to RHS */
  161. gsl_permute_vector_complex (p, x);
  162. /* Solve for c using forward-substitution, L c = P b */
  163. gsl_blas_ztrsv (CblasLower, CblasNoTrans, CblasUnit, LU, x);
  164. /* Perform back-substitution, U x = c */
  165. gsl_blas_ztrsv (CblasUpper, CblasNoTrans, CblasNonUnit, LU, x);
  166. return GSL_SUCCESS;
  167. }
  168. }
  169. int
  170. gsl_linalg_complex_LU_refine (const gsl_matrix_complex * A, const gsl_matrix_complex * LU, const gsl_permutation * p, const gsl_vector_complex * b, gsl_vector_complex * x, gsl_vector_complex * residual)
  171. {
  172. if (A->size1 != A->size2)
  173. {
  174. GSL_ERROR ("matrix a must be square", GSL_ENOTSQR);
  175. }
  176. if (LU->size1 != LU->size2)
  177. {
  178. GSL_ERROR ("LU matrix must be square", GSL_ENOTSQR);
  179. }
  180. else if (A->size1 != LU->size2)
  181. {
  182. GSL_ERROR ("LU matrix must be decomposition of a", GSL_ENOTSQR);
  183. }
  184. else if (LU->size1 != p->size)
  185. {
  186. GSL_ERROR ("permutation length must match matrix size", GSL_EBADLEN);
  187. }
  188. else if (LU->size1 != b->size)
  189. {
  190. GSL_ERROR ("matrix size must match b size", GSL_EBADLEN);
  191. }
  192. else if (LU->size1 != x->size)
  193. {
  194. GSL_ERROR ("matrix size must match solution size", GSL_EBADLEN);
  195. }
  196. else
  197. {
  198. /* Compute residual, residual = (A * x - b) */
  199. gsl_vector_complex_memcpy (residual, b);
  200. {
  201. gsl_complex one = GSL_COMPLEX_ONE;
  202. gsl_complex negone = GSL_COMPLEX_NEGONE;
  203. gsl_blas_zgemv (CblasNoTrans, one, A, x, negone, residual);
  204. }
  205. /* Find correction, delta = - (A^-1) * residual, and apply it */
  206. gsl_linalg_complex_LU_svx (LU, p, residual);
  207. {
  208. gsl_complex negone= GSL_COMPLEX_NEGONE;
  209. gsl_blas_zaxpy (negone, residual, x);
  210. }
  211. return GSL_SUCCESS;
  212. }
  213. }
  214. int
  215. gsl_linalg_complex_LU_invert (const gsl_matrix_complex * LU, const gsl_permutation * p, gsl_matrix_complex * inverse)
  216. {
  217. size_t i, n = LU->size1;
  218. int status = GSL_SUCCESS;
  219. gsl_matrix_complex_set_identity (inverse);
  220. for (i = 0; i < n; i++)
  221. {
  222. gsl_vector_complex_view c = gsl_matrix_complex_column (inverse, i);
  223. int status_i = gsl_linalg_complex_LU_svx (LU, p, &(c.vector));
  224. if (status_i)
  225. status = status_i;
  226. }
  227. return status;
  228. }
  229. gsl_complex
  230. gsl_linalg_complex_LU_det (gsl_matrix_complex * LU, int signum)
  231. {
  232. size_t i, n = LU->size1;
  233. gsl_complex det = gsl_complex_rect((double) signum, 0.0);
  234. for (i = 0; i < n; i++)
  235. {
  236. gsl_complex zi = gsl_matrix_complex_get (LU, i, i);
  237. det = gsl_complex_mul (det, zi);
  238. }
  239. return det;
  240. }
  241. double
  242. gsl_linalg_complex_LU_lndet (gsl_matrix_complex * LU)
  243. {
  244. size_t i, n = LU->size1;
  245. double lndet = 0.0;
  246. for (i = 0; i < n; i++)
  247. {
  248. gsl_complex z = gsl_matrix_complex_get (LU, i, i);
  249. lndet += log (gsl_complex_abs (z));
  250. }
  251. return lndet;
  252. }
  253. gsl_complex
  254. gsl_linalg_complex_LU_sgndet (gsl_matrix_complex * LU, int signum)
  255. {
  256. size_t i, n = LU->size1;
  257. gsl_complex phase = gsl_complex_rect((double) signum, 0.0);
  258. for (i = 0; i < n; i++)
  259. {
  260. gsl_complex z = gsl_matrix_complex_get (LU, i, i);
  261. double r = gsl_complex_abs(z);
  262. if (r == 0)
  263. {
  264. phase = gsl_complex_rect(0.0, 0.0);
  265. break;
  266. }
  267. else
  268. {
  269. z = gsl_complex_div_real(z, r);
  270. phase = gsl_complex_mul(phase, z);
  271. }
  272. }
  273. return phase;
  274. }