gsl_multiroots__fdfsolver.c 4.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177
  1. /* multiroots/fdfsolver.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 <stdlib.h>
  21. #include <string.h>
  22. #include "gsl_errno.h"
  23. #include "gsl_multiroots.h"
  24. gsl_multiroot_fdfsolver *
  25. gsl_multiroot_fdfsolver_alloc (const gsl_multiroot_fdfsolver_type * T,
  26. size_t n)
  27. {
  28. int status;
  29. gsl_multiroot_fdfsolver * s;
  30. s = (gsl_multiroot_fdfsolver *) malloc (sizeof (gsl_multiroot_fdfsolver));
  31. if (s == 0)
  32. {
  33. GSL_ERROR_VAL ("failed to allocate space for multiroot solver struct",
  34. GSL_ENOMEM, 0);
  35. }
  36. s->x = gsl_vector_calloc (n);
  37. if (s->x == 0)
  38. {
  39. free (s);
  40. GSL_ERROR_VAL ("failed to allocate space for x", GSL_ENOMEM, 0);
  41. }
  42. s->f = gsl_vector_calloc (n);
  43. if (s->f == 0)
  44. {
  45. gsl_vector_free (s->x);
  46. free (s);
  47. GSL_ERROR_VAL ("failed to allocate space for f", GSL_ENOMEM, 0);
  48. }
  49. s->J = gsl_matrix_calloc (n,n);
  50. if (s->J == 0)
  51. {
  52. gsl_vector_free (s->x);
  53. gsl_vector_free (s->f);
  54. free (s);
  55. GSL_ERROR_VAL ("failed to allocate space for g", GSL_ENOMEM, 0);
  56. }
  57. s->dx = gsl_vector_calloc (n);
  58. if (s->dx == 0)
  59. {
  60. gsl_matrix_free (s->J);
  61. gsl_vector_free (s->x);
  62. gsl_vector_free (s->f);
  63. free (s);
  64. GSL_ERROR_VAL ("failed to allocate space for dx", GSL_ENOMEM, 0);
  65. }
  66. s->state = malloc (T->size);
  67. if (s->state == 0)
  68. {
  69. gsl_vector_free (s->dx);
  70. gsl_vector_free (s->x);
  71. gsl_vector_free (s->f);
  72. gsl_matrix_free (s->J);
  73. free (s); /* exception in constructor, avoid memory leak */
  74. GSL_ERROR_VAL ("failed to allocate space for multiroot solver state",
  75. GSL_ENOMEM, 0);
  76. }
  77. s->type = T ;
  78. status = (s->type->alloc)(s->state, n);
  79. if (status != GSL_SUCCESS)
  80. {
  81. free (s->state);
  82. gsl_vector_free (s->dx);
  83. gsl_vector_free (s->x);
  84. gsl_vector_free (s->f);
  85. gsl_matrix_free (s->J);
  86. free (s); /* exception in constructor, avoid memory leak */
  87. GSL_ERROR_VAL ("failed to set solver", status, 0);
  88. }
  89. s->fdf = NULL;
  90. return s;
  91. }
  92. int
  93. gsl_multiroot_fdfsolver_set (gsl_multiroot_fdfsolver * s,
  94. gsl_multiroot_function_fdf * f,
  95. const gsl_vector * x)
  96. {
  97. if (s->x->size != f->n)
  98. {
  99. GSL_ERROR ("function incompatible with solver size", GSL_EBADLEN);
  100. }
  101. if (x->size != f->n)
  102. {
  103. GSL_ERROR ("vector length not compatible with function", GSL_EBADLEN);
  104. }
  105. s->fdf = f;
  106. gsl_vector_memcpy(s->x,x);
  107. return (s->type->set) (s->state, s->fdf, s->x, s->f, s->J, s->dx);
  108. }
  109. int
  110. gsl_multiroot_fdfsolver_iterate (gsl_multiroot_fdfsolver * s)
  111. {
  112. return (s->type->iterate) (s->state, s->fdf, s->x, s->f, s->J, s->dx);
  113. }
  114. void
  115. gsl_multiroot_fdfsolver_free (gsl_multiroot_fdfsolver * s)
  116. {
  117. (s->type->free) (s->state);
  118. free (s->state);
  119. gsl_vector_free (s->dx);
  120. gsl_vector_free (s->x);
  121. gsl_vector_free (s->f);
  122. gsl_matrix_free (s->J);
  123. free (s);
  124. }
  125. const char *
  126. gsl_multiroot_fdfsolver_name (const gsl_multiroot_fdfsolver * s)
  127. {
  128. return s->type->name;
  129. }
  130. gsl_vector *
  131. gsl_multiroot_fdfsolver_root (const gsl_multiroot_fdfsolver * s)
  132. {
  133. return s->x;
  134. }
  135. gsl_vector *
  136. gsl_multiroot_fdfsolver_dx (const gsl_multiroot_fdfsolver * s)
  137. {
  138. return s->dx;
  139. }
  140. gsl_vector *
  141. gsl_multiroot_fdfsolver_f (const gsl_multiroot_fdfsolver * s)
  142. {
  143. return s->f;
  144. }