File indexing completed on 2025-02-21 10:03:51
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023 #ifndef __GSL_MULTIMIN_H__
0024 #define __GSL_MULTIMIN_H__
0025
0026 #include <stdlib.h>
0027 #include <gsl/gsl_types.h>
0028 #include <gsl/gsl_math.h>
0029 #include <gsl/gsl_vector.h>
0030 #include <gsl/gsl_matrix.h>
0031 #include <gsl/gsl_min.h>
0032
0033 #undef __BEGIN_DECLS
0034 #undef __END_DECLS
0035 #ifdef __cplusplus
0036 # define __BEGIN_DECLS extern "C" {
0037 # define __END_DECLS }
0038 #else
0039 # define __BEGIN_DECLS
0040 # define __END_DECLS
0041 #endif
0042
0043 __BEGIN_DECLS
0044
0045
0046
0047 struct gsl_multimin_function_struct
0048 {
0049 double (* f) (const gsl_vector * x, void * params);
0050 size_t n;
0051 void * params;
0052 };
0053
0054 typedef struct gsl_multimin_function_struct gsl_multimin_function;
0055
0056 #define GSL_MULTIMIN_FN_EVAL(F,x) (*((F)->f))(x,(F)->params)
0057
0058
0059
0060 struct gsl_multimin_function_fdf_struct
0061 {
0062 double (* f) (const gsl_vector * x, void * params);
0063 void (* df) (const gsl_vector * x, void * params,gsl_vector * df);
0064 void (* fdf) (const gsl_vector * x, void * params,double *f,gsl_vector * df);
0065 size_t n;
0066 void * params;
0067 };
0068
0069 typedef struct gsl_multimin_function_fdf_struct gsl_multimin_function_fdf;
0070
0071 #define GSL_MULTIMIN_FN_EVAL_F(F,x) (*((F)->f))(x,(F)->params)
0072 #define GSL_MULTIMIN_FN_EVAL_DF(F,x,g) (*((F)->df))(x,(F)->params,(g))
0073 #define GSL_MULTIMIN_FN_EVAL_F_DF(F,x,y,g) (*((F)->fdf))(x,(F)->params,(y),(g))
0074
0075 int gsl_multimin_diff (const gsl_multimin_function * f,
0076 const gsl_vector * x, gsl_vector * g);
0077
0078
0079
0080 typedef struct
0081 {
0082 const char *name;
0083 size_t size;
0084 int (*alloc) (void *state, size_t n);
0085 int (*set) (void *state, gsl_multimin_function * f,
0086 const gsl_vector * x,
0087 double * size,
0088 const gsl_vector * step_size);
0089 int (*iterate) (void *state, gsl_multimin_function * f,
0090 gsl_vector * x,
0091 double * size,
0092 double * fval);
0093 void (*free) (void *state);
0094 }
0095 gsl_multimin_fminimizer_type;
0096
0097 typedef struct
0098 {
0099
0100 const gsl_multimin_fminimizer_type *type;
0101 gsl_multimin_function *f;
0102
0103 double fval;
0104 gsl_vector * x;
0105
0106 double size;
0107
0108 void *state;
0109 }
0110 gsl_multimin_fminimizer;
0111
0112 gsl_multimin_fminimizer *
0113 gsl_multimin_fminimizer_alloc(const gsl_multimin_fminimizer_type *T,
0114 size_t n);
0115
0116 int
0117 gsl_multimin_fminimizer_set (gsl_multimin_fminimizer * s,
0118 gsl_multimin_function * f,
0119 const gsl_vector * x,
0120 const gsl_vector * step_size);
0121
0122 void
0123 gsl_multimin_fminimizer_free(gsl_multimin_fminimizer *s);
0124
0125 const char *
0126 gsl_multimin_fminimizer_name (const gsl_multimin_fminimizer * s);
0127
0128 int
0129 gsl_multimin_fminimizer_iterate(gsl_multimin_fminimizer *s);
0130
0131 gsl_vector *
0132 gsl_multimin_fminimizer_x (const gsl_multimin_fminimizer * s);
0133
0134 double
0135 gsl_multimin_fminimizer_minimum (const gsl_multimin_fminimizer * s);
0136
0137 double
0138 gsl_multimin_fminimizer_size (const gsl_multimin_fminimizer * s);
0139
0140
0141
0142 int
0143 gsl_multimin_test_gradient(const gsl_vector * g, double epsabs);
0144
0145 int
0146 gsl_multimin_test_size(const double size, double epsabs);
0147
0148
0149
0150 typedef struct
0151 {
0152 const char *name;
0153 size_t size;
0154 int (*alloc) (void *state, size_t n);
0155 int (*set) (void *state, gsl_multimin_function_fdf * fdf,
0156 const gsl_vector * x, double * f,
0157 gsl_vector * gradient, double step_size, double tol);
0158 int (*iterate) (void *state,gsl_multimin_function_fdf * fdf,
0159 gsl_vector * x, double * f,
0160 gsl_vector * gradient, gsl_vector * dx);
0161 int (*restart) (void *state);
0162 void (*free) (void *state);
0163 }
0164 gsl_multimin_fdfminimizer_type;
0165
0166 typedef struct
0167 {
0168
0169 const gsl_multimin_fdfminimizer_type *type;
0170 gsl_multimin_function_fdf *fdf;
0171
0172 double f;
0173 gsl_vector * x;
0174 gsl_vector * gradient;
0175 gsl_vector * dx;
0176
0177 void *state;
0178 }
0179 gsl_multimin_fdfminimizer;
0180
0181 gsl_multimin_fdfminimizer *
0182 gsl_multimin_fdfminimizer_alloc(const gsl_multimin_fdfminimizer_type *T,
0183 size_t n);
0184
0185 int
0186 gsl_multimin_fdfminimizer_set (gsl_multimin_fdfminimizer * s,
0187 gsl_multimin_function_fdf *fdf,
0188 const gsl_vector * x,
0189 double step_size, double tol);
0190
0191 void
0192 gsl_multimin_fdfminimizer_free(gsl_multimin_fdfminimizer *s);
0193
0194 const char *
0195 gsl_multimin_fdfminimizer_name (const gsl_multimin_fdfminimizer * s);
0196
0197 int
0198 gsl_multimin_fdfminimizer_iterate(gsl_multimin_fdfminimizer *s);
0199
0200 int
0201 gsl_multimin_fdfminimizer_restart(gsl_multimin_fdfminimizer *s);
0202
0203 gsl_vector *
0204 gsl_multimin_fdfminimizer_x (const gsl_multimin_fdfminimizer * s);
0205
0206 gsl_vector *
0207 gsl_multimin_fdfminimizer_dx (const gsl_multimin_fdfminimizer * s);
0208
0209 gsl_vector *
0210 gsl_multimin_fdfminimizer_gradient (const gsl_multimin_fdfminimizer * s);
0211
0212 double
0213 gsl_multimin_fdfminimizer_minimum (const gsl_multimin_fdfminimizer * s);
0214
0215 GSL_VAR const gsl_multimin_fdfminimizer_type *gsl_multimin_fdfminimizer_steepest_descent;
0216 GSL_VAR const gsl_multimin_fdfminimizer_type *gsl_multimin_fdfminimizer_conjugate_pr;
0217 GSL_VAR const gsl_multimin_fdfminimizer_type *gsl_multimin_fdfminimizer_conjugate_fr;
0218 GSL_VAR const gsl_multimin_fdfminimizer_type *gsl_multimin_fdfminimizer_vector_bfgs;
0219 GSL_VAR const gsl_multimin_fdfminimizer_type *gsl_multimin_fdfminimizer_vector_bfgs2;
0220 GSL_VAR const gsl_multimin_fminimizer_type *gsl_multimin_fminimizer_nmsimplex;
0221 GSL_VAR const gsl_multimin_fminimizer_type *gsl_multimin_fminimizer_nmsimplex2;
0222 GSL_VAR const gsl_multimin_fminimizer_type *gsl_multimin_fminimizer_nmsimplex2rand;
0223
0224 __END_DECLS
0225
0226 #endif