00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, 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, class vec_t=ovector_view, 00053 class mat_t=omatrix_view> 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 00067 This structure is given so that the user can specify the 00068 minimizer to use. 00069 */ 00070 typedef struct { 00071 /// The fitting function 00072 func_t &f; 00073 /// The user-specified parameter 00074 param_t &vp; 00075 /// The number of data 00076 int ndat; 00077 /// The x values 00078 vec_t *xdat; 00079 /// The y values 00080 vec_t *ydat; 00081 /// The y uncertainties 00082 vec_t *yerr; 00083 /// The number of fitting parameters 00084 int npar; 00085 } func_par; 00086 00087 /** \brief Fit the data specified in (xdat,ydat) to 00088 the function \c fitfun with the parameters in \c par. 00089 00090 The covariance matrix for the parameters is returned in \c 00091 covar and the value of \f$ \chi^2 \f$ is returned in \c chi2. 00092 */ 00093 virtual int fit(size_t ndat, vec_t &xdat, vec_t &ydat, vec_t &yerr, 00094 size_t npar, vec_t &par, mat_t &covar, double &chi2, 00095 param_t &pa, func_t &fitfun) { 00096 00097 int status, i, iter=0; 00098 double dtemp; 00099 00100 // --------------------------------------------------- 00101 // First minimize with the specified multi_min object 00102 // --------------------------------------------------- 00103 00104 func_par fpar={fitfun,pa,ndat,&xdat,&ydat,&yerr,npar}; 00105 00106 multi_funct_mfptr_noerr<min_fit,func_par *> mmf(this,&min_fit::min_func); 00107 func_par *fp=&fpar; 00108 mmp->mmin(npar,par,dtemp,fp,mmf); 00109 00110 gsl_multifit_fdfsolver *s=gsl_multifit_fdfsolver_alloc 00111 (gsl_multifit_fdfsolver_lmsder,ndat,npar); 00112 00113 // ------------------------------------------------------------- 00114 // Now use the GSL routines to calculate the covariance matrix, 00115 // and chi-squared. 00116 // ------------------------------------------------------------- 00117 00118 gsl_vector *x=(gsl_vector *)(&par); 00119 gsl_multifit_function_fdf f; 00120 00121 f.f=func; 00122 f.df=dfunc; 00123 f.fdf=fdfunc; 00124 f.n=ndat; 00125 f.p=npar; 00126 f.params=&fpar; 00127 00128 s=gsl_multifit_fdfsolver_alloc 00129 (gsl_multifit_fdfsolver_lmsder,ndat,npar); 00130 gsl_multifit_fdfsolver_set(s,&f,x); 00131 00132 gsl_matrix *gcovar=gsl_matrix_alloc(npar,npar); 00133 gsl_multifit_covar(s->J,0.0,gcovar); 00134 for(i=0;i<((int)npar);i++) { 00135 par[i]=gsl_vector_get(s->x,i); 00136 for(int j=0;j<((int)npar);j++) { 00137 covar[i][j]=gsl_matrix_get(gcovar,i,j); 00138 } 00139 } 00140 gsl_matrix_free(gcovar); 00141 00142 chi2=gsl_blas_dnrm2(s->f); 00143 chi2*=chi2; 00144 00145 gsl_multifit_fdfsolver_free(s); 00146 00147 return 0; 00148 } 00149 00150 /** \brief Set the multi_min object to use (default is gsl_mmin_nmsimplex) 00151 */ 00152 int set_multi_min(multi_min<func_par *, 00153 multi_funct<func_par *,vec_t> > &mm) { 00154 mmp=&mm; 00155 min_set=true; 00156 return 0; 00157 } 00158 00159 /// The default minimizer 00160 gsl_mmin_simp<func_par *,multi_funct<func_par *,vec_t> > def_multi_min; 00161 00162 /// Return string denoting type ("min_fit") 00163 virtual const char *type() { return "min_fit"; } 00164 00165 #ifndef DOXYGEN_INTERNAL 00166 00167 protected: 00168 00169 /// The minimizer 00170 multi_min<func_par *,multi_funct<func_par *,vec_t> > *mmp; 00171 00172 /// True if the minimizer has been set by the user 00173 bool min_set; 00174 00175 /// Evaluate the function 00176 static int func(const gsl_vector *x, void *pa, gsl_vector *f) { 00177 func_par *fp=(func_par *)pa; 00178 int i, j; 00179 double yi; 00180 ovector xp(fp->npar); 00181 for(j=0;j<fp->npar;j++) xp[j]=gsl_vector_get(x,j); 00182 for(i=0;i<fp->ndat;i++) { 00183 fp->f(fp->npar,xp,(*fp->xdat)[i],yi,fp->vp); 00184 gsl_vector_set(f,i,(yi-(*fp->ydat)[i])/(*fp->yerr)[i]); 00185 } 00186 00187 return 0; 00188 } 00189 00190 /// Evaluate the jacobian 00191 static int dfunc(const gsl_vector *x, void *pa, gsl_matrix *jac) { 00192 func_par *fp=(func_par *)pa; 00193 int i, j; 00194 double yi, ylo, yhi, xtemp, eps=1.0e-4; 00195 ovector xp(fp->npar); 00196 00197 for(j=0;j<fp->npar;j++) xp[j]=gsl_vector_get(x,j); 00198 for(j=0;j<fp->npar;j++) { 00199 for(i=0;i<fp->ndat;i++) { 00200 xtemp=xp[j]; 00201 xp[j]+=eps; 00202 fp->f(fp->npar,xp,(*fp->xdat)[i],yi,fp->vp); 00203 yhi=(yi-(*fp->ydat)[i])/(*fp->yerr)[i]; 00204 xp[j]=xtemp; 00205 fp->f(fp->npar,xp,(*fp->xdat)[i],yi,fp->vp); 00206 ylo=(yi-(*fp->ydat)[i])/(*fp->yerr)[i]; 00207 gsl_matrix_set(jac,i,j,(yhi-ylo)/eps); 00208 } 00209 } 00210 00211 return 0; 00212 } 00213 00214 /// Evaluate the function and the jacobian 00215 static int fdfunc(const gsl_vector *x, void *pa, gsl_vector *f, 00216 gsl_matrix *jac) { 00217 func(x,pa,f); 00218 dfunc(x,pa,jac); 00219 return 0; 00220 } 00221 00222 /// The function to minimize 00223 double min_func(size_t np, const vec_t &xp, func_par *&fp) { 00224 int i; 00225 double ret=0.0; 00226 double yi; 00227 00228 ovector xpnoc(np); 00229 for(size_t ii=0;ii<np;ii++) xpnoc[ii]=xp[ii]; 00230 00231 for(i=0;i<fp->ndat;i++) { 00232 fp->f(np,xpnoc,(*fp->xdat)[i],yi,fp->vp); 00233 ret+=pow((yi-(*fp->ydat)[i])/(*fp->yerr)[i],2.0); 00234 } 00235 00236 return ret; 00237 } 00238 00239 #endif 00240 00241 }; 00242 00243 #ifndef DOXYGENP 00244 } 00245 #endif 00246 00247 #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