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