1 /* multiroots/gsl_multiroots.h
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.
20 #ifndef __GSL_MULTIROOTS_H__
21 #define __GSL_MULTIROOTS_H__
24 #include <gsl/gsl_types.h>
25 #include <gsl/gsl_math.h>
26 #include <gsl/gsl_vector.h>
27 #include <gsl/gsl_matrix.h>
32 # define __BEGIN_DECLS extern "C" {
33 # define __END_DECLS }
35 # define __BEGIN_DECLS /* empty */
36 # define __END_DECLS /* empty */
41 /* Definition of vector-valued functions with parameters based on gsl_vector */
43 struct gsl_multiroot_function_struct
45 int (* f) (const gsl_vector * x, void * params, gsl_vector * f);
50 typedef struct gsl_multiroot_function_struct gsl_multiroot_function ;
52 #define GSL_MULTIROOT_FN_EVAL(F,x,y) (*((F)->f))(x,(F)->params,(y))
54 int gsl_multiroot_fdjacobian (gsl_multiroot_function * F,
55 const gsl_vector * x, const gsl_vector * f,
56 double epsrel, gsl_matrix * jacobian);
63 int (*alloc) (void *state, size_t n);
64 int (*set) (void *state, gsl_multiroot_function * function, gsl_vector * x, gsl_vector * f, gsl_vector * dx);
65 int (*iterate) (void *state, gsl_multiroot_function * function, gsl_vector * x, gsl_vector * f, gsl_vector * dx);
66 void (*free) (void *state);
68 gsl_multiroot_fsolver_type;
72 const gsl_multiroot_fsolver_type * type;
73 gsl_multiroot_function * function ;
79 gsl_multiroot_fsolver;
81 gsl_multiroot_fsolver *
82 gsl_multiroot_fsolver_alloc (const gsl_multiroot_fsolver_type * T,
85 void gsl_multiroot_fsolver_free (gsl_multiroot_fsolver * s);
87 int gsl_multiroot_fsolver_set (gsl_multiroot_fsolver * s,
88 gsl_multiroot_function * f,
89 const gsl_vector * x);
91 int gsl_multiroot_fsolver_iterate (gsl_multiroot_fsolver * s);
93 const char * gsl_multiroot_fsolver_name (const gsl_multiroot_fsolver * s);
94 gsl_vector * gsl_multiroot_fsolver_root (const gsl_multiroot_fsolver * s);
95 gsl_vector * gsl_multiroot_fsolver_dx (const gsl_multiroot_fsolver * s);
96 gsl_vector * gsl_multiroot_fsolver_f (const gsl_multiroot_fsolver * s);
98 /* Definition of vector-valued functions and gradient with parameters
99 based on gsl_vector */
101 struct gsl_multiroot_function_fdf_struct
103 int (* f) (const gsl_vector * x, void * params, gsl_vector * f);
104 int (* df) (const gsl_vector * x, void * params, gsl_matrix * df);
105 int (* fdf) (const gsl_vector * x, void * params, gsl_vector * f, gsl_matrix *df);
110 typedef struct gsl_multiroot_function_fdf_struct gsl_multiroot_function_fdf ;
112 #define GSL_MULTIROOT_FN_EVAL_F(F,x,y) ((*((F)->f))(x,(F)->params,(y)))
113 #define GSL_MULTIROOT_FN_EVAL_DF(F,x,dy) ((*((F)->df))(x,(F)->params,(dy)))
114 #define GSL_MULTIROOT_FN_EVAL_F_DF(F,x,y,dy) ((*((F)->fdf))(x,(F)->params,(y),(dy)))
120 int (*alloc) (void *state, size_t n);
121 int (*set) (void *state, gsl_multiroot_function_fdf * fdf, gsl_vector * x, gsl_vector * f, gsl_matrix * J, gsl_vector * dx);
122 int (*iterate) (void *state, gsl_multiroot_function_fdf * fdf, gsl_vector * x, gsl_vector * f, gsl_matrix * J, gsl_vector * dx);
123 void (*free) (void *state);
125 gsl_multiroot_fdfsolver_type;
129 const gsl_multiroot_fdfsolver_type * type;
130 gsl_multiroot_function_fdf * fdf ;
137 gsl_multiroot_fdfsolver;
139 gsl_multiroot_fdfsolver *
140 gsl_multiroot_fdfsolver_alloc (const gsl_multiroot_fdfsolver_type * T,
144 gsl_multiroot_fdfsolver_set (gsl_multiroot_fdfsolver * s,
145 gsl_multiroot_function_fdf * fdf,
146 const gsl_vector * x);
149 gsl_multiroot_fdfsolver_iterate (gsl_multiroot_fdfsolver * s);
152 gsl_multiroot_fdfsolver_free (gsl_multiroot_fdfsolver * s);
154 const char * gsl_multiroot_fdfsolver_name (const gsl_multiroot_fdfsolver * s);
155 gsl_vector * gsl_multiroot_fdfsolver_root (const gsl_multiroot_fdfsolver * s);
156 gsl_vector * gsl_multiroot_fdfsolver_dx (const gsl_multiroot_fdfsolver * s);
157 gsl_vector * gsl_multiroot_fdfsolver_f (const gsl_multiroot_fdfsolver * s);
159 int gsl_multiroot_test_delta (const gsl_vector * dx, const gsl_vector * x,
160 double epsabs, double epsrel);
162 int gsl_multiroot_test_residual (const gsl_vector * f, double epsabs);
164 GSL_VAR const gsl_multiroot_fsolver_type * gsl_multiroot_fsolver_dnewton;
165 GSL_VAR const gsl_multiroot_fsolver_type * gsl_multiroot_fsolver_broyden;
166 GSL_VAR const gsl_multiroot_fsolver_type * gsl_multiroot_fsolver_hybrid;
167 GSL_VAR const gsl_multiroot_fsolver_type * gsl_multiroot_fsolver_hybrids;
169 GSL_VAR const gsl_multiroot_fdfsolver_type * gsl_multiroot_fdfsolver_newton;
170 GSL_VAR const gsl_multiroot_fdfsolver_type * gsl_multiroot_fdfsolver_gnewton;
171 GSL_VAR const gsl_multiroot_fdfsolver_type * gsl_multiroot_fdfsolver_hybridj;
172 GSL_VAR const gsl_multiroot_fdfsolver_type * gsl_multiroot_fdfsolver_hybridsj;
177 #endif /* __GSL_MULTIROOTS_H__ */