gsl_multiroots__gnewton.c 5.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227
  1. /* multiroots/gnewton.c
  2. *
  3. * Copyright (C) 1996, 1997, 1998, 1999, 2000, 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 <stddef.h>
  21. #include <stdlib.h>
  22. #include <stdio.h>
  23. #include <math.h>
  24. #include <float.h>
  25. #include "gsl_math.h"
  26. #include "gsl_errno.h"
  27. #include "gsl_multiroots.h"
  28. #include "gsl_linalg.h"
  29. #include "gsl_multiroots__enorm.c"
  30. /* Simple globally convergent Newton method (rejects uphill steps) */
  31. typedef struct
  32. {
  33. double phi;
  34. gsl_vector * x_trial;
  35. gsl_vector * d;
  36. gsl_matrix * lu;
  37. gsl_permutation * permutation;
  38. }
  39. gnewton_state_t;
  40. static int gnewton_alloc (void * vstate, size_t n);
  41. static int gnewton_set (void * vstate, gsl_multiroot_function_fdf * fdf, gsl_vector * x, gsl_vector * f, gsl_matrix * J, gsl_vector * dx);
  42. static int gnewton_iterate (void * vstate, gsl_multiroot_function_fdf * fdf, gsl_vector * x, gsl_vector * f, gsl_matrix * J, gsl_vector * dx);
  43. static void gnewton_free (void * vstate);
  44. static int
  45. gnewton_alloc (void * vstate, size_t n)
  46. {
  47. gnewton_state_t * state = (gnewton_state_t *) vstate;
  48. gsl_vector * d, * x_trial ;
  49. gsl_permutation * p;
  50. gsl_matrix * m;
  51. m = gsl_matrix_calloc (n,n);
  52. if (m == 0)
  53. {
  54. GSL_ERROR ("failed to allocate space for lu", GSL_ENOMEM);
  55. }
  56. state->lu = m ;
  57. p = gsl_permutation_calloc (n);
  58. if (p == 0)
  59. {
  60. gsl_matrix_free(m);
  61. GSL_ERROR ("failed to allocate space for permutation", GSL_ENOMEM);
  62. }
  63. state->permutation = p ;
  64. d = gsl_vector_calloc (n);
  65. if (d == 0)
  66. {
  67. gsl_permutation_free(p);
  68. gsl_matrix_free(m);
  69. GSL_ERROR ("failed to allocate space for d", GSL_ENOMEM);
  70. }
  71. state->d = d;
  72. x_trial = gsl_vector_calloc (n);
  73. if (x_trial == 0)
  74. {
  75. gsl_vector_free(d);
  76. gsl_permutation_free(p);
  77. gsl_matrix_free(m);
  78. GSL_ERROR ("failed to allocate space for x_trial", GSL_ENOMEM);
  79. }
  80. state->x_trial = x_trial;
  81. return GSL_SUCCESS;
  82. }
  83. static int
  84. gnewton_set (void * vstate, gsl_multiroot_function_fdf * FDF, gsl_vector * x, gsl_vector * f, gsl_matrix * J, gsl_vector * dx)
  85. {
  86. gnewton_state_t * state = (gnewton_state_t *) vstate;
  87. size_t i, n = FDF->n ;
  88. GSL_MULTIROOT_FN_EVAL_F_DF (FDF, x, f, J);
  89. for (i = 0; i < n; i++)
  90. {
  91. gsl_vector_set (dx, i, 0.0);
  92. }
  93. state->phi = enorm(f);
  94. return GSL_SUCCESS;
  95. }
  96. static int
  97. gnewton_iterate (void * vstate, gsl_multiroot_function_fdf * fdf, gsl_vector * x, gsl_vector * f, gsl_matrix * J, gsl_vector * dx)
  98. {
  99. gnewton_state_t * state = (gnewton_state_t *) vstate;
  100. int signum ;
  101. double t, phi0, phi1;
  102. size_t i;
  103. size_t n = fdf->n ;
  104. gsl_matrix_memcpy (state->lu, J);
  105. gsl_linalg_LU_decomp (state->lu, state->permutation, &signum);
  106. gsl_linalg_LU_solve (state->lu, state->permutation, f, state->d);
  107. t = 1;
  108. phi0 = state->phi;
  109. new_step:
  110. for (i = 0; i < n; i++)
  111. {
  112. double di = gsl_vector_get (state->d, i);
  113. double xi = gsl_vector_get (x, i);
  114. gsl_vector_set (state->x_trial, i, xi - t*di);
  115. }
  116. {
  117. int status = GSL_MULTIROOT_FN_EVAL_F (fdf, state->x_trial, f);
  118. if (status != GSL_SUCCESS)
  119. {
  120. return GSL_EBADFUNC;
  121. }
  122. }
  123. phi1 = enorm (f);
  124. if (phi1 > phi0 && t > GSL_DBL_EPSILON)
  125. {
  126. /* full step goes uphill, take a reduced step instead */
  127. double theta = phi1 / phi0;
  128. double u = (sqrt(1.0 + 6.0 * theta) - 1.0) / (3.0 * theta);
  129. t *= u ;
  130. goto new_step;
  131. }
  132. /* copy x_trial into x */
  133. gsl_vector_memcpy (x, state->x_trial);
  134. for (i = 0; i < n; i++)
  135. {
  136. double di = gsl_vector_get (state->d, i);
  137. gsl_vector_set (dx, i, -t*di);
  138. }
  139. {
  140. int status = GSL_MULTIROOT_FN_EVAL_DF (fdf, x, J);
  141. if (status != GSL_SUCCESS)
  142. {
  143. return GSL_EBADFUNC;
  144. }
  145. }
  146. state->phi = phi1;
  147. return GSL_SUCCESS;
  148. }
  149. static void
  150. gnewton_free (void * vstate)
  151. {
  152. gnewton_state_t * state = (gnewton_state_t *) vstate;
  153. gsl_vector_free(state->d);
  154. gsl_vector_free(state->x_trial);
  155. gsl_matrix_free(state->lu);
  156. gsl_permutation_free(state->permutation);
  157. }
  158. static const gsl_multiroot_fdfsolver_type gnewton_type =
  159. {"gnewton", /* name */
  160. sizeof (gnewton_state_t),
  161. &gnewton_alloc,
  162. &gnewton_set,
  163. &gnewton_iterate,
  164. &gnewton_free};
  165. const gsl_multiroot_fdfsolver_type * gsl_multiroot_fdfsolver_gnewton = &gnewton_type;