3 * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2007 Brian Gough
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.
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.
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.
22 static void compute_diag (const gsl_matrix * J, gsl_vector * diag);
23 static void update_diag (const gsl_matrix * J, gsl_vector * diag);
24 static double compute_delta (gsl_vector * diag, gsl_vector * x);
25 static void compute_df (const gsl_vector * f_trial, const gsl_vector * f, gsl_vector * df);
26 static void compute_wv (const gsl_vector * qtdf, const gsl_vector *rdx, const gsl_vector *dx, const gsl_vector *diag, double pnorm, gsl_vector * w, gsl_vector * v);
28 static double scaled_enorm (const gsl_vector * d, const gsl_vector * f);
30 static double scaled_enorm (const gsl_vector * d, const gsl_vector * f) {
32 size_t i, n = f->size ;
33 for (i = 0; i < n ; i++) {
34 double fi= gsl_vector_get(f, i);
35 double di= gsl_vector_get(d, i);
42 static double enorm_sum (const gsl_vector * a, const gsl_vector * b);
44 static double enorm_sum (const gsl_vector * a, const gsl_vector * b) {
46 size_t i, n = a->size ;
47 for (i = 0; i < n ; i++) {
48 double ai= gsl_vector_get(a, i);
49 double bi= gsl_vector_get(b, i);
57 compute_wv (const gsl_vector * qtdf, const gsl_vector *rdx, const gsl_vector *dx, const gsl_vector *diag, double pnorm, gsl_vector * w, gsl_vector * v)
59 size_t i, n = qtdf->size;
61 for (i = 0; i < n; i++)
63 double qtdfi = gsl_vector_get (qtdf, i);
64 double rdxi = gsl_vector_get (rdx, i);
65 double dxi = gsl_vector_get (dx, i);
66 double diagi = gsl_vector_get (diag, i);
68 gsl_vector_set (w, i, (qtdfi - rdxi) / pnorm);
69 gsl_vector_set (v, i, diagi * diagi * dxi / pnorm);
75 compute_df (const gsl_vector * f_trial, const gsl_vector * f, gsl_vector * df)
77 size_t i, n = f->size;
79 for (i = 0; i < n; i++)
81 double dfi = gsl_vector_get (f_trial, i) - gsl_vector_get (f, i);
82 gsl_vector_set (df, i, dfi);
87 compute_diag (const gsl_matrix * J, gsl_vector * diag)
89 size_t i, j, n = diag->size;
91 for (j = 0; j < n; j++)
94 for (i = 0; i < n; i++)
96 double Jij = gsl_matrix_get (J, i, j);
102 gsl_vector_set (diag, j, sqrt (sum));
107 update_diag (const gsl_matrix * J, gsl_vector * diag)
109 size_t i, j, n = diag->size;
111 for (j = 0; j < n; j++)
113 double cnorm, diagj, sum = 0;
114 for (i = 0; i < n; i++)
116 double Jij = gsl_matrix_get (J, i, j);
123 diagj = gsl_vector_get (diag, j);
126 gsl_vector_set (diag, j, cnorm);
131 compute_delta (gsl_vector * diag, gsl_vector * x)
133 double Dx = scaled_enorm (diag, x);
136 return (Dx > 0) ? factor * Dx : factor;
140 compute_actual_reduction (double fnorm, double fnorm1)
146 double u = fnorm1 / fnorm;
158 compute_predicted_reduction (double fnorm, double fnorm1)
164 double u = fnorm1 / fnorm;
176 compute_qtf (const gsl_matrix * q, const gsl_vector * f, gsl_vector * qtf)
178 size_t i, j, N = f->size ;
180 for (j = 0; j < N; j++)
183 for (i = 0; i < N; i++)
184 sum += gsl_matrix_get (q, i, j) * gsl_vector_get (f, i);
186 gsl_vector_set (qtf, j, sum);
191 compute_rdx (const gsl_matrix * r, const gsl_vector * dx, gsl_vector * rdx)
193 size_t i, j, N = dx->size ;
195 for (i = 0; i < N; i++)
199 for (j = i; j < N; j++)
201 sum += gsl_matrix_get (r, i, j) * gsl_vector_get (dx, j);
204 gsl_vector_set (rdx, i, sum);
210 compute_trial_step (gsl_vector *x, gsl_vector * dx, gsl_vector * x_trial)
212 size_t i, N = x->size;
214 for (i = 0; i < N; i++)
216 double pi = gsl_vector_get (dx, i);
217 double xi = gsl_vector_get (x, i);
218 gsl_vector_set (x_trial, i, xi + pi);
223 newton_direction (const gsl_matrix * r, const gsl_vector * qtf, gsl_vector * p)
225 const size_t N = r->size2;
229 status = gsl_linalg_R_solve (r, qtf, p);
232 printf("rsolve status = %d\n", status);
235 for (i = 0; i < N; i++)
237 double pi = gsl_vector_get (p, i);
238 gsl_vector_set (p, i, -pi);
245 gradient_direction (const gsl_matrix * r, const gsl_vector * qtf,
246 const gsl_vector * diag, gsl_vector * g)
248 const size_t M = r->size1;
249 const size_t N = r->size2;
253 for (j = 0; j < M; j++)
258 for (i = 0; i < N; i++)
260 sum += gsl_matrix_get (r, i, j) * gsl_vector_get (qtf, i);
263 dj = gsl_vector_get (diag, j);
264 gsl_vector_set (g, j, -sum / dj);
269 minimum_step (double gnorm, const gsl_vector * diag, gsl_vector * g)
271 const size_t N = g->size;
274 for (i = 0; i < N; i++)
276 double gi = gsl_vector_get (g, i);
277 double di = gsl_vector_get (diag, i);
278 gsl_vector_set (g, i, (gi / gnorm) / di);
283 compute_Rg (const gsl_matrix * r, const gsl_vector * gradient, gsl_vector * Rg)
285 const size_t N = r->size2;
288 for (i = 0; i < N; i++)
292 for (j = i; j < N; j++)
294 double gj = gsl_vector_get (gradient, j);
295 double rij = gsl_matrix_get (r, i, j);
299 gsl_vector_set (Rg, i, sum);
304 scaled_addition (double alpha, gsl_vector * newton, double beta, gsl_vector * gradient, gsl_vector * p)
306 const size_t N = p->size;
309 for (i = 0; i < N; i++)
311 double ni = gsl_vector_get (newton, i);
312 double gi = gsl_vector_get (gradient, i);
313 gsl_vector_set (p, i, alpha * ni + beta * gi);
318 dogleg (const gsl_matrix * r, const gsl_vector * qtf,
319 const gsl_vector * diag, double delta,
320 gsl_vector * newton, gsl_vector * gradient, gsl_vector * p)
322 double qnorm, gnorm, sgnorm, bnorm, temp;
324 newton_direction (r, qtf, newton);
327 printf("newton = "); gsl_vector_fprintf(stdout, newton, "%g"); printf("\n");
330 qnorm = scaled_enorm (diag, newton);
334 gsl_vector_memcpy (p, newton);
336 printf("took newton (qnorm = %g <= delta = %g)\n", qnorm, delta);
341 gradient_direction (r, qtf, diag, gradient);
344 printf("grad = "); gsl_vector_fprintf(stdout, gradient, "%g"); printf("\n");
347 gnorm = enorm (gradient);
351 double alpha = delta / qnorm;
353 scaled_addition (alpha, newton, beta, gradient, p);
355 printf("took scaled newton because gnorm = 0\n");
360 minimum_step (gnorm, diag, gradient);
362 compute_Rg (r, gradient, p); /* Use p as temporary space to compute Rg */
365 printf("mingrad = "); gsl_vector_fprintf(stdout, gradient, "%g"); printf("\n");
366 printf("Rg = "); gsl_vector_fprintf(stdout, p, "%g"); printf("\n");
370 sgnorm = (gnorm / temp) / temp;
376 scaled_addition (alpha, newton, beta, gradient, p);
378 printf("took gradient\n");
386 double bg = bnorm / gnorm;
387 double bq = bnorm / qnorm;
388 double dq = delta / qnorm;
389 double dq2 = dq * dq;
390 double sd = sgnorm / delta;
391 double sd2 = sd * sd;
393 double t1 = bg * bq * sd;
395 double t2 = t1 - dq * sd2 + sqrt (u * u + (1-dq2) * (1 - sd2));
397 double alpha = dq * (1 - sd2) / t2;
398 double beta = (1 - alpha) * sgnorm;
401 printf("bnorm = %g\n", bnorm);
402 printf("gnorm = %g\n", gnorm);
403 printf("qnorm = %g\n", qnorm);
404 printf("delta = %g\n", delta);
405 printf("alpha = %g beta = %g\n", alpha, beta);
406 printf("took scaled combination of newton and gradient\n");
409 scaled_addition (alpha, newton, beta, gradient, p);