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 /* multifit/fsolver.c 00024 * 00025 * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2007 Brian Gough 00026 * 00027 * This program is free software; you can redistribute it and/or modify 00028 * it under the terms of the GNU General Public License as published by 00029 * the Free Software Foundation; either version 3 of the License, or (at 00030 * your option) any later version. 00031 * 00032 * This program is distributed in the hope that it will be useful, but 00033 * WITHOUT ANY WARRANTY; without even the implied warranty of 00034 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00035 * General Public License for more details. 00036 * 00037 * You should have received a copy of the GNU General Public License 00038 * along with this program; if not, write to the Free Software 00039 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 00040 * 02110-1301, USA. 00041 */ 00042 00043 #ifndef O2SCL_GSL_FIT_H 00044 #define O2SCL_GSL_FIT_H 00045 00046 #include <o2scl/fit_base.h> 00047 #include <gsl/gsl_rng.h> 00048 #include <gsl/gsl_vector.h> 00049 #include <gsl/gsl_blas.h> 00050 #include <gsl/gsl_multifit_nlin.h> 00051 #include <o2scl/omatrix_tlate.h> 00052 00053 #ifndef DOXYGENP 00054 namespace o2scl { 00055 #endif 00056 00057 /** 00058 \brief Non-linear least-squares fitting class (GSL) 00059 00060 The GSL-based fitting class using a Levenberg-Marquardt type 00061 algorithm. The algorithm stops when 00062 \f[ 00063 |dx_i| < \mathrm{epsabs}+\mathrm{epsrel}\times|x_i| 00064 \f] 00065 where \f$dx\f$ is the last step and \f$x\f$ is the current 00066 position. If test_gradient is true, then additionally fit() 00067 requires that 00068 \f[ 00069 \sum_i |g_i| < \mathrm{epsabs} 00070 \f] 00071 where \f$g_i\f$ is the \f$i\f$-th component of the gradient of 00072 the function \f$\Phi(x)\f$ where 00073 \f[ 00074 \Phi(x) = || F(x) ||^2 00075 \f] 00076 00077 \todo Properly generalize other vector types than 00078 ovector_base 00079 \todo Allow the user to specify the derivatives 00080 \todo Fix so that the user can specify automatic 00081 scaling of the fitting parameters, where the initial 00082 guess are used for scaling so that the fitting parameters 00083 are near unity. 00084 */ 00085 template<class param_t, class func_t=fit_funct<param_t>, 00086 class vec_t=ovector_base, 00087 class mat_t=omatrix_base, class bool_vec_t=bool *> class gsl_fit : 00088 public fit_base<param_t,func_t,vec_t,mat_t> { 00089 00090 public: 00091 00092 gsl_fit() { 00093 gsl_set_error_handler(err_hnd->gsl_hnd); 00094 max_iter=500; 00095 epsabs=1.0e-4; 00096 epsrel=1.0e-4; 00097 use_scaled=true; 00098 test_gradient=false; 00099 } 00100 00101 virtual ~gsl_fit() {} 00102 00103 /** \brief Fit the data specified in (xdat,ydat) to 00104 the function \c fitfun with the parameters in \c par. 00105 00106 The covariance matrix for the parameters is returned in \c covar 00107 and the value of \f$ \chi^2 \f$ is returned in \c chi2. 00108 */ 00109 virtual int fit(size_t ndat, vec_t &xdat, vec_t &ydat, vec_t &yerr, 00110 size_t npar, vec_t &par, mat_t &covar, double &chi2, 00111 param_t &pa, func_t &fitfun) { 00112 00113 int status, iter=0; 00114 size_t i, j; 00115 00116 gsl_vector *x=gsl_vector_alloc(npar); 00117 for(i=0;i<npar;i++) gsl_vector_set(x,i,par[i]); 00118 00119 // Create fitting function and parameters 00120 gsl_multifit_function_fdf f; 00121 gsl_vector *grad=gsl_vector_alloc(npar); 00122 func_par fpar={fitfun,&pa,ndat,&xdat,&ydat,&yerr,npar}; 00123 f.f=func; 00124 f.df=dfunc; 00125 f.fdf=fdfunc; 00126 f.n=ndat; 00127 f.p=npar; 00128 f.params=&fpar; 00129 00130 // Create fitting object 00131 gsl_multifit_fdfsolver *s; 00132 if (use_scaled) { 00133 s=gsl_multifit_fdfsolver_alloc 00134 (gsl_multifit_fdfsolver_lmsder,ndat,npar); 00135 } else { 00136 s=gsl_multifit_fdfsolver_alloc 00137 (gsl_multifit_fdfsolver_lmder,ndat,npar); 00138 } 00139 gsl_multifit_fdfsolver_set(s,&f,x); 00140 00141 // Perform the fit 00142 do { 00143 iter++; 00144 status=gsl_multifit_fdfsolver_iterate(s); 00145 00146 if (status) { 00147 break; 00148 } 00149 00150 status=gsl_multifit_test_delta(s->dx,s->x,epsabs,epsrel); 00151 if (test_gradient && status==gsl_continue) { 00152 gsl_multifit_gradient(s->J,s->x,grad); 00153 status=gsl_multifit_test_gradient(grad,epsabs); 00154 } 00155 00156 print_iter(npar,s->dx,s->x,iter,epsabs,epsrel); 00157 00158 } while (status == gsl_continue && iter < max_iter); 00159 00160 // Create covariance matrix and copy results to the 00161 // user-specified vector 00162 gsl_matrix *gcovar=gsl_matrix_alloc(npar,npar); 00163 gsl_multifit_covar(s->J,0.0,gcovar); 00164 for(i=0;i<npar;i++) { 00165 par[i]=gsl_vector_get(s->x,i); 00166 for(j=0;j<npar;j++) { 00167 covar.set(i,j,gsl_matrix_get(gcovar,i,j)); 00168 } 00169 } 00170 gsl_matrix_free(gcovar); 00171 00172 // Compute chi squared 00173 chi2=gsl_blas_dnrm2(s->f); 00174 chi2*=chi2; 00175 00176 // Free memory 00177 gsl_multifit_fdfsolver_free(s); 00178 gsl_vector_free(grad); 00179 00180 return 0; 00181 } 00182 00183 /// (default 500) 00184 int max_iter; 00185 00186 /// (default 1.0e-4) 00187 double epsabs; 00188 00189 /// (default 1.0e-4) 00190 double epsrel; 00191 00192 /// If true, test the gradient also (default false) 00193 bool test_gradient; 00194 00195 /** \brief Use the scaled routine if true (default true) 00196 */ 00197 bool use_scaled; 00198 00199 /// Return string denoting type ("gsl_fit") 00200 virtual const char *type() { return "gsl_fit"; } 00201 00202 #ifndef DOXYGEN_INTERNAL 00203 00204 protected: 00205 00206 /// Print the progress in the current iteration 00207 virtual int print_iter(int nv, gsl_vector *x, gsl_vector *dx, int iter, 00208 double l_epsabs, double l_epsrel) { 00209 if (this->verbose<=0) return 0; 00210 00211 int i; 00212 char ch; 00213 double dxmax, epsmin; 00214 00215 std::cout << "Iteration: " << iter << std::endl; 00216 text_out_file outs(&std::cout,79); 00217 outs.word_out("x:"); 00218 dxmax=gsl_vector_get(dx,0); 00219 epsmin=l_epsabs+l_epsrel*gsl_vector_get(x,0); 00220 for(i=0;i<nv;i++) { 00221 if (dxmax<gsl_vector_get(dx,i)) dxmax=gsl_vector_get(dx,i); 00222 if (epsmin>l_epsabs+l_epsrel*gsl_vector_get(x,i)) { 00223 epsmin=l_epsabs+l_epsrel*gsl_vector_get(x,i); 00224 } 00225 outs.double_out(gsl_vector_get(x,i)); 00226 } 00227 outs.end_line(); 00228 std::cout << " Val: " << dxmax << " Lim: " << epsmin << std::endl; 00229 if (this->verbose>1) { 00230 std::cout << "Press a key and type enter to continue. "; 00231 std::cin >> ch; 00232 } 00233 00234 return 0; 00235 } 00236 00237 /** \brief A structure for passing to the functions func(), dfunc(), 00238 and fdfunc() 00239 */ 00240 typedef struct { 00241 /// The function object 00242 func_t &f; 00243 /// The user-specified parameter 00244 param_t *vp; 00245 /// The number 00246 int ndat; 00247 /// The x values 00248 vec_t *xdat; 00249 /// The y values 00250 vec_t *ydat; 00251 /// The y uncertainties 00252 vec_t *yerr; 00253 /// The number of parameters 00254 int npar; 00255 } func_par; 00256 00257 /// Evaluate the function 00258 static int func(const gsl_vector *x, void *pa, gsl_vector *f) { 00259 func_par *fp=(func_par *)pa; 00260 int i, j; 00261 double yi; 00262 ovector xp(fp->npar); 00263 for(j=0;j<fp->npar;j++) xp[j]=gsl_vector_get(x,j); 00264 for(i=0;i<fp->ndat;i++) { 00265 fp->f(fp->npar,xp,(*fp->xdat)[i],yi,*(fp->vp)); 00266 gsl_vector_set(f,i,(yi-(*fp->ydat)[i])/(*fp->yerr)[i]); 00267 } 00268 00269 return 0; 00270 } 00271 00272 /// Evaluate the jacobian 00273 static int dfunc(const gsl_vector *x, void *pa, gsl_matrix *jac) { 00274 func_par *fp=(func_par *)pa; 00275 int i, j; 00276 double yi, ylo, yhi, xtemp, eps=1.0e-4; 00277 ovector xp(fp->npar); 00278 00279 for(j=0;j<fp->npar;j++) xp[j]=gsl_vector_get(x,j); 00280 for(j=0;j<fp->npar;j++) { 00281 for(i=0;i<fp->ndat;i++) { 00282 xtemp=xp[j]; 00283 xp[j]+=eps; 00284 fp->f(fp->npar,xp,(*fp->xdat)[i],yi,*(fp->vp)); 00285 yhi=(yi-(*fp->ydat)[i])/(*fp->yerr)[i]; 00286 xp[j]=xtemp; 00287 fp->f(fp->npar,xp,(*fp->xdat)[i],yi,*(fp->vp)); 00288 ylo=(yi-(*fp->ydat)[i])/(*fp->yerr)[i]; 00289 gsl_matrix_set(jac,i,j,(yhi-ylo)/eps); 00290 } 00291 } 00292 00293 return 0; 00294 } 00295 00296 /// Evaluate the function and the jacobian 00297 static int fdfunc(const gsl_vector *x, void *pa, gsl_vector *f, 00298 gsl_matrix *jac) { 00299 func(x,pa,f); 00300 dfunc(x,pa,jac); 00301 return 0; 00302 } 00303 00304 #endif 00305 00306 }; 00307 00308 #ifndef DOXYGENP 00309 } 00310 #endif 00311 00312 #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