gsl_linalg__tridiag.c 15 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581
  1. /* linalg/tridiag.c
  2. *
  3. * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2002, 2004, 2007 Gerard Jungman, Brian Gough, David Necas
  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. /* Author: G. Jungman */
  20. #include "gsl__config.h"
  21. #include <stdlib.h>
  22. #include <math.h>
  23. #include "gsl_errno.h"
  24. #include "gsl_linalg__tridiag.h"
  25. #include "gsl_linalg.h"
  26. /* for description of method see [Engeln-Mullges + Uhlig, p. 92]
  27. *
  28. * diag[0] offdiag[0] 0 .....
  29. * offdiag[0] diag[1] offdiag[1] .....
  30. * 0 offdiag[1] diag[2]
  31. * 0 0 offdiag[2] .....
  32. */
  33. static
  34. int
  35. solve_tridiag(
  36. const double diag[], size_t d_stride,
  37. const double offdiag[], size_t o_stride,
  38. const double b[], size_t b_stride,
  39. double x[], size_t x_stride,
  40. size_t N)
  41. {
  42. int status = GSL_SUCCESS;
  43. double *gamma = (double *) malloc (N * sizeof (double));
  44. double *alpha = (double *) malloc (N * sizeof (double));
  45. double *c = (double *) malloc (N * sizeof (double));
  46. double *z = (double *) malloc (N * sizeof (double));
  47. if (gamma == 0 || alpha == 0 || c == 0 || z == 0)
  48. {
  49. GSL_ERROR("failed to allocate working space", GSL_ENOMEM);
  50. }
  51. else
  52. {
  53. size_t i, j;
  54. /* Cholesky decomposition
  55. A = L.D.L^t
  56. lower_diag(L) = gamma
  57. diag(D) = alpha
  58. */
  59. alpha[0] = diag[0];
  60. gamma[0] = offdiag[0] / alpha[0];
  61. if (alpha[0] == 0) {
  62. status = GSL_EZERODIV;
  63. }
  64. for (i = 1; i < N - 1; i++)
  65. {
  66. alpha[i] = diag[d_stride * i] - offdiag[o_stride*(i - 1)] * gamma[i - 1];
  67. gamma[i] = offdiag[o_stride * i] / alpha[i];
  68. if (alpha[i] == 0) {
  69. status = GSL_EZERODIV;
  70. }
  71. }
  72. if (N > 1)
  73. {
  74. alpha[N - 1] = diag[d_stride * (N - 1)] - offdiag[o_stride*(N - 2)] * gamma[N - 2];
  75. }
  76. /* update RHS */
  77. z[0] = b[0];
  78. for (i = 1; i < N; i++)
  79. {
  80. z[i] = b[b_stride * i] - gamma[i - 1] * z[i - 1];
  81. }
  82. for (i = 0; i < N; i++)
  83. {
  84. c[i] = z[i] / alpha[i];
  85. }
  86. /* backsubstitution */
  87. x[x_stride * (N - 1)] = c[N - 1];
  88. if (N >= 2)
  89. {
  90. for (i = N - 2, j = 0; j <= N - 2; j++, i--)
  91. {
  92. x[x_stride * i] = c[i] - gamma[i] * x[x_stride * (i + 1)];
  93. }
  94. }
  95. }
  96. if (z != 0)
  97. free (z);
  98. if (c != 0)
  99. free (c);
  100. if (alpha != 0)
  101. free (alpha);
  102. if (gamma != 0)
  103. free (gamma);
  104. if (status == GSL_EZERODIV) {
  105. GSL_ERROR ("matrix must be positive definite", status);
  106. }
  107. return status;
  108. }
  109. /* plain gauss elimination, only not bothering with the zeroes
  110. *
  111. * diag[0] abovediag[0] 0 .....
  112. * belowdiag[0] diag[1] abovediag[1] .....
  113. * 0 belowdiag[1] diag[2]
  114. * 0 0 belowdiag[2] .....
  115. */
  116. static
  117. int
  118. solve_tridiag_nonsym(
  119. const double diag[], size_t d_stride,
  120. const double abovediag[], size_t a_stride,
  121. const double belowdiag[], size_t b_stride,
  122. const double rhs[], size_t r_stride,
  123. double x[], size_t x_stride,
  124. size_t N)
  125. {
  126. int status = GSL_SUCCESS;
  127. double *alpha = (double *) malloc (N * sizeof (double));
  128. double *z = (double *) malloc (N * sizeof (double));
  129. if (alpha == 0 || z == 0)
  130. {
  131. GSL_ERROR("failed to allocate working space", GSL_ENOMEM);
  132. }
  133. else
  134. {
  135. size_t i, j;
  136. /* Bidiagonalization (eliminating belowdiag)
  137. & rhs update
  138. diag' = alpha
  139. rhs' = z
  140. */
  141. alpha[0] = diag[0];
  142. z[0] = rhs[0];
  143. if (alpha[0] == 0) {
  144. status = GSL_EZERODIV;
  145. }
  146. for (i = 1; i < N; i++)
  147. {
  148. const double t = belowdiag[b_stride*(i - 1)]/alpha[i-1];
  149. alpha[i] = diag[d_stride*i] - t*abovediag[a_stride*(i - 1)];
  150. z[i] = rhs[r_stride*i] - t*z[i-1];
  151. if (alpha[i] == 0) {
  152. status = GSL_EZERODIV;
  153. }
  154. }
  155. /* backsubstitution */
  156. x[x_stride * (N - 1)] = z[N - 1]/alpha[N - 1];
  157. if (N >= 2)
  158. {
  159. for (i = N - 2, j = 0; j <= N - 2; j++, i--)
  160. {
  161. x[x_stride * i] = (z[i] - abovediag[a_stride*i] * x[x_stride * (i + 1)])/alpha[i];
  162. }
  163. }
  164. }
  165. if (z != 0)
  166. free (z);
  167. if (alpha != 0)
  168. free (alpha);
  169. if (status == GSL_EZERODIV) {
  170. GSL_ERROR ("matrix must be positive definite", status);
  171. }
  172. return status;
  173. }
  174. /* for description of method see [Engeln-Mullges + Uhlig, p. 96]
  175. *
  176. * diag[0] offdiag[0] 0 ..... offdiag[N-1]
  177. * offdiag[0] diag[1] offdiag[1] .....
  178. * 0 offdiag[1] diag[2]
  179. * 0 0 offdiag[2] .....
  180. * ... ...
  181. * offdiag[N-1] ...
  182. *
  183. */
  184. static
  185. int
  186. solve_cyc_tridiag(
  187. const double diag[], size_t d_stride,
  188. const double offdiag[], size_t o_stride,
  189. const double b[], size_t b_stride,
  190. double x[], size_t x_stride,
  191. size_t N)
  192. {
  193. int status = GSL_SUCCESS;
  194. double * delta = (double *) malloc (N * sizeof (double));
  195. double * gamma = (double *) malloc (N * sizeof (double));
  196. double * alpha = (double *) malloc (N * sizeof (double));
  197. double * c = (double *) malloc (N * sizeof (double));
  198. double * z = (double *) malloc (N * sizeof (double));
  199. if (delta == 0 || gamma == 0 || alpha == 0 || c == 0 || z == 0)
  200. {
  201. GSL_ERROR("failed to allocate working space", GSL_ENOMEM);
  202. }
  203. else
  204. {
  205. size_t i, j;
  206. double sum = 0.0;
  207. /* factor */
  208. if (N == 1)
  209. {
  210. x[0] = b[0] / diag[0];
  211. return GSL_SUCCESS;
  212. }
  213. alpha[0] = diag[0];
  214. gamma[0] = offdiag[0] / alpha[0];
  215. delta[0] = offdiag[o_stride * (N-1)] / alpha[0];
  216. if (alpha[0] == 0) {
  217. status = GSL_EZERODIV;
  218. }
  219. for (i = 1; i < N - 2; i++)
  220. {
  221. alpha[i] = diag[d_stride * i] - offdiag[o_stride * (i-1)] * gamma[i - 1];
  222. gamma[i] = offdiag[o_stride * i] / alpha[i];
  223. delta[i] = -delta[i - 1] * offdiag[o_stride * (i-1)] / alpha[i];
  224. if (alpha[i] == 0) {
  225. status = GSL_EZERODIV;
  226. }
  227. }
  228. for (i = 0; i < N - 2; i++)
  229. {
  230. sum += alpha[i] * delta[i] * delta[i];
  231. }
  232. alpha[N - 2] = diag[d_stride * (N - 2)] - offdiag[o_stride * (N - 3)] * gamma[N - 3];
  233. gamma[N - 2] = (offdiag[o_stride * (N - 2)] - offdiag[o_stride * (N - 3)] * delta[N - 3]) / alpha[N - 2];
  234. alpha[N - 1] = diag[d_stride * (N - 1)] - sum - alpha[(N - 2)] * gamma[N - 2] * gamma[N - 2];
  235. /* update */
  236. z[0] = b[0];
  237. for (i = 1; i < N - 1; i++)
  238. {
  239. z[i] = b[b_stride * i] - z[i - 1] * gamma[i - 1];
  240. }
  241. sum = 0.0;
  242. for (i = 0; i < N - 2; i++)
  243. {
  244. sum += delta[i] * z[i];
  245. }
  246. z[N - 1] = b[b_stride * (N - 1)] - sum - gamma[N - 2] * z[N - 2];
  247. for (i = 0; i < N; i++)
  248. {
  249. c[i] = z[i] / alpha[i];
  250. }
  251. /* backsubstitution */
  252. x[x_stride * (N - 1)] = c[N - 1];
  253. x[x_stride * (N - 2)] = c[N - 2] - gamma[N - 2] * x[x_stride * (N - 1)];
  254. if (N >= 3)
  255. {
  256. for (i = N - 3, j = 0; j <= N - 3; j++, i--)
  257. {
  258. x[x_stride * i] = c[i] - gamma[i] * x[x_stride * (i + 1)] - delta[i] * x[x_stride * (N - 1)];
  259. }
  260. }
  261. }
  262. if (z != 0)
  263. free (z);
  264. if (c != 0)
  265. free (c);
  266. if (alpha != 0)
  267. free (alpha);
  268. if (gamma != 0)
  269. free (gamma);
  270. if (delta != 0)
  271. free (delta);
  272. if (status == GSL_EZERODIV) {
  273. GSL_ERROR ("matrix must be positive definite", status);
  274. }
  275. return status;
  276. }
  277. /* solve following system w/o the corner elements and then use
  278. * Sherman-Morrison formula to compensate for them
  279. *
  280. * diag[0] abovediag[0] 0 ..... belowdiag[N-1]
  281. * belowdiag[0] diag[1] abovediag[1] .....
  282. * 0 belowdiag[1] diag[2]
  283. * 0 0 belowdiag[2] .....
  284. * ... ...
  285. * abovediag[N-1] ...
  286. */
  287. static
  288. int solve_cyc_tridiag_nonsym(
  289. const double diag[], size_t d_stride,
  290. const double abovediag[], size_t a_stride,
  291. const double belowdiag[], size_t b_stride,
  292. const double rhs[], size_t r_stride,
  293. double x[], size_t x_stride,
  294. size_t N)
  295. {
  296. int status = GSL_SUCCESS;
  297. double *alpha = (double *) malloc (N * sizeof (double));
  298. double *zb = (double *) malloc (N * sizeof (double));
  299. double *zu = (double *) malloc (N * sizeof (double));
  300. double *w = (double *) malloc (N * sizeof (double));
  301. if (alpha == 0 || zb == 0 || zu == 0 || w == 0)
  302. {
  303. GSL_ERROR("failed to allocate working space", GSL_ENOMEM);
  304. }
  305. else
  306. {
  307. double beta;
  308. /* Bidiagonalization (eliminating belowdiag)
  309. & rhs update
  310. diag' = alpha
  311. rhs' = zb
  312. rhs' for Aq=u is zu
  313. */
  314. zb[0] = rhs[0];
  315. if (diag[0] != 0) beta = -diag[0]; else beta = 1;
  316. {
  317. const double q = 1 - abovediag[0]*belowdiag[0]/(diag[0]*diag[d_stride]);
  318. if (fabs(q/beta) > 0.5 && fabs(q/beta) < 2) {
  319. beta *= (fabs(q/beta) < 1) ? 0.5 : 2;
  320. }
  321. }
  322. zu[0] = beta;
  323. alpha[0] = diag[0] - beta;
  324. if (alpha[0] == 0) {
  325. status = GSL_EZERODIV;
  326. }
  327. {
  328. size_t i;
  329. for (i = 1; i+1 < N; i++)
  330. {
  331. const double t = belowdiag[b_stride*(i - 1)]/alpha[i-1];
  332. alpha[i] = diag[d_stride*i] - t*abovediag[a_stride*(i - 1)];
  333. zb[i] = rhs[r_stride*i] - t*zb[i-1];
  334. zu[i] = -t*zu[i-1];
  335. /* FIXME!!! */
  336. if (alpha[i] == 0) {
  337. status = GSL_EZERODIV;
  338. }
  339. }
  340. }
  341. {
  342. const size_t i = N-1;
  343. const double t = belowdiag[b_stride*(i - 1)]/alpha[i-1];
  344. alpha[i] = diag[d_stride*i]
  345. - abovediag[a_stride*i]*belowdiag[b_stride*i]/beta
  346. - t*abovediag[a_stride*(i - 1)];
  347. zb[i] = rhs[r_stride*i] - t*zb[i-1];
  348. zu[i] = abovediag[a_stride*i] - t*zu[i-1];
  349. /* FIXME!!! */
  350. if (alpha[i] == 0) {
  351. status = GSL_EZERODIV;
  352. }
  353. }
  354. /* backsubstitution */
  355. {
  356. size_t i, j;
  357. w[N-1] = zu[N-1]/alpha[N-1];
  358. x[N-1] = zb[N-1]/alpha[N-1];
  359. for (i = N - 2, j = 0; j <= N - 2; j++, i--)
  360. {
  361. w[i] = (zu[i] - abovediag[a_stride*i] * w[i+1])/alpha[i];
  362. x[i*x_stride] = (zb[i] - abovediag[a_stride*i] * x[x_stride*(i + 1)])/alpha[i];
  363. }
  364. }
  365. /* Sherman-Morrison */
  366. {
  367. const double vw = w[0] + belowdiag[b_stride*(N - 1)]/beta * w[N-1];
  368. const double vx = x[0] + belowdiag[b_stride*(N - 1)]/beta * x[x_stride*(N - 1)];
  369. /* FIXME!!! */
  370. if (vw + 1 == 0) {
  371. status = GSL_EZERODIV;
  372. }
  373. {
  374. size_t i;
  375. for (i = 0; i < N; i++)
  376. x[i] -= vx/(1 + vw)*w[i];
  377. }
  378. }
  379. }
  380. if (zb != 0)
  381. free (zb);
  382. if (zu != 0)
  383. free (zu);
  384. if (w != 0)
  385. free (w);
  386. if (alpha != 0)
  387. free (alpha);
  388. if (status == GSL_EZERODIV) {
  389. GSL_ERROR ("matrix must be positive definite", status);
  390. }
  391. return status;
  392. }
  393. int
  394. gsl_linalg_solve_symm_tridiag(
  395. const gsl_vector * diag,
  396. const gsl_vector * offdiag,
  397. const gsl_vector * rhs,
  398. gsl_vector * solution)
  399. {
  400. if(diag->size != rhs->size)
  401. {
  402. GSL_ERROR ("size of diag must match rhs", GSL_EBADLEN);
  403. }
  404. else if (offdiag->size != rhs->size-1)
  405. {
  406. GSL_ERROR ("size of offdiag must match rhs-1", GSL_EBADLEN);
  407. }
  408. else if (solution->size != rhs->size)
  409. {
  410. GSL_ERROR ("size of solution must match rhs", GSL_EBADLEN);
  411. }
  412. else
  413. {
  414. return solve_tridiag(diag->data, diag->stride,
  415. offdiag->data, offdiag->stride,
  416. rhs->data, rhs->stride,
  417. solution->data, solution->stride,
  418. diag->size);
  419. }
  420. }
  421. int
  422. gsl_linalg_solve_tridiag(
  423. const gsl_vector * diag,
  424. const gsl_vector * abovediag,
  425. const gsl_vector * belowdiag,
  426. const gsl_vector * rhs,
  427. gsl_vector * solution)
  428. {
  429. if(diag->size != rhs->size)
  430. {
  431. GSL_ERROR ("size of diag must match rhs", GSL_EBADLEN);
  432. }
  433. else if (abovediag->size != rhs->size-1)
  434. {
  435. GSL_ERROR ("size of abovediag must match rhs-1", GSL_EBADLEN);
  436. }
  437. else if (belowdiag->size != rhs->size-1)
  438. {
  439. GSL_ERROR ("size of belowdiag must match rhs-1", GSL_EBADLEN);
  440. }
  441. else if (solution->size != rhs->size)
  442. {
  443. GSL_ERROR ("size of solution must match rhs", GSL_EBADLEN);
  444. }
  445. else
  446. {
  447. return solve_tridiag_nonsym(diag->data, diag->stride,
  448. abovediag->data, abovediag->stride,
  449. belowdiag->data, belowdiag->stride,
  450. rhs->data, rhs->stride,
  451. solution->data, solution->stride,
  452. diag->size);
  453. }
  454. }
  455. int
  456. gsl_linalg_solve_symm_cyc_tridiag(
  457. const gsl_vector * diag,
  458. const gsl_vector * offdiag,
  459. const gsl_vector * rhs,
  460. gsl_vector * solution)
  461. {
  462. if(diag->size != rhs->size)
  463. {
  464. GSL_ERROR ("size of diag must match rhs", GSL_EBADLEN);
  465. }
  466. else if (offdiag->size != rhs->size)
  467. {
  468. GSL_ERROR ("size of offdiag must match rhs", GSL_EBADLEN);
  469. }
  470. else if (solution->size != rhs->size)
  471. {
  472. GSL_ERROR ("size of solution must match rhs", GSL_EBADLEN);
  473. }
  474. else if (diag->size < 3)
  475. {
  476. GSL_ERROR ("size of cyclic system must be 3 or more", GSL_EBADLEN);
  477. }
  478. else
  479. {
  480. return solve_cyc_tridiag(diag->data, diag->stride,
  481. offdiag->data, offdiag->stride,
  482. rhs->data, rhs->stride,
  483. solution->data, solution->stride,
  484. diag->size);
  485. }
  486. }
  487. int
  488. gsl_linalg_solve_cyc_tridiag(
  489. const gsl_vector * diag,
  490. const gsl_vector * abovediag,
  491. const gsl_vector * belowdiag,
  492. const gsl_vector * rhs,
  493. gsl_vector * solution)
  494. {
  495. if(diag->size != rhs->size)
  496. {
  497. GSL_ERROR ("size of diag must match rhs", GSL_EBADLEN);
  498. }
  499. else if (abovediag->size != rhs->size)
  500. {
  501. GSL_ERROR ("size of abovediag must match rhs", GSL_EBADLEN);
  502. }
  503. else if (belowdiag->size != rhs->size)
  504. {
  505. GSL_ERROR ("size of belowdiag must match rhs", GSL_EBADLEN);
  506. }
  507. else if (solution->size != rhs->size)
  508. {
  509. GSL_ERROR ("size of solution must match rhs", GSL_EBADLEN);
  510. }
  511. else if (diag->size < 3)
  512. {
  513. GSL_ERROR ("size of cyclic system must be 3 or more", GSL_EBADLEN);
  514. }
  515. else
  516. {
  517. return solve_cyc_tridiag_nonsym(diag->data, diag->stride,
  518. abovediag->data, abovediag->stride,
  519. belowdiag->data, belowdiag->stride,
  520. rhs->data, rhs->stride,
  521. solution->data, solution->stride,
  522. diag->size);
  523. }
  524. }