gsl_multiroots.h 6.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178
  1. /* multiroots/gsl_multiroots.h
  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. #ifndef __GSL_MULTIROOTS_H__
  20. #define __GSL_MULTIROOTS_H__
  21. #include <stdlib.h>
  22. #include "gsl_types.h"
  23. #include "gsl_math.h"
  24. #include "gsl_vector.h"
  25. #include "gsl_matrix.h"
  26. #undef __BEGIN_DECLS
  27. #undef __END_DECLS
  28. #ifdef __cplusplus
  29. # define __BEGIN_DECLS extern "C" {
  30. # define __END_DECLS }
  31. #else
  32. # define __BEGIN_DECLS /* empty */
  33. # define __END_DECLS /* empty */
  34. #endif
  35. __BEGIN_DECLS
  36. /* Definition of vector-valued functions with parameters based on gsl_vector */
  37. struct gsl_multiroot_function_struct
  38. {
  39. int (* f) (const gsl_vector * x, void * params, gsl_vector * f);
  40. size_t n;
  41. void * params;
  42. };
  43. typedef struct gsl_multiroot_function_struct gsl_multiroot_function ;
  44. #define GSL_MULTIROOT_FN_EVAL(F,x,y) (*((F)->f))(x,(F)->params,(y))
  45. int gsl_multiroot_fdjacobian (gsl_multiroot_function * F,
  46. const gsl_vector * x, const gsl_vector * f,
  47. double epsrel, gsl_matrix * jacobian);
  48. typedef struct
  49. {
  50. const char *name;
  51. size_t size;
  52. int (*alloc) (void *state, size_t n);
  53. int (*set) (void *state, gsl_multiroot_function * function, gsl_vector * x, gsl_vector * f, gsl_vector * dx);
  54. int (*iterate) (void *state, gsl_multiroot_function * function, gsl_vector * x, gsl_vector * f, gsl_vector * dx);
  55. void (*free) (void *state);
  56. }
  57. gsl_multiroot_fsolver_type;
  58. typedef struct
  59. {
  60. const gsl_multiroot_fsolver_type * type;
  61. gsl_multiroot_function * function ;
  62. gsl_vector * x ;
  63. gsl_vector * f ;
  64. gsl_vector * dx ;
  65. void *state;
  66. }
  67. gsl_multiroot_fsolver;
  68. gsl_multiroot_fsolver *
  69. gsl_multiroot_fsolver_alloc (const gsl_multiroot_fsolver_type * T,
  70. size_t n);
  71. void gsl_multiroot_fsolver_free (gsl_multiroot_fsolver * s);
  72. int gsl_multiroot_fsolver_set (gsl_multiroot_fsolver * s,
  73. gsl_multiroot_function * f,
  74. const gsl_vector * x);
  75. int gsl_multiroot_fsolver_iterate (gsl_multiroot_fsolver * s);
  76. const char * gsl_multiroot_fsolver_name (const gsl_multiroot_fsolver * s);
  77. gsl_vector * gsl_multiroot_fsolver_root (const gsl_multiroot_fsolver * s);
  78. gsl_vector * gsl_multiroot_fsolver_dx (const gsl_multiroot_fsolver * s);
  79. gsl_vector * gsl_multiroot_fsolver_f (const gsl_multiroot_fsolver * s);
  80. /* Definition of vector-valued functions and gradient with parameters
  81. based on gsl_vector */
  82. struct gsl_multiroot_function_fdf_struct
  83. {
  84. int (* f) (const gsl_vector * x, void * params, gsl_vector * f);
  85. int (* df) (const gsl_vector * x, void * params, gsl_matrix * df);
  86. int (* fdf) (const gsl_vector * x, void * params, gsl_vector * f, gsl_matrix *df);
  87. size_t n;
  88. void * params;
  89. };
  90. typedef struct gsl_multiroot_function_fdf_struct gsl_multiroot_function_fdf ;
  91. #define GSL_MULTIROOT_FN_EVAL_F(F,x,y) ((*((F)->f))(x,(F)->params,(y)))
  92. #define GSL_MULTIROOT_FN_EVAL_DF(F,x,dy) ((*((F)->df))(x,(F)->params,(dy)))
  93. #define GSL_MULTIROOT_FN_EVAL_F_DF(F,x,y,dy) ((*((F)->fdf))(x,(F)->params,(y),(dy)))
  94. typedef struct
  95. {
  96. const char *name;
  97. size_t size;
  98. int (*alloc) (void *state, size_t n);
  99. int (*set) (void *state, gsl_multiroot_function_fdf * fdf, gsl_vector * x, gsl_vector * f, gsl_matrix * J, gsl_vector * dx);
  100. int (*iterate) (void *state, gsl_multiroot_function_fdf * fdf, gsl_vector * x, gsl_vector * f, gsl_matrix * J, gsl_vector * dx);
  101. void (*free) (void *state);
  102. }
  103. gsl_multiroot_fdfsolver_type;
  104. typedef struct
  105. {
  106. const gsl_multiroot_fdfsolver_type * type;
  107. gsl_multiroot_function_fdf * fdf ;
  108. gsl_vector * x;
  109. gsl_vector * f;
  110. gsl_matrix * J;
  111. gsl_vector * dx;
  112. void *state;
  113. }
  114. gsl_multiroot_fdfsolver;
  115. gsl_multiroot_fdfsolver *
  116. gsl_multiroot_fdfsolver_alloc (const gsl_multiroot_fdfsolver_type * T,
  117. size_t n);
  118. int
  119. gsl_multiroot_fdfsolver_set (gsl_multiroot_fdfsolver * s,
  120. gsl_multiroot_function_fdf * fdf,
  121. const gsl_vector * x);
  122. int
  123. gsl_multiroot_fdfsolver_iterate (gsl_multiroot_fdfsolver * s);
  124. void
  125. gsl_multiroot_fdfsolver_free (gsl_multiroot_fdfsolver * s);
  126. const char * gsl_multiroot_fdfsolver_name (const gsl_multiroot_fdfsolver * s);
  127. gsl_vector * gsl_multiroot_fdfsolver_root (const gsl_multiroot_fdfsolver * s);
  128. gsl_vector * gsl_multiroot_fdfsolver_dx (const gsl_multiroot_fdfsolver * s);
  129. gsl_vector * gsl_multiroot_fdfsolver_f (const gsl_multiroot_fdfsolver * s);
  130. int gsl_multiroot_test_delta (const gsl_vector * dx, const gsl_vector * x,
  131. double epsabs, double epsrel);
  132. int gsl_multiroot_test_residual (const gsl_vector * f, double epsabs);
  133. GSL_VAR const gsl_multiroot_fsolver_type * gsl_multiroot_fsolver_dnewton;
  134. GSL_VAR const gsl_multiroot_fsolver_type * gsl_multiroot_fsolver_broyden;
  135. GSL_VAR const gsl_multiroot_fsolver_type * gsl_multiroot_fsolver_hybrid;
  136. GSL_VAR const gsl_multiroot_fsolver_type * gsl_multiroot_fsolver_hybrids;
  137. GSL_VAR const gsl_multiroot_fdfsolver_type * gsl_multiroot_fdfsolver_newton;
  138. GSL_VAR const gsl_multiroot_fdfsolver_type * gsl_multiroot_fdfsolver_gnewton;
  139. GSL_VAR const gsl_multiroot_fdfsolver_type * gsl_multiroot_fdfsolver_hybridj;
  140. GSL_VAR const gsl_multiroot_fdfsolver_type * gsl_multiroot_fdfsolver_hybridsj;
  141. __END_DECLS
  142. #endif /* __GSL_MULTIROOTS_H__ */