![]() |
Object-oriented Scientific Computing Library: Version 0.910
|
00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006-2012, 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 /** \brief Non-linear least-squares fitting class with generic minimizer 00040 00041 This minimizes a generic fitting function using any multi_min 00042 object, and then uses the GSL routines to calculate the 00043 uncertainties in the parameters and the covariance matrix. 00044 00045 This can be useful for fitting problems which might be better 00046 handled by more complex minimizers than those that are used in 00047 gsl_fit. For problems with many local minima near the global 00048 minimum, using a sim_anneal object with this class can sometimes 00049 produce better results than gsl_fit. 00050 00051 Default template arguments 00052 - \c func_t - \ref fit_funct<> 00053 - \c vec_t - \ref ovector_base 00054 - \c alloc_vec_t - \ref ovector 00055 */ 00056 template<class func_t=fit_funct<>, class vec_t=ovector_base, 00057 class mat_t=omatrix_base> class min_fit : 00058 public fit_base<func_t,vec_t,mat_t> { 00059 00060 public: 00061 00062 min_fit() { 00063 min_set=false; 00064 mmp=&def_multi_min; 00065 } 00066 00067 virtual ~min_fit() {} 00068 00069 /** \brief A structure for passing information to the GSL functions 00070 for the \ref min_fit class 00071 00072 This structure is given so that the user can specify the 00073 minimizer to use. 00074 */ 00075 typedef struct { 00076 /// The fitting function 00077 func_t &f; 00078 /// The number of data 00079 int ndat; 00080 /// The x values 00081 vec_t *xdat; 00082 /// The y values 00083 vec_t *ydat; 00084 /// The y uncertainties 00085 vec_t *yerr; 00086 /// The number of fitting parameters 00087 int npar; 00088 } func_par; 00089 00090 /** \brief Fit the data specified in (xdat,ydat) to 00091 the function \c fitfun with the parameters in \c par. 00092 00093 The covariance matrix for the parameters is returned in \c 00094 covar and the value of \f$ \chi^2 \f$ is returned in \c chi2. 00095 */ 00096 virtual int fit(size_t ndat, vec_t &xdat, vec_t &ydat, vec_t &yerr, 00097 size_t npar, vec_t &par, mat_t &covar, double &chi2, 00098 func_t &fitfun) { 00099 00100 int status, i, iter=0; 00101 double dtemp; 00102 00103 // --------------------------------------------------- 00104 // First minimize with the specified multi_min object 00105 // --------------------------------------------------- 00106 00107 func_par fpar={fitfun,ndat,&xdat,&ydat,&yerr,npar}; 00108 00109 func_par *fp=&fpar; 00110 multi_funct_mfptr_param<min_fit,func_par *> 00111 mmf(this,&min_fit::min_func,fp); 00112 mmp->mmin(npar,par,dtemp,mmf); 00113 00114 gsl_multifit_fdfsolver *s=gsl_multifit_fdfsolver_alloc 00115 (gsl_multifit_fdfsolver_lmsder,ndat,npar); 00116 00117 // ------------------------------------------------------------- 00118 // Now use the GSL routines to calculate the covariance matrix, 00119 // and chi-squared. 00120 // ------------------------------------------------------------- 00121 00122 gsl_vector *x=(gsl_vector *)(&par); 00123 gsl_multifit_function_fdf f; 00124 00125 f.f=func; 00126 f.df=dfunc; 00127 f.fdf=fdfunc; 00128 f.n=ndat; 00129 f.p=npar; 00130 f.params=&fpar; 00131 00132 s=gsl_multifit_fdfsolver_alloc 00133 (gsl_multifit_fdfsolver_lmsder,ndat,npar); 00134 gsl_multifit_fdfsolver_set(s,&f,x); 00135 00136 gsl_matrix *gcovar=gsl_matrix_alloc(npar,npar); 00137 gsl_multifit_covar(s->J,0.0,gcovar); 00138 for(i=0;i<((int)npar);i++) { 00139 par[i]=gsl_vector_get(s->x,i); 00140 for(int j=0;j<((int)npar);j++) { 00141 covar[i][j]=gsl_matrix_get(gcovar,i,j); 00142 } 00143 } 00144 gsl_matrix_free(gcovar); 00145 00146 chi2=gsl_blas_dnrm2(s->f); 00147 chi2*=chi2; 00148 00149 gsl_multifit_fdfsolver_free(s); 00150 00151 return 0; 00152 } 00153 00154 /** \brief Set the multi_min object to use (default is gsl_mmin_nmsimplex) 00155 */ 00156 int set_multi_min(multi_min<multi_funct<vec_t> > &mm) { 00157 mmp=&mm; 00158 min_set=true; 00159 return 0; 00160 } 00161 00162 /// The default minimizer 00163 gsl_mmin_simp<multi_funct<vec_t> > def_multi_min; 00164 00165 /// Return string denoting type ("min_fit") 00166 virtual const char *type() { return "min_fit"; } 00167 00168 #ifndef DOXYGEN_INTERNAL 00169 00170 protected: 00171 00172 /// The minimizer 00173 multi_min<multi_funct<vec_t> > *mmp; 00174 00175 /// True if the minimizer has been set by the user 00176 bool min_set; 00177 00178 /// Evaluate the function 00179 static int func(const gsl_vector *x, void *pa, gsl_vector *f) { 00180 func_par *fp=(func_par *)pa; 00181 int i, j; 00182 double yi; 00183 ovector xp(fp->npar); 00184 for(j=0;j<fp->npar;j++) xp[j]=gsl_vector_get(x,j); 00185 for(i=0;i<fp->ndat;i++) { 00186 yi=fp->f(fp->npar,xp,(*fp->xdat)[i]); 00187 gsl_vector_set(f,i,(yi-(*fp->ydat)[i])/(*fp->yerr)[i]); 00188 } 00189 00190 return 0; 00191 } 00192 00193 /// Evaluate the jacobian 00194 static int dfunc(const gsl_vector *x, void *pa, gsl_matrix *jac) { 00195 func_par *fp=(func_par *)pa; 00196 int i, j; 00197 double yi, ylo, yhi, xtemp, eps=1.0e-4; 00198 ovector xp(fp->npar); 00199 00200 for(j=0;j<fp->npar;j++) xp[j]=gsl_vector_get(x,j); 00201 for(j=0;j<fp->npar;j++) { 00202 for(i=0;i<fp->ndat;i++) { 00203 xtemp=xp[j]; 00204 xp[j]+=eps; 00205 yi=fp->f(fp->npar,xp,(*fp->xdat)[i]); 00206 yhi=(yi-(*fp->ydat)[i])/(*fp->yerr)[i]; 00207 xp[j]=xtemp; 00208 yi=fp->f(fp->npar,xp,(*fp->xdat)[i]); 00209 ylo=(yi-(*fp->ydat)[i])/(*fp->yerr)[i]; 00210 gsl_matrix_set(jac,i,j,(yhi-ylo)/eps); 00211 } 00212 } 00213 00214 return 0; 00215 } 00216 00217 /// Evaluate the function and the jacobian 00218 static int fdfunc(const gsl_vector *x, void *pa, gsl_vector *f, 00219 gsl_matrix *jac) { 00220 func(x,pa,f); 00221 dfunc(x,pa,jac); 00222 return 0; 00223 } 00224 00225 /// The function to minimize 00226 double min_func(size_t np, const vec_t &xp, func_par *&fp) { 00227 int i; 00228 double ret=0.0; 00229 double yi; 00230 00231 ovector xpnoc(np); 00232 for(size_t ii=0;ii<np;ii++) xpnoc[ii]=xp[ii]; 00233 00234 for(i=0;i<fp->ndat;i++) { 00235 yi=fp->f(np,xpnoc,(*fp->xdat)[i]); 00236 ret+=pow((yi-(*fp->ydat)[i])/(*fp->yerr)[i],2.0); 00237 } 00238 00239 return ret; 00240 } 00241 00242 #endif 00243 00244 }; 00245 00246 #ifndef DOXYGENP 00247 } 00248 #endif 00249 00250 #endif
Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).