00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, 2008, 2009, Andrew W. Steiner 00005 00006 This file is part of O2scl. 00007 00008 O2scl is free software; you can redistribute it and/or modify 00009 it under the terms of the GNU General Public License as published by 00010 the Free Software Foundation; either version 3 of the License, or 00011 (at your option) any later version. 00012 00013 O2scl is distributed in the hope that it will be useful, 00014 but WITHOUT ANY WARRANTY; without even the implied warranty of 00015 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00016 GNU General Public License for more details. 00017 00018 You should have received a copy of the GNU General Public License 00019 along with O2scl. If not, see <http://www.gnu.org/licenses/>. 00020 00021 ------------------------------------------------------------------- 00022 */ 00023 #ifndef O2SCL_MIN_FIT_H 00024 #define O2SCL_MIN_FIT_H 00025 00026 #include <gsl/gsl_rng.h> 00027 #include <gsl/gsl_vector.h> 00028 #include <gsl/gsl_blas.h> 00029 #include <gsl/gsl_multifit_nlin.h> 00030 #include <o2scl/multi_min.h> 00031 #include <o2scl/multi_funct.h> 00032 #include <o2scl/gsl_mmin_simp.h> 00033 #include <o2scl/fit_base.h> 00034 00035 #ifndef DOXYGENP 00036 namespace o2scl { 00037 #endif 00038 00039 /** 00040 \brief Non-linear least-squares fitting class with generic minimizer 00041 00042 This minimizes a generic fitting function using any multi_min 00043 object, and then uses the GSL routines to calculate the 00044 uncertainties in the parameters and the covariance matrix. 00045 00046 This can be useful for fitting problems which might be better 00047 handled by more complex minimizers than those that are used in 00048 gsl_fit. For problems with many local minima near the global 00049 minimum, using a sim_anneal object with this class can sometimes 00050 produce better results than gsl_fit. 00051 */ 00052 template<class param_t, class func_t=fit_funct<param_t>, 00053 class vec_t=ovector_base, class mat_t=omatrix_base> class min_fit : 00054 public fit_base<param_t,func_t,vec_t,mat_t> { 00055 public: 00056 00057 min_fit() { 00058 min_set=false; 00059 mmp=&def_multi_min; 00060 } 00061 00062 virtual ~min_fit() {} 00063 00064 /** 00065 \brief A structure for passing information to the GSL functions 00066 for the \ref min_fit class 00067 00068 This structure is given so that the user can specify the 00069 minimizer to use. 00070 */ 00071 typedef struct { 00072 /// The fitting function 00073 func_t &f; 00074 /// The user-specified parameter 00075 param_t &vp; 00076 /// The number of data 00077 int ndat; 00078 /// The x values 00079 vec_t *xdat; 00080 /// The y values 00081 vec_t *ydat; 00082 /// The y uncertainties 00083 vec_t *yerr; 00084 /// The number of fitting parameters 00085 int npar; 00086 } func_par; 00087 00088 /** \brief Fit the data specified in (xdat,ydat) to 00089 the function \c fitfun with the parameters in \c par. 00090 00091 The covariance matrix for the parameters is returned in \c 00092 covar and the value of \f$ \chi^2 \f$ is returned in \c chi2. 00093 */ 00094 virtual int fit(size_t ndat, vec_t &xdat, vec_t &ydat, vec_t &yerr, 00095 size_t npar, vec_t &par, mat_t &covar, double &chi2, 00096 param_t &pa, func_t &fitfun) { 00097 00098 int status, i, iter=0; 00099 double dtemp; 00100 00101 // --------------------------------------------------- 00102 // First minimize with the specified multi_min object 00103 // --------------------------------------------------- 00104 00105 func_par fpar={fitfun,pa,ndat,&xdat,&ydat,&yerr,npar}; 00106 00107 multi_funct_mfptr_noerr<min_fit,func_par *> mmf(this,&min_fit::min_func); 00108 func_par *fp=&fpar; 00109 mmp->mmin(npar,par,dtemp,fp,mmf); 00110 00111 gsl_multifit_fdfsolver *s=gsl_multifit_fdfsolver_alloc 00112 (gsl_multifit_fdfsolver_lmsder,ndat,npar); 00113 00114 // ------------------------------------------------------------- 00115 // Now use the GSL routines to calculate the covariance matrix, 00116 // and chi-squared. 00117 // ------------------------------------------------------------- 00118 00119 gsl_vector *x=(gsl_vector *)(&par); 00120 gsl_multifit_function_fdf f; 00121 00122 f.f=func; 00123 f.df=dfunc; 00124 f.fdf=fdfunc; 00125 f.n=ndat; 00126 f.p=npar; 00127 f.params=&fpar; 00128 00129 s=gsl_multifit_fdfsolver_alloc 00130 (gsl_multifit_fdfsolver_lmsder,ndat,npar); 00131 gsl_multifit_fdfsolver_set(s,&f,x); 00132 00133 gsl_matrix *gcovar=gsl_matrix_alloc(npar,npar); 00134 gsl_multifit_covar(s->J,0.0,gcovar); 00135 for(i=0;i<((int)npar);i++) { 00136 par[i]=gsl_vector_get(s->x,i); 00137 for(int j=0;j<((int)npar);j++) { 00138 covar[i][j]=gsl_matrix_get(gcovar,i,j); 00139 } 00140 } 00141 gsl_matrix_free(gcovar); 00142 00143 chi2=gsl_blas_dnrm2(s->f); 00144 chi2*=chi2; 00145 00146 gsl_multifit_fdfsolver_free(s); 00147 00148 return 0; 00149 } 00150 00151 /** \brief Set the multi_min object to use (default is gsl_mmin_nmsimplex) 00152 */ 00153 int set_multi_min(multi_min<func_par *, 00154 multi_funct<func_par *,vec_t> > &mm) { 00155 mmp=&mm; 00156 min_set=true; 00157 return 0; 00158 } 00159 00160 /// The default minimizer 00161 gsl_mmin_simp<func_par *,multi_funct<func_par *,vec_t> > def_multi_min; 00162 00163 /// Return string denoting type ("min_fit") 00164 virtual const char *type() { return "min_fit"; } 00165 00166 #ifndef DOXYGEN_INTERNAL 00167 00168 protected: 00169 00170 /// The minimizer 00171 multi_min<func_par *,multi_funct<func_par *,vec_t> > *mmp; 00172 00173 /// True if the minimizer has been set by the user 00174 bool min_set; 00175 00176 /// Evaluate the function 00177 static int func(const gsl_vector *x, void *pa, gsl_vector *f) { 00178 func_par *fp=(func_par *)pa; 00179 int i, j; 00180 double yi; 00181 ovector xp(fp->npar); 00182 for(j=0;j<fp->npar;j++) xp[j]=gsl_vector_get(x,j); 00183 for(i=0;i<fp->ndat;i++) { 00184 fp->f(fp->npar,xp,(*fp->xdat)[i],yi,fp->vp); 00185 gsl_vector_set(f,i,(yi-(*fp->ydat)[i])/(*fp->yerr)[i]); 00186 } 00187 00188 return 0; 00189 } 00190 00191 /// Evaluate the jacobian 00192 static int dfunc(const gsl_vector *x, void *pa, gsl_matrix *jac) { 00193 func_par *fp=(func_par *)pa; 00194 int i, j; 00195 double yi, ylo, yhi, xtemp, eps=1.0e-4; 00196 ovector xp(fp->npar); 00197 00198 for(j=0;j<fp->npar;j++) xp[j]=gsl_vector_get(x,j); 00199 for(j=0;j<fp->npar;j++) { 00200 for(i=0;i<fp->ndat;i++) { 00201 xtemp=xp[j]; 00202 xp[j]+=eps; 00203 fp->f(fp->npar,xp,(*fp->xdat)[i],yi,fp->vp); 00204 yhi=(yi-(*fp->ydat)[i])/(*fp->yerr)[i]; 00205 xp[j]=xtemp; 00206 fp->f(fp->npar,xp,(*fp->xdat)[i],yi,fp->vp); 00207 ylo=(yi-(*fp->ydat)[i])/(*fp->yerr)[i]; 00208 gsl_matrix_set(jac,i,j,(yhi-ylo)/eps); 00209 } 00210 } 00211 00212 return 0; 00213 } 00214 00215 /// Evaluate the function and the jacobian 00216 static int fdfunc(const gsl_vector *x, void *pa, gsl_vector *f, 00217 gsl_matrix *jac) { 00218 func(x,pa,f); 00219 dfunc(x,pa,jac); 00220 return 0; 00221 } 00222 00223 /// The function to minimize 00224 double min_func(size_t np, const vec_t &xp, func_par *&fp) { 00225 int i; 00226 double ret=0.0; 00227 double yi; 00228 00229 ovector xpnoc(np); 00230 for(size_t ii=0;ii<np;ii++) xpnoc[ii]=xp[ii]; 00231 00232 for(i=0;i<fp->ndat;i++) { 00233 fp->f(np,xpnoc,(*fp->xdat)[i],yi,fp->vp); 00234 ret+=pow((yi-(*fp->ydat)[i])/(*fp->yerr)[i],2.0); 00235 } 00236 00237 return ret; 00238 } 00239 00240 #endif 00241 00242 }; 00243 00244 #ifndef DOXYGENP 00245 } 00246 #endif 00247 00248 #endif
Documentation generated with Doxygen and provided under the GNU Free Documentation License. See License Information for details.
Project hosting provided by
,
O2scl Sourceforge Project Page