gsl_multiroots__fdjac.c 2.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112
  1. /* multiroots/fdjac.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 "gsl_multiroots.h"
  21. int
  22. gsl_multiroot_fdjacobian (gsl_multiroot_function * F,
  23. const gsl_vector * x, const gsl_vector * f,
  24. double epsrel, gsl_matrix * jacobian)
  25. {
  26. const size_t n = x->size;
  27. const size_t m = f->size;
  28. const size_t n1 = jacobian->size1;
  29. const size_t n2 = jacobian->size2;
  30. int status = 0;
  31. if (m != n1 || n != n2)
  32. {
  33. GSL_ERROR ("function and jacobian are not conformant", GSL_EBADLEN);
  34. }
  35. {
  36. size_t i,j;
  37. gsl_vector *x1, *f1;
  38. x1 = gsl_vector_alloc (n);
  39. if (x1 == 0)
  40. {
  41. GSL_ERROR ("failed to allocate space for x1 workspace", GSL_ENOMEM);
  42. }
  43. f1 = gsl_vector_alloc (m);
  44. if (f1 == 0)
  45. {
  46. gsl_vector_free (x1);
  47. GSL_ERROR ("failed to allocate space for f1 workspace", GSL_ENOMEM);
  48. }
  49. gsl_vector_memcpy (x1, x); /* copy x into x1 */
  50. for (j = 0; j < n; j++)
  51. {
  52. double xj = gsl_vector_get (x, j);
  53. double dx = epsrel * fabs (xj);
  54. if (dx == 0)
  55. {
  56. dx = epsrel;
  57. }
  58. gsl_vector_set (x1, j, xj + dx);
  59. {
  60. int f_stat = GSL_MULTIROOT_FN_EVAL (F, x1, f1);
  61. if (f_stat != GSL_SUCCESS)
  62. {
  63. status = GSL_EBADFUNC;
  64. break; /* n.b. avoid memory leak for x1,f1 */
  65. }
  66. }
  67. gsl_vector_set (x1, j, xj);
  68. for (i = 0; i < m; i++)
  69. {
  70. double g1 = gsl_vector_get (f1, i);
  71. double g0 = gsl_vector_get (f, i);
  72. gsl_matrix_set (jacobian, i, j, (g1 - g0) / dx);
  73. }
  74. {
  75. gsl_vector_view col = gsl_matrix_column (jacobian, j);
  76. int null_col = gsl_vector_isnull (&col.vector);
  77. /* if column is null, return an error - this may be due to
  78. dx being too small. Try increasing epsrel */
  79. if (null_col) {
  80. status = GSL_ESING;
  81. }
  82. }
  83. }
  84. gsl_vector_free (x1);
  85. gsl_vector_free (f1);
  86. }
  87. if (status)
  88. return status;
  89. else
  90. return GSL_SUCCESS;
  91. }