00001
00008 #include <stdio.h>
00009 #include <time.h>
00010 #include <gsl/gsl_blas.h>
00011 #include <gsl/gsl_deriv.h>
00012 #include <gsl/gsl_errno.h>
00013 #include <gsl/gsl_linalg.h>
00014 #include <gsl/gsl_math.h>
00015 #include <gsl/gsl_matrix.h>
00016 #include <gsl/gsl_multifit_nlin.h>
00017 #include <gsl/gsl_multimin.h>
00018 #include <gsl/gsl_randist.h>
00019 #include <gsl/gsl_rng.h>
00020 #include <gsl/gsl_siman.h>
00021 #include <gsl/gsl_vector.h>
00022
00023
00024
00025
00026
00027 enum _method {
00028 MULTIDIM_MINIMIZAION = 0,
00029 LEAST_SQUARES = 1,
00030 SIMULATED_ANNEAL = 2,
00031 };
00032 enum _algorithm {
00033 CONJUGATE_FR = 0,
00034 CONJUGATE_PR = 1,
00035 VECTOR_BFGS = 2,
00036 STEEPEST_DESCENT = 3,
00037 LMSDR = 4,
00038 LMDR = 5,
00039 };
00040 enum _type {
00041 NOTHING = 0,
00042 SINGLE = 1,
00043 MULTIPLE = 2,
00044 };
00045
00046
00052 struct _data {
00053 unsigned int idx;
00054 unsigned int num;
00055 double *xpt;
00056 double *ypt;
00057 double *sig;
00058 };
00059
00060
00065 struct _meth {
00066 const gsl_multimin_fdfminimizer_type *mdmin_type; \
00068 const gsl_multifit_fdfsolver_type *lstsq_type;
00070 int (*solver_function_call)(void *);
00072 };
00073
00074
00080 struct _func {
00081 unsigned int idx; \
00083 unsigned int num; \
00085 double (*fun_s)(double, double *); \
00087 double (*fun_m)(double, double *, unsigned int, unsigned int); \
00089 };
00090
00091
00097 struct _para {
00098 unsigned int idx;
00099 unsigned int num;
00100 double *par;
00101 unsigned int *fix;
00102 double *err;
00103 double *tmp;
00104 };
00105
00106
00126 typedef struct _fitf{
00127 struct _data *data;
00128 struct _func *func;
00129 unsigned int numd;
00130 unsigned int numf;
00131 struct _para para;
00132 struct _meth meth;
00133 double chis;
00134 unsigned int indx;
00135 unsigned int maxi;
00136 unsigned int type;
00137 } gsl_fit;
00138
00143
00157 gsl_fit *gsl_fit_alloc(unsigned int method, \
00158 unsigned int algotrithm, \
00159 unsigned int max_iter, \
00160 unsigned int type);
00161
00174 int gsl_fit_add_dataset(gsl_fit *gft, \
00175 double *xpt, double *ypt, double *sig, \
00176 unsigned int num);
00177
00189 int gsl_fit_add_model_single(gsl_fit *gft, \
00190 double (*fun)(double, double *), \
00191 double *par, unsigned int *fix, unsigned int num);
00192
00204 int gsl_fit_add_model_multiple(gsl_fit *gft, \
00205 double (*fun)(double, double *, \
00206 unsigned int, unsigned int), \
00207 double *par, unsigned int *fix, unsigned int num);
00208
00217 int gsl_fit_print(gsl_fit *gft, char *filename);
00218
00226 int gsl_fit_do(gsl_fit *gft);
00227
00239 int gsl_fit_get_result(gsl_fit *gft, \
00240 unsigned int function, \
00241 double *par, double *err, unsigned int num);
00242
00250 void gsl_fit_free(gsl_fit *gft);
00251
00252
00253
00254
00289 int lstsq_main(void *void_gsl_fit_struct);
00290
00302 int lstsq_fdf(const gsl_vector *parameters, \
00303 void *void_gsl_fit_struct, \
00304 gsl_vector *r, gsl_matrix *J);
00305
00320 int lstsq_f(const gsl_vector *parameters, \
00321 void *void_gsl_fit_struct, \
00322 gsl_vector *r);
00323
00346 int lstsq_df(const gsl_vector *parameters, \
00347 void *void_gsl_fit_struct, \
00348 gsl_matrix *J);
00349
00363 double lstsq_derivatives(double p, void *void_gsl_fit_struct);
00364
00376 int lstsq_status(gsl_fit *gft, \
00377 int iteration, \
00378 gsl_multifit_fdfsolver *s, \
00379 gsl_matrix *c);
00380
00381
00382
00383
00409 int mdmin_main(void *void_gsl_fit_struct);
00410
00422 void mdmin_fdf(const gsl_vector *parameters, \
00423 void *void_gsl_fit_struct, \
00424 double *f, gsl_vector *df);
00425
00437 double mdmin_f(const gsl_vector *parameters, \
00438 void *void_gsl_fit_struct);
00439
00456 void mdmin_df(const gsl_vector *parameters, \
00457 void *void_gsl_fit_struct, \
00458 gsl_vector *df);
00459
00471 double mdmin_derivatives(double p, void *void_gsl_fit_struct);
00472
00481 double mdmin_chisq(gsl_fit *gft);
00482
00493 void mdmin_status(gsl_fit *gft, \
00494 int iteration, \
00495 gsl_multimin_fdfminimizer *s);
00496
00497
00498
00499
00522 int siman_main(void *void_gsl_fit_struct);
00523
00532 double siman_energy(void *void_gsl_fit_struct);
00533
00544 void siman_step(const gsl_rng *rng, void *void_gsl_fit_struct, \
00545 double step_size);
00546
00556 void siman_copy(void *source, void *dest);
00557
00566 void *siman_construct(void *void_gsl_fit_struct);
00567
00576 void siman_destroy(void *void_gsl_fit_struct);
00577
00586 double siman_metric(void *void_gsl_fit_struct1, void *void_gsl_fit_struct2);
00587
00603 void siman_print(void *void_gsl_fit_struct);