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 protected: 00091 00092 #ifdef O2SCL_NEVER_DEFINED 00093 00094 int covar_new(size_t n, const mat_t &J, double epsrel, 00095 mat_t &covar) { 00096 00097 double tolr; 00098 00099 size_t i, j, k; 00100 size_t kmax = 0; 00101 00102 alloc_mat_t r; 00103 alloc_vec_t tau, norm; 00104 permutation perm; 00105 00106 am.allocate(r,n,n); 00107 av.allocate(tau,n); 00108 av.allocate(norm,n); 00109 perm.allocate(n); 00110 00111 int signum = 0; 00112 gsl_matrix_memcpy (r, J); 00113 gsl_linalg_QRPT_decomp (r, tau, perm, &signum, norm); 00114 00115 /* Form the inverse of R in the full upper triangle of R */ 00116 00117 tolr = epsrel * fabs(gsl_matrix_get(r, 0, 0)); 00118 00119 for (k = 0 ; k < n ; k++) { 00120 00121 double rkk = r[k][k]; 00122 00123 if (fabs(rkk) <= tolr) { 00124 break; 00125 } 00126 00127 r[k][k]=1.0/rkk; 00128 00129 for (j = 0; j < k ; j++) { 00130 00131 double t = r[j][k] / rkk; 00132 r[j][k]=0.0; 00133 00134 for (i = 0; i <= j; i++) { 00135 00136 double rik = r[i][k]; 00137 double rij = r[i][k]; 00138 00139 r[i][k]=rik - t * rij; 00140 } 00141 } 00142 kmax = k; 00143 } 00144 00145 /* Form the full upper triangle of the inverse of R^T R in the full 00146 upper triangle of R */ 00147 00148 for (k = 0; k <= kmax ; k++) { 00149 00150 for (j = 0; j < k; j++) { 00151 00152 double rjk = r[j][k]; 00153 00154 for (i = 0; i <= j ; i++) { 00155 00156 double rij = r[i][j]; 00157 double rik = r[i][k]; 00158 00159 r[i][j]= rij + rjk * rik; 00160 } 00161 } 00162 00163 { 00164 double t = r[k][k]; 00165 00166 for (i = 0; i <= k; i++) { 00167 00168 double rik = r[i][k]; 00169 00170 r[i][k]*=t; 00171 }; 00172 } 00173 } 00174 00175 /* Form the full lower triangle of the covariance matrix in the 00176 strict lower triangle of R and in w */ 00177 00178 for (j = 0 ; j < n ; j++) { 00179 00180 size_t pj = gsl_permutation_get (perm, j); 00181 00182 for (i = 0; i <= j; i++) { 00183 00184 size_t pi = gsl_permutation_get (perm, i); 00185 00186 double rij; 00187 00188 if (j > kmax) { 00189 00190 gsl_matrix_set (r, i, j, 0.0); 00191 rij = 0.0 ; 00192 00193 } else { 00194 00195 rij = gsl_matrix_get (r, i, j); 00196 } 00197 00198 if (pi > pj) { 00199 gsl_matrix_set (r, pi, pj, rij); 00200 00201 } else if (pi < pj) { 00202 00203 gsl_matrix_set (r, pj, pi, rij); 00204 } 00205 00206 } 00207 00208 { 00209 double rjj = r[j][j]; 00210 covar[pj][pj]=rjj; 00211 } 00212 } 00213 00214 00215 /* symmetrize the covariance matrix */ 00216 00217 for (j = 0 ; j < n ; j++) { 00218 00219 for (i = 0; i < j ; i++) { 00220 00221 double rji = gsl_matrix_get (r, j, i); 00222 00223 gsl_matrix_set (covar, j, i, rji); 00224 gsl_matrix_set (covar, i, j, rji); 00225 } 00226 } 00227 00228 gsl_matrix_free (r); 00229 gsl_permutation_free (perm); 00230 gsl_vector_free (tau); 00231 gsl_vector_free (norm); 00232 00233 return GSL_SUCCESS; 00234 } 00235 00236 #endif 00237 00238 public: 00239 00240 gsl_fit() { 00241 gsl_set_error_handler(err_hnd->gsl_hnd); 00242 max_iter=500; 00243 epsabs=1.0e-4; 00244 epsrel=1.0e-4; 00245 use_scaled=true; 00246 test_gradient=false; 00247 } 00248 00249 virtual ~gsl_fit() {} 00250 00251 /** \brief Fit the data specified in (xdat,ydat) to 00252 the function \c fitfun with the parameters in \c par. 00253 00254 The covariance matrix for the parameters is returned in \c covar 00255 and the value of \f$ \chi^2 \f$ is returned in \c chi2. 00256 */ 00257 virtual int fit(size_t ndat, vec_t &xdat, vec_t &ydat, vec_t &yerr, 00258 size_t npar, vec_t &par, mat_t &covar, double &chi2, 00259 param_t &pa, func_t &fitfun) { 00260 00261 int status, iter=0; 00262 size_t i, j; 00263 00264 gsl_vector *x=gsl_vector_alloc(npar); 00265 for(i=0;i<npar;i++) gsl_vector_set(x,i,par[i]); 00266 00267 // Create fitting function and parameters 00268 gsl_multifit_function_fdf f; 00269 gsl_vector *grad=gsl_vector_alloc(npar); 00270 func_par fpar={fitfun,&pa,ndat,&xdat,&ydat,&yerr,npar}; 00271 f.f=func; 00272 f.df=dfunc; 00273 f.fdf=fdfunc; 00274 f.n=ndat; 00275 f.p=npar; 00276 f.params=&fpar; 00277 00278 // Create fitting object 00279 gsl_multifit_fdfsolver *s; 00280 if (use_scaled) { 00281 s=gsl_multifit_fdfsolver_alloc 00282 (gsl_multifit_fdfsolver_lmsder,ndat,npar); 00283 } else { 00284 s=gsl_multifit_fdfsolver_alloc 00285 (gsl_multifit_fdfsolver_lmder,ndat,npar); 00286 } 00287 gsl_multifit_fdfsolver_set(s,&f,x); 00288 00289 // Perform the fit 00290 do { 00291 iter++; 00292 status=gsl_multifit_fdfsolver_iterate(s); 00293 00294 if (status) { 00295 break; 00296 } 00297 00298 status=gsl_multifit_test_delta(s->dx,s->x,epsabs,epsrel); 00299 if (test_gradient && status==gsl_continue) { 00300 gsl_multifit_gradient(s->J,s->x,grad); 00301 status=gsl_multifit_test_gradient(grad,epsabs); 00302 } 00303 00304 print_iter(npar,s->dx,s->x,iter,epsabs,epsrel); 00305 00306 } while (status == gsl_continue && iter < max_iter); 00307 00308 // Create covariance matrix and copy results to the 00309 // user-specified vector 00310 gsl_matrix *gcovar=gsl_matrix_alloc(npar,npar); 00311 gsl_multifit_covar(s->J,0.0,gcovar); 00312 for(i=0;i<npar;i++) { 00313 par[i]=gsl_vector_get(s->x,i); 00314 for(j=0;j<npar;j++) { 00315 covar.set(i,j,gsl_matrix_get(gcovar,i,j)); 00316 } 00317 } 00318 gsl_matrix_free(gcovar); 00319 00320 // Compute chi squared 00321 chi2=gsl_blas_dnrm2(s->f); 00322 chi2*=chi2; 00323 00324 // Free memory 00325 gsl_multifit_fdfsolver_free(s); 00326 gsl_vector_free(grad); 00327 00328 return 0; 00329 } 00330 00331 /// (default 500) 00332 int max_iter; 00333 00334 /// (default 1.0e-4) 00335 double epsabs; 00336 00337 /// (default 1.0e-4) 00338 double epsrel; 00339 00340 /// If true, test the gradient also (default false) 00341 bool test_gradient; 00342 00343 /** \brief Use the scaled routine if true (default true) 00344 */ 00345 bool use_scaled; 00346 00347 /// Return string denoting type ("gsl_fit") 00348 virtual const char *type() { return "gsl_fit"; } 00349 00350 #ifndef DOXYGEN_INTERNAL 00351 00352 protected: 00353 00354 /// Print the progress in the current iteration 00355 virtual int print_iter(int nv, gsl_vector *x, gsl_vector *dx, int iter, 00356 double l_epsabs, double l_epsrel) { 00357 if (this->verbose<=0) return 0; 00358 00359 int i; 00360 char ch; 00361 double dxmax, epsmin; 00362 00363 std::cout << "Iteration: " << iter << std::endl; 00364 text_out_file outs(&std::cout,79); 00365 outs.word_out("x:"); 00366 dxmax=gsl_vector_get(dx,0); 00367 epsmin=l_epsabs+l_epsrel*gsl_vector_get(x,0); 00368 for(i=0;i<nv;i++) { 00369 if (dxmax<gsl_vector_get(dx,i)) dxmax=gsl_vector_get(dx,i); 00370 if (epsmin>l_epsabs+l_epsrel*gsl_vector_get(x,i)) { 00371 epsmin=l_epsabs+l_epsrel*gsl_vector_get(x,i); 00372 } 00373 outs.double_out(gsl_vector_get(x,i)); 00374 } 00375 outs.end_line(); 00376 std::cout << " Val: " << dxmax << " Lim: " << epsmin << std::endl; 00377 if (this->verbose>1) { 00378 std::cout << "Press a key and type enter to continue. "; 00379 std::cin >> ch; 00380 } 00381 00382 return 0; 00383 } 00384 00385 /** \brief A structure for passing to the functions func(), dfunc(), 00386 and fdfunc() 00387 */ 00388 typedef struct { 00389 /// The function object 00390 func_t &f; 00391 /// The user-specified parameter 00392 param_t *vp; 00393 /// The number 00394 int ndat; 00395 /// The x values 00396 vec_t *xdat; 00397 /// The y values 00398 vec_t *ydat; 00399 /// The y uncertainties 00400 vec_t *yerr; 00401 /// The number of parameters 00402 int npar; 00403 } func_par; 00404 00405 /// Evaluate the function 00406 static int func(const gsl_vector *x, void *pa, gsl_vector *f) { 00407 func_par *fp=(func_par *)pa; 00408 int i, j; 00409 double yi; 00410 ovector xp(fp->npar); 00411 for(j=0;j<fp->npar;j++) xp[j]=gsl_vector_get(x,j); 00412 for(i=0;i<fp->ndat;i++) { 00413 fp->f(fp->npar,xp,(*fp->xdat)[i],yi,*(fp->vp)); 00414 gsl_vector_set(f,i,(yi-(*fp->ydat)[i])/(*fp->yerr)[i]); 00415 } 00416 00417 return 0; 00418 } 00419 00420 /// Evaluate the jacobian 00421 static int dfunc(const gsl_vector *x, void *pa, gsl_matrix *jac) { 00422 func_par *fp=(func_par *)pa; 00423 int i, j; 00424 double yi, ylo, yhi, xtemp, eps=1.0e-4; 00425 ovector xp(fp->npar); 00426 00427 for(j=0;j<fp->npar;j++) xp[j]=gsl_vector_get(x,j); 00428 for(j=0;j<fp->npar;j++) { 00429 for(i=0;i<fp->ndat;i++) { 00430 xtemp=xp[j]; 00431 xp[j]+=eps; 00432 fp->f(fp->npar,xp,(*fp->xdat)[i],yi,*(fp->vp)); 00433 yhi=(yi-(*fp->ydat)[i])/(*fp->yerr)[i]; 00434 xp[j]=xtemp; 00435 fp->f(fp->npar,xp,(*fp->xdat)[i],yi,*(fp->vp)); 00436 ylo=(yi-(*fp->ydat)[i])/(*fp->yerr)[i]; 00437 gsl_matrix_set(jac,i,j,(yhi-ylo)/eps); 00438 } 00439 } 00440 00441 return 0; 00442 } 00443 00444 /// Evaluate the function and the jacobian 00445 static int fdfunc(const gsl_vector *x, void *pa, gsl_vector *f, 00446 gsl_matrix *jac) { 00447 func(x,pa,f); 00448 dfunc(x,pa,jac); 00449 return 0; 00450 } 00451 00452 #endif 00453 00454 }; 00455 00456 #ifndef DOXYGENP 00457 } 00458 #endif 00459 00460 #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