gsl_multifit__covar.c 4.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193
  1. /* multifit/covar.c
  2. *
  3. * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2004, 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 "gsl_math.h"
  21. #include "gsl_errno.h"
  22. #include "gsl_permutation.h"
  23. #include "gsl_linalg.h"
  24. #include "gsl_multifit_nlin.h"
  25. /* Compute the covariance matrix
  26. cov = inv (J^T J)
  27. by QRP^T decomposition of J
  28. */
  29. int
  30. gsl_multifit_covar (const gsl_matrix * J, double epsrel, gsl_matrix * covar)
  31. {
  32. double tolr;
  33. size_t i, j, k;
  34. size_t kmax = 0;
  35. gsl_matrix * r;
  36. gsl_vector * tau;
  37. gsl_vector * norm;
  38. gsl_permutation * perm;
  39. size_t m = J->size1, n = J->size2 ;
  40. if (m < n)
  41. {
  42. GSL_ERROR ("Jacobian be rectangular M x N with M >= N", GSL_EBADLEN);
  43. }
  44. if (covar->size1 != covar->size2 || covar->size1 != n)
  45. {
  46. GSL_ERROR ("covariance matrix must be square and match second dimension of jacobian", GSL_EBADLEN);
  47. }
  48. r = gsl_matrix_alloc (m, n);
  49. tau = gsl_vector_alloc (n);
  50. perm = gsl_permutation_alloc (n) ;
  51. norm = gsl_vector_alloc (n) ;
  52. {
  53. int signum = 0;
  54. gsl_matrix_memcpy (r, J);
  55. gsl_linalg_QRPT_decomp (r, tau, perm, &signum, norm);
  56. }
  57. /* Form the inverse of R in the full upper triangle of R */
  58. tolr = epsrel * fabs(gsl_matrix_get(r, 0, 0));
  59. for (k = 0 ; k < n ; k++)
  60. {
  61. double rkk = gsl_matrix_get(r, k, k);
  62. if (fabs(rkk) <= tolr)
  63. {
  64. break;
  65. }
  66. gsl_matrix_set(r, k, k, 1.0/rkk);
  67. for (j = 0; j < k ; j++)
  68. {
  69. double t = gsl_matrix_get(r, j, k) / rkk;
  70. gsl_matrix_set (r, j, k, 0.0);
  71. for (i = 0; i <= j; i++)
  72. {
  73. double rik = gsl_matrix_get (r, i, k);
  74. double rij = gsl_matrix_get (r, i, j);
  75. gsl_matrix_set (r, i, k, rik - t * rij);
  76. }
  77. }
  78. kmax = k;
  79. }
  80. /* Form the full upper triangle of the inverse of R^T R in the full
  81. upper triangle of R */
  82. for (k = 0; k <= kmax ; k++)
  83. {
  84. for (j = 0; j < k; j++)
  85. {
  86. double rjk = gsl_matrix_get (r, j, k);
  87. for (i = 0; i <= j ; i++)
  88. {
  89. double rij = gsl_matrix_get (r, i, j);
  90. double rik = gsl_matrix_get (r, i, k);
  91. gsl_matrix_set (r, i, j, rij + rjk * rik);
  92. }
  93. }
  94. {
  95. double t = gsl_matrix_get (r, k, k);
  96. for (i = 0; i <= k; i++)
  97. {
  98. double rik = gsl_matrix_get (r, i, k);
  99. gsl_matrix_set (r, i, k, t * rik);
  100. };
  101. }
  102. }
  103. /* Form the full lower triangle of the covariance matrix in the
  104. strict lower triangle of R and in w */
  105. for (j = 0 ; j < n ; j++)
  106. {
  107. size_t pj = gsl_permutation_get (perm, j);
  108. for (i = 0; i <= j; i++)
  109. {
  110. size_t pi = gsl_permutation_get (perm, i);
  111. double rij;
  112. if (j > kmax)
  113. {
  114. gsl_matrix_set (r, i, j, 0.0);
  115. rij = 0.0 ;
  116. }
  117. else
  118. {
  119. rij = gsl_matrix_get (r, i, j);
  120. }
  121. if (pi > pj)
  122. {
  123. gsl_matrix_set (r, pi, pj, rij);
  124. }
  125. else if (pi < pj)
  126. {
  127. gsl_matrix_set (r, pj, pi, rij);
  128. }
  129. }
  130. {
  131. double rjj = gsl_matrix_get (r, j, j);
  132. gsl_matrix_set (covar, pj, pj, rjj);
  133. }
  134. }
  135. /* symmetrize the covariance matrix */
  136. for (j = 0 ; j < n ; j++)
  137. {
  138. for (i = 0; i < j ; i++)
  139. {
  140. double rji = gsl_matrix_get (r, j, i);
  141. gsl_matrix_set (covar, j, i, rji);
  142. gsl_matrix_set (covar, i, j, rji);
  143. }
  144. }
  145. gsl_matrix_free (r);
  146. gsl_permutation_free (perm);
  147. gsl_vector_free (tau);
  148. gsl_vector_free (norm);
  149. return GSL_SUCCESS;
  150. }