Home  · Classes  · Annotated Classes  · Modules  · Members  · Namespaces  · Related Pages

LevMarqFitter1D.h (Maintainer: Marcel Grunert)

Go to the documentation of this file.
00001 // -*- Mode: C++; tab-width: 2; -*-
00002 // vi: set ts=2:
00003 //
00004 // --------------------------------------------------------------------------
00005 //                   OpenMS Mass Spectrometry Framework
00006 // --------------------------------------------------------------------------
00007 //  Copyright (C) 2003-2008 -- Oliver Kohlbacher, Knut Reinert
00008 //
00009 //  This library is free software; you can redistribute it and/or
00010 //  modify it under the terms of the GNU Lesser General Public
00011 //  License as published by the Free Software Foundation; either
00012 //  version 2.1 of the License, or (at your option) any later version.
00013 //
00014 //  This library is distributed in the hope that it will be useful,
00015 //  but WITHOUT ANY WARRANTY; without even the implied warranty of
00016 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00017 //  Lesser General Public License for more details.
00018 //
00019 //  You should have received a copy of the GNU Lesser General Public
00020 //  License along with this library; if not, write to the Free Software
00021 //  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00022 //
00023 // --------------------------------------------------------------------------
00024 // $Maintainer: Marcel Grunert $
00025 // --------------------------------------------------------------------------
00026 
00027 #ifndef OPENMS_TRANSFORMATIONS_FEATUREFINDER_LEVMARQFITTER1D_H
00028 #define OPENMS_TRANSFORMATIONS_FEATUREFINDER_LEVMARQFITTER1D_H
00029 
00030 #include <OpenMS/TRANSFORMATIONS/FEATUREFINDER/Fitter1D.h>
00031 #include <OpenMS/MATH/MISC/LinearInterpolation.h>
00032 #include <OpenMS/MATH/MISC/MathFunctions.h>
00033 
00034 #include <gsl/gsl_rng.h> // gsl random number generators
00035 #include <gsl/gsl_randist.h> // gsl random number distributions
00036 #include <gsl/gsl_vector.h> // gsl vector and matrix definitions
00037 #include <gsl/gsl_multifit_nlin.h> // gsl multidimensional fitting
00038 #include <gsl/gsl_blas.h> // gsl linear algebra stuff
00039 
00040 namespace OpenMS
00041 {
00042   
00048     class LevMarqFitter1D
00049     : public Fitter1D
00050     {
00051 
00052       public:
00053       
00054         typedef std::vector < double > ContainerType;
00055         
00057         LevMarqFitter1D()
00058         : Fitter1D()
00059         {
00060           this->defaults_.setValue( "max_iteration", 500, "Maximum number of iterations using by Levenberg-Marquardt algorithm.", true );
00061           this->defaults_.setValue( "deltaAbsError", 0.0001, "Absolute error used by the Levenberg-Marquardt algorithm.", true );
00062           this->defaults_.setValue( "deltaRelError", 0.0001, "Relative error used by the Levenberg-Marquardt algorithm.", true );
00063         }
00064   
00066         LevMarqFitter1D(const LevMarqFitter1D& source)
00067         : Fitter1D(source),
00068             max_iteration_(source.max_iteration_),
00069             abs_error_(source.abs_error_),
00070             rel_error_(source.rel_error_)
00071         {
00072         }
00073                           
00075         virtual ~LevMarqFitter1D()
00076         {
00077         }
00078   
00080         virtual LevMarqFitter1D& operator = (const LevMarqFitter1D& source)
00081         {
00082             if (&source ==this) return *this;
00083       
00084             Fitter1D::operator = (source);
00085             max_iteration_ = source.max_iteration_;
00086             abs_error_ = source.abs_error_;
00087             rel_error_ = source.rel_error_;
00088       
00089             return *this;
00090         }
00091            
00092     protected:
00093       
00095         Int gsl_status_;
00097         bool symmetric_;
00099         Int max_iteration_;
00101 
00102         CoordinateType abs_error_;
00104         CoordinateType rel_error_;
00106         Math::BasicStatistics<Real> stat_;
00107        
00111         virtual void printState_(Int iter, gsl_multifit_fdfsolver * s) = 0;  
00112      
00114         const String getGslStatus_()
00115         {
00116           return gsl_strerror( gsl_status_ );
00117         }
00118            
00120         void optimize_(const RawDataArrayType& set, Int num_params, CoordinateType x_init[],
00121                        Int (* residual)(const gsl_vector * x, void * params, gsl_vector * f),
00122                        Int (* jacobian)(const gsl_vector * x, void * params, gsl_matrix * J),
00123                        Int (* evaluate)(const gsl_vector * x, void * params, gsl_vector * f, gsl_matrix * J), void* advanced_params )
00124         {
00125           
00126           const gsl_multifit_fdfsolver_type * T;
00127           gsl_multifit_fdfsolver *s;
00128     
00129           Int status;
00130           Int iter = 0;
00131           const UInt n = set.size();
00132     
00133           // number of parameter to be optimize
00134           UInt p = num_params;
00135           
00136           // gsl always expects N>=p or default gsl error handler invoked, 
00137           // cause Jacobian be rectangular M x N with M>=N
00138           if ( n < p ) throw UnableToFit( __FILE__, __LINE__, __PRETTY_FUNCTION__, "UnableToFit-FinalSet", "Skipping feature, gsl always expects N>=p" );
00139     
00140       // allocate space for a covariance matrix of size p by p
00141           gsl_matrix *covar = gsl_matrix_alloc( p, p );
00142           gsl_multifit_function_fdf f;
00143     
00144           gsl_vector_view x;
00145           x = gsl_vector_view_array( x_init, p );
00146         
00147           const gsl_rng_type * type;
00148           gsl_rng * r;
00149           gsl_rng_env_setup();
00150           type = gsl_rng_default;
00151           r = gsl_rng_alloc ( type );
00152           
00153           // set up the function to be fit
00154           f.f = (residual); // the function of residuals
00155           f.df = (jacobian); // the gradient of this function
00156           f.fdf = (evaluate); // combined function and gradient
00157           f.n = set.size(); // number of points in the data set
00158           f.p = p; // number of parameters in the fit function
00159           f.params = advanced_params; // // structure with the data and error bars
00160           
00161           T = gsl_multifit_fdfsolver_lmsder;
00162           s = gsl_multifit_fdfsolver_alloc( T, n, p );
00163           gsl_multifit_fdfsolver_set( s, &f, &x.vector );
00164 
00165 #ifdef DEBUG_FEATUREFINDER
00166           printState_(iter, s);
00167 #endif
00168   
00169           // this is the loop for fitting
00170           do
00171           {
00172             iter++;
00173             
00174             // perform a single iteration of the fitting routine
00175             status = gsl_multifit_fdfsolver_iterate ( s );
00176           
00177 #ifdef DEBUG_FEATUREFINDER
00178             // customized routine to print out current parameters
00179             printState_(iter, s);
00180 #endif
00181            
00182             /* check if solver is stuck */
00183             if ( status ) break;
00184             
00185             // test for convergence with an absolute and relative error
00186             status = gsl_multifit_test_delta( s->dx, s->x, abs_error_, rel_error_ );
00187           }
00188           while ( status == GSL_CONTINUE && iter < max_iteration_ );
00189         
00190           // This function uses Jacobian matrix J to compute the covariance matrix of the best-fit parameters, covar. 
00191           // The parameter epsrel (0.0) is used to remove linear-dependent columns when J is rank deficient.
00192           gsl_multifit_covar( s->J, 0.0, covar );
00193   
00194 #ifdef DEBUG_FEATUREFINDER
00195           gsl_matrix_fprintf( stdout, covar, "covar %g" );
00196 #endif
00197   
00198 #define FIT(i) gsl_vector_get(s->x, i)
00199 #define ERR(i) sqrt(gsl_matrix_get(covar,i,i))
00200   
00201           // Set GSl status
00202           gsl_status_ = status;
00203         
00204 #ifdef DEBUG_FEATUREFINDER
00205           { 
00206             // chi-squared value
00207             DoubleReal chi = gsl_blas_dnrm2(s->f);
00208             DoubleReal dof = n - p;
00209             DoubleReal c = GSL_MAX_DBL(1, chi / sqrt(dof)); 
00210                                 
00211             printf("chisq/dof = %g\n",  pow(chi, 2.0) / dof);
00212                               
00213             for (UInt i=0; i<p; ++i)
00214             {
00215               std::cout << i;
00216         printf(".Parameter = %.5f +/- %.5f\n", FIT( i ), c*ERR( i ) );
00217             }           
00218           }
00219 #endif
00220         
00221           // set optimized parameter  
00222           for (UInt i = 0; i < p; ++i)
00223           {
00224             x_init[i] = FIT( i );
00225           }
00226           
00227           gsl_multifit_fdfsolver_free (s);
00228           gsl_matrix_free (covar);
00229           gsl_rng_free (r);
00230         }       
00231 
00232         void updateMembers_()
00233         {
00234             Fitter1D::updateMembers_();
00235             max_iteration_ = this->param_.getValue("max_iteration");
00236             abs_error_ = this->param_.getValue("deltaAbsError");
00237             rel_error_ = this->param_.getValue("deltaRelError");
00238         }
00239   };
00240 }
00241 
00242 #endif // OPENMS_TRANSFORMATIONS_FEATUREFINDER_LEVMARQFITTER1D_H

Generated Tue Apr 1 15:36:35 2008 -- using doxygen 1.5.4 OpenMS / TOPP 1.1