00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023 #ifndef O2SCL_GSL_FIT_H
00024 #define O2SCL_GSL_FIT_H
00025
00026 #include <o2scl/fit_base.h>
00027 #include <gsl/gsl_rng.h>
00028 #include <gsl/gsl_vector.h>
00029 #include <gsl/gsl_blas.h>
00030 #include <gsl/gsl_multifit_nlin.h>
00031 #include <o2scl/omatrix_tlate.h>
00032
00033 #ifndef DOXYGENP
00034 namespace o2scl {
00035 #endif
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065 template<class param_t, class func_t, class vec_t=ovector_view,
00066 class mat_t=omatrix_view, class bool_vec_t=bool *> class gsl_fit :
00067 public fit_base<param_t,func_t,vec_t,mat_t> {
00068
00069 public:
00070
00071 gsl_fit() {
00072 gsl_set_error_handler(err_hnd->gsl_hnd);
00073 max_iter=500;
00074 epsabs=1.0e-4;
00075 epsrel=1.0e-4;
00076 use_scaled=true;
00077 test_gradient=false;
00078 }
00079
00080 virtual ~gsl_fit() {}
00081
00082
00083
00084
00085
00086
00087
00088 virtual int fit(size_t ndat, vec_t &xdat, vec_t &ydat, vec_t &yerr,
00089 size_t npar, vec_t &par, mat_t &covar, double &chi2,
00090 param_t &pa, func_t &fitfun) {
00091
00092 int status, iter=0;
00093 size_t i, j;
00094
00095 gsl_vector *x=gsl_vector_alloc(npar);
00096 for(i=0;i<npar;i++) gsl_vector_set(x,i,par[i]);
00097
00098
00099 gsl_multifit_function_fdf f;
00100 gsl_vector *grad=gsl_vector_alloc(npar);
00101 func_par fpar={fitfun,&pa,ndat,&xdat,&ydat,&yerr,npar};
00102 f.f=func;
00103 f.df=dfunc;
00104 f.fdf=fdfunc;
00105 f.n=ndat;
00106 f.p=npar;
00107 f.params=&fpar;
00108
00109
00110 gsl_multifit_fdfsolver *s;
00111 if (use_scaled) {
00112 s=gsl_multifit_fdfsolver_alloc
00113 (gsl_multifit_fdfsolver_lmsder,ndat,npar);
00114 } else {
00115 s=gsl_multifit_fdfsolver_alloc
00116 (gsl_multifit_fdfsolver_lmder,ndat,npar);
00117 }
00118 gsl_multifit_fdfsolver_set(s,&f,x);
00119
00120
00121 do {
00122 iter++;
00123 status=gsl_multifit_fdfsolver_iterate(s);
00124
00125 if (status) {
00126 break;
00127 }
00128
00129 status=gsl_multifit_test_delta(s->dx,s->x,epsabs,epsrel);
00130 if (test_gradient && status==gsl_continue) {
00131 gsl_multifit_gradient(s->J,s->x,grad);
00132 status=gsl_multifit_test_gradient(grad,epsabs);
00133 }
00134
00135 print_iter(npar,s->dx,s->x,iter,epsabs,epsrel);
00136
00137 } while (status == gsl_continue && iter < max_iter);
00138
00139
00140
00141 gsl_matrix *gcovar=gsl_matrix_alloc(npar,npar);
00142 gsl_multifit_covar(s->J,0.0,gcovar);
00143 for(i=0;i<npar;i++) {
00144 par[i]=gsl_vector_get(s->x,i);
00145 for(j=0;j<npar;j++) {
00146 covar.set(i,j,gsl_matrix_get(gcovar,i,j));
00147 }
00148 }
00149 gsl_matrix_free(gcovar);
00150
00151
00152 chi2=gsl_blas_dnrm2(s->f);
00153 chi2*=chi2;
00154
00155
00156 gsl_multifit_fdfsolver_free(s);
00157 gsl_vector_free(grad);
00158
00159 return 0;
00160 }
00161
00162
00163 int max_iter;
00164
00165
00166 double epsabs;
00167
00168
00169 double epsrel;
00170
00171
00172 bool test_gradient;
00173
00174
00175
00176 bool use_scaled;
00177
00178
00179 virtual const char *type() { return "gsl_fit"; }
00180
00181 #ifndef DOXYGEN_INTERNAL
00182
00183 protected:
00184
00185
00186 virtual int print_iter(int nv, gsl_vector *x, gsl_vector *dx, int iter,
00187 double l_epsabs, double l_epsrel) {
00188 if (this->verbose<=0) return 0;
00189
00190 int i;
00191 char ch;
00192 double dxmax, epsmin;
00193
00194 std::cout << "Iteration: " << iter << std::endl;
00195 text_out_file outs(&std::cout,79);
00196 outs.word_out("x:");
00197 dxmax=gsl_vector_get(dx,0);
00198 epsmin=l_epsabs+l_epsrel*gsl_vector_get(x,0);
00199 for(i=0;i<nv;i++) {
00200 if (dxmax<gsl_vector_get(dx,i)) dxmax=gsl_vector_get(dx,i);
00201 if (epsmin>l_epsabs+l_epsrel*gsl_vector_get(x,i)) {
00202 epsmin=l_epsabs+l_epsrel*gsl_vector_get(x,i);
00203 }
00204 outs.double_out(gsl_vector_get(x,i));
00205 }
00206 outs.end_line();
00207 std::cout << " Val: " << dxmax << " Lim: " << epsmin << std::endl;
00208 if (this->verbose>1) {
00209 std::cout << "Press a key and type enter to continue. ";
00210 std::cin >> ch;
00211 }
00212
00213 return 0;
00214 }
00215
00216
00217
00218
00219 typedef struct {
00220
00221 func_t &f;
00222
00223 param_t *vp;
00224
00225 int ndat;
00226
00227 vec_t *xdat;
00228
00229 vec_t *ydat;
00230
00231 vec_t *yerr;
00232
00233 int npar;
00234 } func_par;
00235
00236
00237 static int func(const gsl_vector *x, void *pa, gsl_vector *f) {
00238 func_par *fp=(func_par *)pa;
00239 int i, j;
00240 double yi;
00241 ovector xp(fp->npar);
00242 for(j=0;j<fp->npar;j++) xp[j]=gsl_vector_get(x,j);
00243 for(i=0;i<fp->ndat;i++) {
00244 fp->f(fp->npar,xp,(*fp->xdat)[i],yi,*(fp->vp));
00245 gsl_vector_set(f,i,(yi-(*fp->ydat)[i])/(*fp->yerr)[i]);
00246 }
00247
00248 return 0;
00249 }
00250
00251
00252 static int dfunc(const gsl_vector *x, void *pa, gsl_matrix *jac) {
00253 func_par *fp=(func_par *)pa;
00254 int i, j;
00255 double yi, ylo, yhi, xtemp, eps=1.0e-4;
00256 ovector xp(fp->npar);
00257
00258 for(j=0;j<fp->npar;j++) xp[j]=gsl_vector_get(x,j);
00259 for(j=0;j<fp->npar;j++) {
00260 for(i=0;i<fp->ndat;i++) {
00261 xtemp=xp[j];
00262 xp[j]+=eps;
00263 fp->f(fp->npar,xp,(*fp->xdat)[i],yi,*(fp->vp));
00264 yhi=(yi-(*fp->ydat)[i])/(*fp->yerr)[i];
00265 xp[j]=xtemp;
00266 fp->f(fp->npar,xp,(*fp->xdat)[i],yi,*(fp->vp));
00267 ylo=(yi-(*fp->ydat)[i])/(*fp->yerr)[i];
00268 gsl_matrix_set(jac,i,j,(yhi-ylo)/eps);
00269 }
00270 }
00271
00272 return 0;
00273 }
00274
00275
00276 static int fdfunc(const gsl_vector *x, void *pa, gsl_vector *f,
00277 gsl_matrix *jac) {
00278 func(x,pa,f);
00279 dfunc(x,pa,jac);
00280 return 0;
00281 }
00282
00283 #endif
00284
00285 };
00286
00287 #ifndef DOXYGENP
00288 }
00289 #endif
00290
00291 #endif