1 /* multimin/gsl_multimin.h
3 * Copyright (C) 1996, 1997, 1998, 1999, 2000 Fabrice Rossi
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 /* Modified by Tuomo Keskitalo to include fminimizer and
21 Nelder Mead related lines */
23 #ifndef __GSL_MULTIMIN_H__
24 #define __GSL_MULTIMIN_H__
27 #include <gsl/gsl_types.h>
28 #include <gsl/gsl_math.h>
29 #include <gsl/gsl_vector.h>
30 #include <gsl/gsl_matrix.h>
31 #include <gsl/gsl_min.h>
36 # define __BEGIN_DECLS extern "C" {
37 # define __END_DECLS }
39 # define __BEGIN_DECLS /* empty */
40 # define __END_DECLS /* empty */
45 /* Definition of an arbitrary real-valued function with gsl_vector input and */
47 struct gsl_multimin_function_struct
49 double (* f) (const gsl_vector * x, void * params);
54 typedef struct gsl_multimin_function_struct gsl_multimin_function;
56 #define GSL_MULTIMIN_FN_EVAL(F,x) (*((F)->f))(x,(F)->params)
58 /* Definition of an arbitrary differentiable real-valued function */
59 /* with gsl_vector input and parameters */
60 struct gsl_multimin_function_fdf_struct
62 double (* f) (const gsl_vector * x, void * params);
63 void (* df) (const gsl_vector * x, void * params,gsl_vector * df);
64 void (* fdf) (const gsl_vector * x, void * params,double *f,gsl_vector * df);
69 typedef struct gsl_multimin_function_fdf_struct gsl_multimin_function_fdf;
71 #define GSL_MULTIMIN_FN_EVAL_F(F,x) (*((F)->f))(x,(F)->params)
72 #define GSL_MULTIMIN_FN_EVAL_DF(F,x,g) (*((F)->df))(x,(F)->params,(g))
73 #define GSL_MULTIMIN_FN_EVAL_F_DF(F,x,y,g) (*((F)->fdf))(x,(F)->params,(y),(g))
75 int gsl_multimin_diff (const gsl_multimin_function * f,
76 const gsl_vector * x, gsl_vector * g);
78 /* minimization of non-differentiable functions */
84 int (*alloc) (void *state, size_t n);
85 int (*set) (void *state, gsl_multimin_function * f,
88 const gsl_vector * step_size);
89 int (*iterate) (void *state, gsl_multimin_function * f,
93 void (*free) (void *state);
95 gsl_multimin_fminimizer_type;
99 /* multi dimensional part */
100 const gsl_multimin_fminimizer_type *type;
101 gsl_multimin_function *f;
110 gsl_multimin_fminimizer;
112 gsl_multimin_fminimizer *
113 gsl_multimin_fminimizer_alloc(const gsl_multimin_fminimizer_type *T,
117 gsl_multimin_fminimizer_set (gsl_multimin_fminimizer * s,
118 gsl_multimin_function * f,
119 const gsl_vector * x,
120 const gsl_vector * step_size);
123 gsl_multimin_fminimizer_free(gsl_multimin_fminimizer *s);
126 gsl_multimin_fminimizer_name (const gsl_multimin_fminimizer * s);
129 gsl_multimin_fminimizer_iterate(gsl_multimin_fminimizer *s);
132 gsl_multimin_fminimizer_x (const gsl_multimin_fminimizer * s);
135 gsl_multimin_fminimizer_minimum (const gsl_multimin_fminimizer * s);
138 gsl_multimin_fminimizer_size (const gsl_multimin_fminimizer * s);
140 /* Convergence test functions */
143 gsl_multimin_test_gradient(const gsl_vector * g,double epsabs);
146 gsl_multimin_test_size(const double size ,double epsabs);
148 /* minimisation of differentiable functions */
154 int (*alloc) (void *state, size_t n);
155 int (*set) (void *state, gsl_multimin_function_fdf * fdf,
156 const gsl_vector * x, double * f,
157 gsl_vector * gradient, double step_size, double tol);
158 int (*iterate) (void *state,gsl_multimin_function_fdf * fdf,
159 gsl_vector * x, double * f,
160 gsl_vector * gradient, gsl_vector * dx);
161 int (*restart) (void *state);
162 void (*free) (void *state);
164 gsl_multimin_fdfminimizer_type;
168 /* multi dimensional part */
169 const gsl_multimin_fdfminimizer_type *type;
170 gsl_multimin_function_fdf *fdf;
174 gsl_vector * gradient;
179 gsl_multimin_fdfminimizer;
181 gsl_multimin_fdfminimizer *
182 gsl_multimin_fdfminimizer_alloc(const gsl_multimin_fdfminimizer_type *T,
186 gsl_multimin_fdfminimizer_set (gsl_multimin_fdfminimizer * s,
187 gsl_multimin_function_fdf *fdf,
188 const gsl_vector * x,
189 double step_size, double tol);
192 gsl_multimin_fdfminimizer_free(gsl_multimin_fdfminimizer *s);
195 gsl_multimin_fdfminimizer_name (const gsl_multimin_fdfminimizer * s);
198 gsl_multimin_fdfminimizer_iterate(gsl_multimin_fdfminimizer *s);
201 gsl_multimin_fdfminimizer_restart(gsl_multimin_fdfminimizer *s);
204 gsl_multimin_fdfminimizer_x (const gsl_multimin_fdfminimizer * s);
207 gsl_multimin_fdfminimizer_dx (const gsl_multimin_fdfminimizer * s);
210 gsl_multimin_fdfminimizer_gradient (const gsl_multimin_fdfminimizer * s);
213 gsl_multimin_fdfminimizer_minimum (const gsl_multimin_fdfminimizer * s);
215 GSL_VAR const gsl_multimin_fdfminimizer_type *gsl_multimin_fdfminimizer_steepest_descent;
216 GSL_VAR const gsl_multimin_fdfminimizer_type *gsl_multimin_fdfminimizer_conjugate_pr;
217 GSL_VAR const gsl_multimin_fdfminimizer_type *gsl_multimin_fdfminimizer_conjugate_fr;
218 GSL_VAR const gsl_multimin_fdfminimizer_type *gsl_multimin_fdfminimizer_vector_bfgs;
219 GSL_VAR const gsl_multimin_fdfminimizer_type *gsl_multimin_fdfminimizer_vector_bfgs2;
220 GSL_VAR const gsl_multimin_fminimizer_type *gsl_multimin_fminimizer_nmsimplex;
225 #endif /* __GSL_MULTIMIN_H__ */