FittingFunction Class Reference

#include <FittingFunction.h>

Inherited by FittingFunctionGaussian, FittingFunctionGaussianWithBackground, FittingFunctionLorentzian, FittingFunctionLorentzianWithBackground, FittingFunctionSigmoid, and FittingFunctionSigmoidWithBackground.

Inheritance diagram for FittingFunction:

Inheritance graph
[legend]
Collaboration diagram for FittingFunction:

Collaboration graph
[legend]
List of all members.

Public Member Functions

 FittingFunction ()
virtual ~FittingFunction ()
virtual PtrFunction getFunction (void) const =0
virtual PtrFunctionDerivate getFunctionDerivate (void) const =0
virtual PtrFunctionAndDerivate getFunctionAndDerivate (void) const =0
virtual int doFit (long points, gsl_vector *x, gsl_vector *y, gsl_vector *sigma, gsl_vector *guess, bool bEstimateGuess, int iSearchStoppingMethod, bool bUseScaled)=0
virtual int generateFunctionFit (double dValMinX, double dResolutionX, long lNbPoints, gsl_vector *funcX, gsl_vector *funcY)=0
virtual int generateFunctionFit (long lNbPoints, gsl_vector *x, gsl_vector *funcX, gsl_vector *funcY)=0
double getParameter (std::string sParameterName) const
double getParameter (int iIndiceParameter) const
double * getAllParameters () const
double getParameterError (std::string sParameterName) const
double getParameterError (int iIndiceParameter) const
double * getAllParametersErrors () const
int getNumberOfParameters () const
double getQualityFactor () const
double getSigma () const
double getFHWM () const
double getHWHM () const
double getSearchLimit () const
long getNbMaxIteration () const
int getNbParameters () const
long getNbIterations () const
void setNumberOfParameters (int iNbParameters)
void setSearchLimit (double dSearchValue)
void setNbMaxIteration (long lNbMaxIterationValue)
double computeStandardDeviation ()
double getStandardDeviation ()
virtual double computeValue (double dX, double dPosition, double dWidth, double dHeight, double dBackground)=0
virtual double computeDerivateValue (double dX, double dPosition, double dWidth, double dHeight)=0
virtual void initializeInitialsParameters (gsl_vector *vExperimentalDataX, gsl_vector *vExperimentalDataY, gsl_vector *vInitialGuess, bool bEstimate)=0
virtual void estimateInitialGuess (gsl_vector *vExperimentalDataX, gsl_vector *vExperimentalDataY)=0
double getInitialGuess (int iParameterIndex)
double computeMean (gsl_vector *vVector)
double computeVariance (gsl_vector *vVector)
double computeResidual (gsl_vector *vVectorExp, gsl_vector *vVectorPred)
double computeDeterminationCoefficient (gsl_vector *vExperimentalDataY, gsl_vector *vFittedDataY)
double computeObservedFStatistic (gsl_vector *vExperimentalDataY, gsl_vector *vFittedDataY)
int getIndexValue (double dWantedValue, gsl_vector *vVector)
virtual std::string getEquation () const =0

Protected Attributes

STRING2INT _mFunctionMap
int _iNbParameters
double * _pParameters
double * _pParametersErrors
double _dQualityFactor
double _dSigma
double _dHWHM
double _dFHWM
double _dSearchLimit
long _lNbMaxIteration
double _dStandardDeviation
long _lNbIterations
FittingData_cData
FittingData_cDataFitted
gsl_vector * _vInitialGuess
gsl_vector * _vExperimentalDataX
gsl_vector * _vExperimentalDataY

Detailed Description

Definition at line 58 of file FittingFunction.h.


Constructor & Destructor Documentation

FittingFunction::FittingFunction  ) 
 

Definition at line 17 of file FittingFunction.cpp.

00018 {
00019 }

FittingFunction::~FittingFunction  )  [virtual]
 

Definition at line 22 of file FittingFunction.cpp.

00023 {
00024 
00025 }


Member Function Documentation

virtual double FittingFunction::computeDerivateValue double  dX,
double  dPosition,
double  dWidth,
double  dHeight
[pure virtual]
 

Implemented in FittingFunctionGaussian, FittingFunctionGaussianWithBackground, FittingFunctionLorentzian, FittingFunctionLorentzianWithBackground, FittingFunctionSigmoid, and FittingFunctionSigmoidWithBackground.

double FittingFunction::computeDeterminationCoefficient gsl_vector *  vExperimentalDataY,
gsl_vector *  vFittedDataY
 

Definition at line 156 of file FittingFunction.cpp.

References computeResidual(), and computeVariance().

00157 {
00158         double dDeviation       = computeVariance(vExperimentalDataY);
00159         double dResidual        = computeResidual(vExperimentalDataY,vFittedDataY);
00160 
00161         double dRSQ                     = 1 - dResidual/dDeviation;
00162         return dRSQ;
00163 }

Here is the call graph for this function:

double FittingFunction::computeMean gsl_vector *  vVector  ) 
 

Definition at line 109 of file FittingFunction.cpp.

Referenced by computeVariance().

00110 {
00111         double dMean    = 0.0;
00112         double dSum             = 0.0;
00113         size_t dNbData  = vVector->size;
00114 
00115         for (unsigned int i=0;i<dNbData;i++)
00116         {
00117                 dSum+= gsl_vector_get(vVector,i);
00118         }
00119         dMean = dSum / dNbData;
00120         return dMean;
00121 }

double FittingFunction::computeObservedFStatistic gsl_vector *  vExperimentalDataY,
gsl_vector *  vFittedDataY
 

Definition at line 166 of file FittingFunction.cpp.

References computeResidual(), computeVariance(), and getNbParameters().

00167 {
00168         double dDeviation               = computeVariance(vExperimentalDataY);
00169         double dResidual                = computeResidual(vExperimentalDataY,vFittedDataY);
00170         size_t lNbData                  = vExperimentalDataY->size;
00171 
00172         double dSSR                             = dDeviation - dResidual;
00173         double dNbParameters    = getNbParameters();
00174         double dDegreeOfFreedom = lNbData - 1 - dNbParameters;
00175         
00176         double dMSE                             = dResidual/dDegreeOfFreedom;
00177         double dMSR                             = dSSR/dNbParameters;
00178         double dFStatistic              = dMSR/dMSE;
00179 
00180         return dFStatistic;
00181 }

Here is the call graph for this function:

double FittingFunction::computeResidual gsl_vector *  vVectorExp,
gsl_vector *  vVectorPred
 

Definition at line 140 of file FittingFunction.cpp.

Referenced by computeDeterminationCoefficient(), and computeObservedFStatistic().

00141 {
00142         double dSum                     = 0.0;
00143         double dDiff            = 0.0;
00144         double dResidual        = 0.0;
00145         size_t dNbData  = vVectorExp->size;
00146 
00147         for (unsigned int i=0;i<dNbData;i++)
00148         {
00149                 dDiff=gsl_vector_get(vVectorExp,i) - gsl_vector_get(vVectorPred,i);
00150                 dSum+= dDiff*dDiff;
00151         }
00152         dResidual = dSum;// / dNbData;
00153         return dResidual;
00154 } 

double FittingFunction::computeStandardDeviation  ) 
 

Definition at line 84 of file FittingFunction.cpp.

References _cData, _dStandardDeviation, getStandardDeviation(), FittingData::n, and FittingData::x.

00085 {
00086 
00087         size_t lNbdata  = _cData->n;
00088 
00089         double dSum = 0.0;
00090         double dSumSquare = 0.0;
00091         double dXValue;
00092         for (unsigned int i=0;i<lNbdata;i++)
00093         {
00094                 dXValue         = *(_cData->x+i);
00095                 dSum            += dXValue;
00096                 dSumSquare      += dXValue*dXValue;
00097         }
00098 
00099         double dInvNbData = 1/lNbdata;
00100         
00101         double dSqrt = dInvNbData*dSumSquare - dInvNbData*dInvNbData*dSum*dSum;
00102         
00103         _dStandardDeviation = sqrt(dSqrt);
00104 
00105         return getStandardDeviation();
00106 
00107 }

Here is the call graph for this function:

virtual double FittingFunction::computeValue double  dX,
double  dPosition,
double  dWidth,
double  dHeight,
double  dBackground
[pure virtual]
 

Implemented in FittingFunctionGaussian, FittingFunctionGaussianWithBackground, FittingFunctionLorentzian, FittingFunctionLorentzianWithBackground, FittingFunctionSigmoid, and FittingFunctionSigmoidWithBackground.

double FittingFunction::computeVariance gsl_vector *  vVector  ) 
 

Definition at line 123 of file FittingFunction.cpp.

References computeMean().

Referenced by computeDeterminationCoefficient(), and computeObservedFStatistic().

00124 {
00125         double dMean            = computeMean(vVector);
00126         double dSum                     = 0.0;
00127         double dDiff            = 0.0;
00128         double dDeviation       = 0.0;
00129         size_t dNbData  = vVector->size;
00130 
00131         for (unsigned int i=0;i<dNbData;i++)
00132         {
00133                 dDiff=gsl_vector_get(vVector,i) - dMean;
00134                 dSum+= dDiff*dDiff;
00135         }
00136         dDeviation = dSum;// / dNbData;
00137         return dDeviation;
00138 }

Here is the call graph for this function:

virtual int FittingFunction::doFit long  points,
gsl_vector *  x,
gsl_vector *  y,
gsl_vector *  sigma,
gsl_vector *  guess,
bool  bEstimateGuess,
int  iSearchStoppingMethod,
bool  bUseScaled
[pure virtual]
 

Implemented in FittingFunctionGaussian, FittingFunctionGaussianWithBackground, FittingFunctionLorentzian, FittingFunctionLorentzianWithBackground, FittingFunctionSigmoid, and FittingFunctionSigmoidWithBackground.

virtual void FittingFunction::estimateInitialGuess gsl_vector *  vExperimentalDataX,
gsl_vector *  vExperimentalDataY
[pure virtual]
 

Implemented in FittingFunctionGaussian, FittingFunctionGaussianWithBackground, FittingFunctionLorentzian, FittingFunctionLorentzianWithBackground, FittingFunctionSigmoid, and FittingFunctionSigmoidWithBackground.

virtual int FittingFunction::generateFunctionFit long  lNbPoints,
gsl_vector *  x,
gsl_vector *  funcX,
gsl_vector *  funcY
[pure virtual]
 

Implemented in FittingFunctionGaussian, FittingFunctionGaussianWithBackground, FittingFunctionLorentzian, FittingFunctionLorentzianWithBackground, FittingFunctionSigmoid, and FittingFunctionSigmoidWithBackground.

virtual int FittingFunction::generateFunctionFit double  dValMinX,
double  dResolutionX,
long  lNbPoints,
gsl_vector *  funcX,
gsl_vector *  funcY
[pure virtual]
 

Implemented in FittingFunctionGaussian, FittingFunctionGaussianWithBackground, FittingFunctionLorentzian, FittingFunctionLorentzianWithBackground, FittingFunctionSigmoid, and FittingFunctionSigmoidWithBackground.

double * FittingFunction::getAllParameters  )  const
 

Definition at line 58 of file FittingFunction.cpp.

References _pParameters.

00058                                                {
00059         
00060         return _pParameters;
00061 }

double * FittingFunction::getAllParametersErrors  )  const
 

Definition at line 73 of file FittingFunction.cpp.

References _pParametersErrors.

00073                                                       {
00074         return _pParametersErrors;
00075 }

virtual std::string FittingFunction::getEquation  )  const [pure virtual]
 

Implemented in FittingFunctionGaussian, FittingFunctionGaussianWithBackground, FittingFunctionLorentzian, FittingFunctionLorentzianWithBackground, FittingFunctionSigmoid, and FittingFunctionSigmoidWithBackground.

double FittingFunction::getFHWM  )  const [inline]
 

Definition at line 103 of file FittingFunction.h.

00103 {return _dFHWM;}

virtual PtrFunction FittingFunction::getFunction void   )  const [pure virtual]
 

Implemented in FittingFunctionGaussian, FittingFunctionGaussianWithBackground, FittingFunctionLorentzian, FittingFunctionLorentzianWithBackground, FittingFunctionSigmoid, and FittingFunctionSigmoidWithBackground.

virtual PtrFunctionAndDerivate FittingFunction::getFunctionAndDerivate void   )  const [pure virtual]
 

Implemented in FittingFunctionGaussian, FittingFunctionGaussianWithBackground, FittingFunctionLorentzian, FittingFunctionLorentzianWithBackground, FittingFunctionSigmoid, and FittingFunctionSigmoidWithBackground.

virtual PtrFunctionDerivate FittingFunction::getFunctionDerivate void   )  const [pure virtual]
 

Implemented in FittingFunctionGaussian, FittingFunctionGaussianWithBackground, FittingFunctionLorentzian, FittingFunctionLorentzianWithBackground, FittingFunctionSigmoid, and FittingFunctionSigmoidWithBackground.

double FittingFunction::getHWHM  )  const [inline]
 

Definition at line 104 of file FittingFunction.h.

00104 {return _dHWHM;}

int FittingFunction::getIndexValue double  dWantedValue,
gsl_vector *  vVector
 

Definition at line 184 of file FittingFunction.cpp.

00185 {
00186         size_t  lVectorSize                     = vVector->size;
00187         double  lHalfVectorSize         = 0.5*lVectorSize;
00188 
00189         size_t iIndexMax                        = gsl_vector_max_index(vVector);
00190         double dMaxValue                        = gsl_vector_get(vVector, iIndexMax);
00191   
00192         double dCurrentValue            = dMaxValue;
00193         double dLowValue                        = dMaxValue;
00194         double dHighValue                       = dMaxValue;
00195 
00196         size_t  iIndexHigh                      = iIndexMax;
00197         size_t  iIndexLow                       = iIndexMax;
00198 
00199         size_t iCurrentIndex            = iIndexMax;
00200         size_t iGoodIndex                       = iIndexMax;
00201 
00202         double dLowStep                         = dLowValue;
00203         double dHighStep                        = dHighValue;
00204 
00205 
00206 
00207     if (iIndexMax < lHalfVectorSize)
00208     {
00209                 //- look in increasing direction
00210                 while(iCurrentIndex < lVectorSize - 1 && dCurrentValue > dWantedValue)
00211                 {
00212                         dCurrentValue   = gsl_vector_get(vVector, ++iCurrentIndex);
00213                         iIndexHigh              = iCurrentIndex;
00214                         iIndexLow               = iCurrentIndex - 1;
00215                         dLowValue               = gsl_vector_get(vVector, iIndexLow);
00216                         dHighValue              = gsl_vector_get(vVector, iIndexHigh);
00217         
00218 //                      cout << "icurrentIndex = " << iCurrentIndex << " | dlowvalue = " << dLowValue << " | dhighvalue = " << dHighValue << endl;
00219 
00220                 }
00221                 dLowStep                        = ::fabs(dLowValue - dWantedValue); 
00222                 dHighStep                       = ::fabs(dHighValue - dWantedValue); 
00223                 if (dLowStep < dHighStep) iGoodIndex = iIndexLow;
00224                 else iGoodIndex = iIndexHigh;
00225     }
00226     else
00227     {
00228                 //- look in decreasing direction
00229                 while( iCurrentIndex > 0 && dCurrentValue > dWantedValue)
00230                 {
00231                         dCurrentValue =  gsl_vector_get(vVector, --iCurrentIndex);
00232                         iIndexLow               = iCurrentIndex;
00233                         iIndexHigh              = iCurrentIndex + 1;
00234                         dLowValue               = gsl_vector_get(vVector, iIndexLow);
00235                         dHighValue              = gsl_vector_get(vVector, iIndexHigh);
00236 //                      cout << "icurrentIndex = " << iCurrentIndex << " | dlowvalue = " << dLowValue << " | dhighvalue = " << dHighValue << endl;
00237                 }
00238                 dLowStep                        = ::fabs(dLowValue - dWantedValue); 
00239                 dHighStep                       = ::fabs(dHighValue - dWantedValue); 
00240                 if (dLowStep < dHighStep) iGoodIndex = iIndexLow;
00241                 else iGoodIndex = iIndexHigh;
00242 
00243     }
00244 
00245     if (iCurrentIndex != iIndexMax)
00246     {
00247                 //cout << "\nFind nearest index for " << dWantedValue << " --> " << iGoodIndex << endl;
00248                 return (int)iGoodIndex;
00249                 
00250     }
00251     else
00252     {
00253                 //iGoodIndex = -1;
00254                 cout << "\nFind nearest index for " << dWantedValue << " --> " << -1 << endl;
00255                 return -1;
00256     }
00257 
00258         
00259 
00260         return -1;
00261 }

double FittingFunction::getInitialGuess int  iParameterIndex  )  [inline]
 

Definition at line 130 of file FittingFunction.h.

00130 {return gsl_vector_get(_vInitialGuess,iParameterIndex);};

long FittingFunction::getNbIterations  )  const [inline]
 

Definition at line 111 of file FittingFunction.h.

00111 {return _lNbIterations;}

long FittingFunction::getNbMaxIteration  )  const [inline]
 

Definition at line 108 of file FittingFunction.h.

00108 {return _lNbMaxIteration;}

int FittingFunction::getNbParameters  )  const [inline]
 

Definition at line 110 of file FittingFunction.h.

Referenced by computeObservedFStatistic().

00110 {return _iNbParameters;}

int FittingFunction::getNumberOfParameters  )  const [inline]
 

Definition at line 97 of file FittingFunction.h.

00097 {return _iNbParameters;}

double FittingFunction::getParameter int  iIndiceParameter  )  const
 

Definition at line 49 of file FittingFunction.cpp.

References _iNbParameters, and _pParameters.

00049                                                               {
00050         
00051         if ((0 <= iIndiceParameter) && (iIndiceParameter < _iNbParameters)) {
00052                 return _pParameters[iIndiceParameter];
00053         }
00054         return 0;
00055 }

double FittingFunction::getParameter std::string  sParameterName  )  const
 

Definition at line 29 of file FittingFunction.cpp.

References _mFunctionMap, and _pParameters.

Referenced by FittingFunctionSigmoidWithBackground::estimateSigmoidIntersections(), FittingFunctionSigmoid::estimateSigmoidIntersections(), FittingFunctionSigmoidWithBackground::printResults(), FittingFunctionSigmoid::printResults(), FittingFunctionLorentzianWithBackground::printResults(), FittingFunctionLorentzian::printResults(), FittingFunctionGaussianWithBackground::printResults(), and FittingFunctionGaussian::printResults().

00029                                                                   {
00030 
00031         STRING2INT::const_iterator MyIterator;
00032         if((MyIterator = _mFunctionMap.find(sParameterName)) != _mFunctionMap.end()) {
00033                 return _pParameters[(*MyIterator).second];
00034         }
00035         else return 0;
00036 }

double FittingFunction::getParameterError int  iIndiceParameter  )  const
 

Definition at line 64 of file FittingFunction.cpp.

References _iNbParameters, and _pParametersErrors.

00064                                                                      {
00065         
00066         if ((0 <= iIndiceParameter) && (iIndiceParameter < _iNbParameters)) {
00067                 return _pParametersErrors[iIndiceParameter];
00068         }
00069         return 0;
00070 }

double FittingFunction::getParameterError std::string  sParameterName  )  const
 

Definition at line 39 of file FittingFunction.cpp.

References _mFunctionMap, and _pParametersErrors.

00039                                                                         {
00040 
00041         STRING2INT::const_iterator MyIterator;
00042         if((MyIterator = _mFunctionMap.find(sParameterName)) != _mFunctionMap.end()) {
00043                 return _pParametersErrors[(*MyIterator).second];
00044         }
00045         else return 0;
00046 }

double FittingFunction::getQualityFactor  )  const [inline]
 

Definition at line 99 of file FittingFunction.h.

00099 {return _dQualityFactor;}

double FittingFunction::getSearchLimit  )  const [inline]
 

Definition at line 106 of file FittingFunction.h.

00106 {return _dSearchLimit;}

double FittingFunction::getSigma  )  const [inline]
 

Definition at line 101 of file FittingFunction.h.

00101 {return _dSigma;}

double FittingFunction::getStandardDeviation  ) 
 

Definition at line 78 of file FittingFunction.cpp.

References _dStandardDeviation.

Referenced by computeStandardDeviation().

00079 {
00080         return _dStandardDeviation;
00081 }

virtual void FittingFunction::initializeInitialsParameters gsl_vector *  vExperimentalDataX,
gsl_vector *  vExperimentalDataY,
gsl_vector *  vInitialGuess,
bool  bEstimate
[pure virtual]
 

Implemented in FittingFunctionGaussian, FittingFunctionGaussianWithBackground, FittingFunctionLorentzian, FittingFunctionLorentzianWithBackground, FittingFunctionSigmoid, and FittingFunctionSigmoidWithBackground.

void FittingFunction::setNbMaxIteration long  lNbMaxIterationValue  )  [inline]
 

Definition at line 117 of file FittingFunction.h.

00117 {_lNbMaxIteration = lNbMaxIterationValue;}

void FittingFunction::setNumberOfParameters int  iNbParameters  )  [inline]
 

Definition at line 113 of file FittingFunction.h.

00113 {_iNbParameters = iNbParameters;}

void FittingFunction::setSearchLimit double  dSearchValue  )  [inline]
 

Definition at line 115 of file FittingFunction.h.

00115 {_dSearchLimit = dSearchValue;}


Member Data Documentation

FittingData* FittingFunction::_cData [protected]
 

Definition at line 177 of file FittingFunction.h.

Referenced by computeStandardDeviation().

FittingData* FittingFunction::_cDataFitted [protected]
 

Definition at line 179 of file FittingFunction.h.

double FittingFunction::_dFHWM [protected]
 

Definition at line 165 of file FittingFunction.h.

Referenced by FittingFunctionSigmoidWithBackground::printResults(), FittingFunctionSigmoid::printResults(), FittingFunctionLorentzianWithBackground::printResults(), FittingFunctionLorentzian::printResults(), FittingFunctionGaussianWithBackground::printResults(), and FittingFunctionGaussian::printResults().

double FittingFunction::_dHWHM [protected]
 

Definition at line 164 of file FittingFunction.h.

Referenced by FittingFunctionSigmoidWithBackground::printResults(), FittingFunctionSigmoid::printResults(), FittingFunctionLorentzianWithBackground::printResults(), FittingFunctionLorentzian::printResults(), FittingFunctionGaussianWithBackground::printResults(), and FittingFunctionGaussian::printResults().

double FittingFunction::_dQualityFactor [protected]
 

Definition at line 160 of file FittingFunction.h.

Referenced by FittingFunctionSigmoidWithBackground::printResults(), FittingFunctionSigmoid::printResults(), FittingFunctionLorentzianWithBackground::printResults(), FittingFunctionLorentzian::printResults(), FittingFunctionGaussianWithBackground::printResults(), and FittingFunctionGaussian::printResults().

double FittingFunction::_dSearchLimit [protected]
 

Definition at line 168 of file FittingFunction.h.

Referenced by FittingFunctionSigmoidWithBackground::doFit(), FittingFunctionSigmoid::doFit(), FittingFunctionLorentzianWithBackground::doFit(), FittingFunctionLorentzian::doFit(), FittingFunctionGaussianWithBackground::doFit(), FittingFunctionGaussian::doFit(), FittingFunctionGaussian::FittingFunctionGaussian(), FittingFunctionGaussianWithBackground::FittingFunctionGaussianWithBackground(), FittingFunctionLorentzian::FittingFunctionLorentzian(), FittingFunctionLorentzianWithBackground::FittingFunctionLorentzianWithBackground(), FittingFunctionSigmoid::FittingFunctionSigmoid(), and FittingFunctionSigmoidWithBackground::FittingFunctionSigmoidWithBackground().

double FittingFunction::_dSigma [protected]
 

Definition at line 162 of file FittingFunction.h.

Referenced by FittingFunctionSigmoidWithBackground::printResults(), FittingFunctionSigmoid::printResults(), FittingFunctionLorentzianWithBackground::printResults(), FittingFunctionLorentzian::printResults(), FittingFunctionGaussianWithBackground::printResults(), and FittingFunctionGaussian::printResults().

double FittingFunction::_dStandardDeviation [protected]
 

Definition at line 172 of file FittingFunction.h.

Referenced by computeStandardDeviation(), and getStandardDeviation().

int FittingFunction::_iNbParameters [protected]
 

Definition at line 153 of file FittingFunction.h.

Referenced by FittingFunctionSigmoidWithBackground::doFit(), FittingFunctionSigmoid::doFit(), FittingFunctionLorentzianWithBackground::doFit(), FittingFunctionLorentzian::doFit(), FittingFunctionGaussianWithBackground::doFit(), FittingFunctionGaussian::doFit(), FittingFunctionGaussian::FittingFunctionGaussian(), FittingFunctionGaussianWithBackground::FittingFunctionGaussianWithBackground(), FittingFunctionLorentzian::FittingFunctionLorentzian(), FittingFunctionLorentzianWithBackground::FittingFunctionLorentzianWithBackground(), FittingFunctionSigmoid::FittingFunctionSigmoid(), FittingFunctionSigmoidWithBackground::FittingFunctionSigmoidWithBackground(), getParameter(), and getParameterError().

long FittingFunction::_lNbIterations [protected]
 

Definition at line 174 of file FittingFunction.h.

long FittingFunction::_lNbMaxIteration [protected]
 

Definition at line 170 of file FittingFunction.h.

Referenced by FittingFunctionGaussian::FittingFunctionGaussian(), FittingFunctionGaussianWithBackground::FittingFunctionGaussianWithBackground(), FittingFunctionLorentzian::FittingFunctionLorentzian(), FittingFunctionLorentzianWithBackground::FittingFunctionLorentzianWithBackground(), FittingFunctionSigmoid::FittingFunctionSigmoid(), and FittingFunctionSigmoidWithBackground::FittingFunctionSigmoidWithBackground().

STRING2INT FittingFunction::_mFunctionMap [protected]
 

Definition at line 151 of file FittingFunction.h.

Referenced by FittingFunctionGaussian::FittingFunctionGaussian(), FittingFunctionGaussianWithBackground::FittingFunctionGaussianWithBackground(), FittingFunctionLorentzian::FittingFunctionLorentzian(), FittingFunctionLorentzianWithBackground::FittingFunctionLorentzianWithBackground(), FittingFunctionSigmoid::FittingFunctionSigmoid(), FittingFunctionSigmoidWithBackground::FittingFunctionSigmoidWithBackground(), getParameter(), and getParameterError().

double* FittingFunction::_pParameters [protected]
 

Definition at line 155 of file FittingFunction.h.

Referenced by FittingFunctionGaussian::FittingFunctionGaussian(), FittingFunctionGaussianWithBackground::FittingFunctionGaussianWithBackground(), FittingFunctionLorentzian::FittingFunctionLorentzian(), FittingFunctionLorentzianWithBackground::FittingFunctionLorentzianWithBackground(), FittingFunctionSigmoid::FittingFunctionSigmoid(), FittingFunctionSigmoidWithBackground::FittingFunctionSigmoidWithBackground(), FittingFunctionSigmoidWithBackground::generateFunctionFit(), FittingFunctionSigmoid::generateFunctionFit(), FittingFunctionLorentzianWithBackground::generateFunctionFit(), FittingFunctionLorentzian::generateFunctionFit(), FittingFunctionGaussianWithBackground::generateFunctionFit(), FittingFunctionGaussian::generateFunctionFit(), getAllParameters(), getParameter(), FittingFunctionGaussian::~FittingFunctionGaussian(), FittingFunctionGaussianWithBackground::~FittingFunctionGaussianWithBackground(), FittingFunctionLorentzian::~FittingFunctionLorentzian(), FittingFunctionLorentzianWithBackground::~FittingFunctionLorentzianWithBackground(), FittingFunctionSigmoid::~FittingFunctionSigmoid(), and FittingFunctionSigmoidWithBackground::~FittingFunctionSigmoidWithBackground().

double* FittingFunction::_pParametersErrors [protected]
 

Definition at line 157 of file FittingFunction.h.

Referenced by FittingFunctionGaussian::FittingFunctionGaussian(), FittingFunctionGaussianWithBackground::FittingFunctionGaussianWithBackground(), FittingFunctionLorentzian::FittingFunctionLorentzian(), FittingFunctionLorentzianWithBackground::FittingFunctionLorentzianWithBackground(), FittingFunctionSigmoid::FittingFunctionSigmoid(), FittingFunctionSigmoidWithBackground::FittingFunctionSigmoidWithBackground(), getAllParametersErrors(), getParameterError(), FittingFunctionSigmoidWithBackground::printResults(), FittingFunctionSigmoid::printResults(), FittingFunctionLorentzianWithBackground::printResults(), FittingFunctionLorentzian::printResults(), FittingFunctionGaussianWithBackground::printResults(), FittingFunctionGaussian::printResults(), FittingFunctionGaussian::~FittingFunctionGaussian(), FittingFunctionGaussianWithBackground::~FittingFunctionGaussianWithBackground(), FittingFunctionLorentzian::~FittingFunctionLorentzian(), FittingFunctionLorentzianWithBackground::~FittingFunctionLorentzianWithBackground(), FittingFunctionSigmoid::~FittingFunctionSigmoid(), and FittingFunctionSigmoidWithBackground::~FittingFunctionSigmoidWithBackground().

gsl_vector* FittingFunction::_vExperimentalDataX [protected]
 

Definition at line 183 of file FittingFunction.h.

gsl_vector* FittingFunction::_vExperimentalDataY [protected]
 

Definition at line 184 of file FittingFunction.h.

gsl_vector* FittingFunction::_vInitialGuess [protected]
 

Definition at line 181 of file FittingFunction.h.

Referenced by FittingFunctionSigmoidWithBackground::doFit(), FittingFunctionSigmoid::doFit(), FittingFunctionLorentzianWithBackground::doFit(), FittingFunctionLorentzian::doFit(), FittingFunctionGaussianWithBackground::doFit(), FittingFunctionGaussian::doFit(), FittingFunctionGaussianWithBackground::estimateInitialGuess(), FittingFunctionGaussian::estimateInitialGuess(), FittingFunctionGaussian::FittingFunctionGaussian(), FittingFunctionGaussianWithBackground::FittingFunctionGaussianWithBackground(), FittingFunctionLorentzian::FittingFunctionLorentzian(), FittingFunctionLorentzianWithBackground::FittingFunctionLorentzianWithBackground(), FittingFunctionSigmoid::FittingFunctionSigmoid(), FittingFunctionSigmoidWithBackground::FittingFunctionSigmoidWithBackground(), FittingFunctionSigmoidWithBackground::initializeInitialsParameters(), FittingFunctionSigmoid::initializeInitialsParameters(), FittingFunctionLorentzianWithBackground::initializeInitialsParameters(), FittingFunctionLorentzian::initializeInitialsParameters(), FittingFunctionGaussianWithBackground::initializeInitialsParameters(), FittingFunctionGaussian::initializeInitialsParameters(), FittingFunctionGaussian::~FittingFunctionGaussian(), FittingFunctionGaussianWithBackground::~FittingFunctionGaussianWithBackground(), FittingFunctionLorentzian::~FittingFunctionLorentzian(), FittingFunctionLorentzianWithBackground::~FittingFunctionLorentzianWithBackground(), FittingFunctionSigmoid::~FittingFunctionSigmoid(), and FittingFunctionSigmoidWithBackground::~FittingFunctionSigmoidWithBackground().


The documentation for this class was generated from the following files:
Generated on Tue Apr 14 09:34:43 2009 for Data Fitting Library by  doxygen 1.4.5