gsl_multifit__lmiterate.c 5.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230
  1. static int
  2. iterate (void *vstate, gsl_multifit_function_fdf * fdf, gsl_vector * x, gsl_vector * f, gsl_matrix * J, gsl_vector * dx, int scale)
  3. {
  4. lmder_state_t *state = (lmder_state_t *) vstate;
  5. gsl_matrix *r = state->r;
  6. gsl_vector *tau = state->tau;
  7. gsl_vector *diag = state->diag;
  8. gsl_vector *qtf = state->qtf;
  9. gsl_vector *x_trial = state->x_trial;
  10. gsl_vector *f_trial = state->f_trial;
  11. gsl_vector *rptdx = state->rptdx;
  12. gsl_vector *newton = state->newton;
  13. gsl_vector *gradient = state->gradient;
  14. gsl_vector *sdiag = state->sdiag;
  15. gsl_vector *w = state->w;
  16. gsl_vector *work1 = state->work1;
  17. gsl_permutation *perm = state->perm;
  18. double prered, actred;
  19. double pnorm, fnorm1, fnorm1p, gnorm;
  20. double ratio;
  21. double dirder;
  22. int iter = 0;
  23. double p1 = 0.1, p25 = 0.25, p5 = 0.5, p75 = 0.75, p0001 = 0.0001;
  24. if (state->fnorm == 0.0)
  25. {
  26. return GSL_SUCCESS;
  27. }
  28. /* Compute qtf = Q^T f */
  29. gsl_vector_memcpy (qtf, f);
  30. gsl_linalg_QR_QTvec (r, tau, qtf);
  31. /* Compute norm of scaled gradient */
  32. compute_gradient_direction (r, perm, qtf, diag, gradient);
  33. {
  34. size_t iamax = gsl_blas_idamax (gradient);
  35. gnorm = fabs(gsl_vector_get (gradient, iamax) / state->fnorm);
  36. }
  37. /* Determine the Levenberg-Marquardt parameter */
  38. lm_iteration:
  39. iter++ ;
  40. {
  41. int status = lmpar (r, perm, qtf, diag, state->delta, &(state->par), newton, gradient, sdiag, dx, w);
  42. if (status)
  43. return status;
  44. }
  45. /* Take a trial step */
  46. gsl_vector_scale (dx, -1.0); /* reverse the step to go downhill */
  47. compute_trial_step (x, dx, state->x_trial);
  48. pnorm = scaled_enorm (diag, dx);
  49. if (state->iter == 1)
  50. {
  51. if (pnorm < state->delta)
  52. {
  53. #ifdef DEBUG
  54. printf("set delta = pnorm = %g\n" , pnorm);
  55. #endif
  56. state->delta = pnorm;
  57. }
  58. }
  59. /* Evaluate function at x + p */
  60. /* return immediately if evaluation raised error */
  61. {
  62. int status = GSL_MULTIFIT_FN_EVAL_F (fdf, x_trial, f_trial);
  63. if (status)
  64. return status;
  65. }
  66. fnorm1 = enorm (f_trial);
  67. /* Compute the scaled actual reduction */
  68. actred = compute_actual_reduction (state->fnorm, fnorm1);
  69. #ifdef DEBUG
  70. printf("lmiterate: fnorm = %g fnorm1 = %g actred = %g\n", state->fnorm, fnorm1, actred);
  71. printf("r = "); gsl_matrix_fprintf(stdout, r, "%g");
  72. printf("perm = "); gsl_permutation_fprintf(stdout, perm, "%d");
  73. printf("dx = "); gsl_vector_fprintf(stdout, dx, "%g");
  74. #endif
  75. /* Compute rptdx = R P^T dx, noting that |J dx| = |R P^T dx| */
  76. compute_rptdx (r, perm, dx, rptdx);
  77. #ifdef DEBUG
  78. printf("rptdx = "); gsl_vector_fprintf(stdout, rptdx, "%g");
  79. #endif
  80. fnorm1p = enorm (rptdx);
  81. /* Compute the scaled predicted reduction = |J dx|^2 + 2 par |D dx|^2 */
  82. {
  83. double t1 = fnorm1p / state->fnorm;
  84. double t2 = (sqrt(state->par) * pnorm) / state->fnorm;
  85. prered = t1 * t1 + t2 * t2 / p5;
  86. dirder = -(t1 * t1 + t2 * t2);
  87. }
  88. /* compute the ratio of the actual to predicted reduction */
  89. if (prered > 0)
  90. {
  91. ratio = actred / prered;
  92. }
  93. else
  94. {
  95. ratio = 0;
  96. }
  97. #ifdef DEBUG
  98. printf("lmiterate: prered = %g dirder = %g ratio = %g\n", prered, dirder,ratio);
  99. #endif
  100. /* update the step bound */
  101. if (ratio > p25)
  102. {
  103. #ifdef DEBUG
  104. printf("ratio > p25\n");
  105. #endif
  106. if (state->par == 0 || ratio >= p75)
  107. {
  108. state->delta = pnorm / p5;
  109. state->par *= p5;
  110. #ifdef DEBUG
  111. printf("updated step bounds: delta = %g, par = %g\n", state->delta, state->par);
  112. #endif
  113. }
  114. }
  115. else
  116. {
  117. double temp = (actred >= 0) ? p5 : p5*dirder / (dirder + p5 * actred);
  118. #ifdef DEBUG
  119. printf("ratio < p25\n");
  120. #endif
  121. if (p1 * fnorm1 >= state->fnorm || temp < p1 )
  122. {
  123. temp = p1;
  124. }
  125. state->delta = temp * GSL_MIN_DBL (state->delta, pnorm/p1);
  126. state->par /= temp;
  127. #ifdef DEBUG
  128. printf("updated step bounds: delta = %g, par = %g\n", state->delta, state->par);
  129. #endif
  130. }
  131. /* test for successful iteration, termination and stringent tolerances */
  132. if (ratio >= p0001)
  133. {
  134. gsl_vector_memcpy (x, x_trial);
  135. gsl_vector_memcpy (f, f_trial);
  136. /* return immediately if evaluation raised error */
  137. {
  138. int status = GSL_MULTIFIT_FN_EVAL_DF (fdf, x_trial, J);
  139. if (status)
  140. return status;
  141. }
  142. /* wa2_j = diag_j * x_j */
  143. state->xnorm = scaled_enorm(diag, x);
  144. state->fnorm = fnorm1;
  145. state->iter++;
  146. /* Rescale if necessary */
  147. if (scale)
  148. {
  149. update_diag (J, diag);
  150. }
  151. {
  152. int signum;
  153. gsl_matrix_memcpy (r, J);
  154. gsl_linalg_QRPT_decomp (r, tau, perm, &signum, work1);
  155. }
  156. return GSL_SUCCESS;
  157. }
  158. else if (fabs(actred) <= GSL_DBL_EPSILON && prered <= GSL_DBL_EPSILON
  159. && p5 * ratio <= 1.0)
  160. {
  161. return GSL_ETOLF ;
  162. }
  163. else if (state->delta <= GSL_DBL_EPSILON * state->xnorm)
  164. {
  165. return GSL_ETOLX;
  166. }
  167. else if (gnorm <= GSL_DBL_EPSILON)
  168. {
  169. return GSL_ETOLG;
  170. }
  171. else if (iter < 10)
  172. {
  173. /* Repeat inner loop if unsuccessful */
  174. goto lm_iteration;
  175. }
  176. return GSL_CONTINUE;
  177. }