ROOT  6.06/09
Reference Guide
TSpectrumFit.h
Go to the documentation of this file.
1 // @(#)root/spectrum:$Id$
2 // Author: Miroslav Morhac 25/09/06
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2006, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 #ifndef ROOT_TSpectrumFit
12 #define ROOT_TSpectrumFit
13 
14 //////////////////////////////////////////////////////////////////////////
15 // //
16 // TSpectrumFit //
17 // //
18 // Class for fitting 1D spectra using AWMI (algorithm without matrix //
19 // inversion) and conjugate gradient algorithms for symmetrical //
20 // matrices (Stiefel-Hestens method). AWMI method allows to fit //
21 // simulaneously 100s up to 1000s peaks. Stiefel method is very stable, //
22 // it converges faster, but is more time consuming //
23 // //
24 //////////////////////////////////////////////////////////////////////////
25 
26 #ifndef ROOT_TNamed
27 #include "TNamed.h"
28 #endif
29 
30 class TH1;
31 
32 class TSpectrumFit : public TNamed {
33 protected:
34  Int_t fNPeaks; //number of peaks present in fit, input parameter, it should be > 0
35  Int_t fNumberIterations; //number of iterations in fitting procedure, input parameter, it should be > 0
36  Int_t fXmin; //first fitted channel
37  Int_t fXmax; //last fitted channel
38  Int_t fStatisticType; //type of statistics, possible values kFitOptimChiCounts (chi square statistics with counts as weighting coefficients), kFitOptimChiFuncValues (chi square statistics with function values as weighting coefficients),kFitOptimMaxLikelihood
39  Int_t fAlphaOptim; //optimization of convergence algorithm, possible values kFitAlphaHalving, kFitAlphaOptimal
40  Int_t fPower; //possible values kFitPower2,4,6,8,10,12, for details see references. It applies only for Awmi fitting function.
41  Int_t fFitTaylor; //order of Taylor expansion, possible values kFitTaylorOrderFirst, kFitTaylorOrderSecond. It applies only for Awmi fitting function.
42  Double_t fAlpha; //convergence coefficient, input parameter, it should be positive number and <=1, for details see references
43  Double_t fChi; //here the fitting functions return resulting chi square
44  Double_t *fPositionInit; //[fNPeaks] array of initial values of peaks positions, input parameters
45  Double_t *fPositionCalc; //[fNPeaks] array of calculated values of fitted positions, output parameters
46  Double_t *fPositionErr; //[fNPeaks] array of position errors
47  Double_t *fAmpInit; //[fNPeaks] array of initial values of peaks amplitudes, input parameters
48  Double_t *fAmpCalc; //[fNPeaks] array of calculated values of fitted amplitudes, output parameters
49  Double_t *fAmpErr; //[fNPeaks] array of amplitude errors
50  Double_t *fArea; //[fNPeaks] array of calculated areas of peaks
51  Double_t *fAreaErr; //[fNPeaks] array of errors of peak areas
52  Double_t fSigmaInit; //initial value of sigma parameter
53  Double_t fSigmaCalc; //calculated value of sigma parameter
54  Double_t fSigmaErr; //error value of sigma parameter
55  Double_t fTInit; //initial value of t parameter (relative amplitude of tail), for details see html manual and references
56  Double_t fTCalc; //calculated value of t parameter
57  Double_t fTErr; //error value of t parameter
58  Double_t fBInit; //initial value of b parameter (slope), for details see html manual and references
59  Double_t fBCalc; //calculated value of b parameter
60  Double_t fBErr; //error value of b parameter
61  Double_t fSInit; //initial value of s parameter (relative amplitude of step), for details see html manual and references
62  Double_t fSCalc; //calculated value of s parameter
63  Double_t fSErr; //error value of s parameter
64  Double_t fA0Init; //initial value of background a0 parameter(backgroud is estimated as a0+a1*x+a2*x*x)
65  Double_t fA0Calc; //calculated value of background a0 parameter
66  Double_t fA0Err; //error value of background a0 parameter
67  Double_t fA1Init; //initial value of background a1 parameter(backgroud is estimated as a0+a1*x+a2*x*x)
68  Double_t fA1Calc; //calculated value of background a1 parameter
69  Double_t fA1Err; //error value of background a1 parameter
70  Double_t fA2Init; //initial value of background a2 parameter(backgroud is estimated as a0+a1*x+a2*x*x)
71  Double_t fA2Calc; //calculated value of background a2 parameter
72  Double_t fA2Err; //error value of background a2 parameter
73  Bool_t *fFixPosition; //[fNPeaks] array of logical values which allow to fix appropriate positions (not fit). However they are present in the estimated functional
74  Bool_t *fFixAmp; //[fNPeaks] array of logical values which allow to fix appropriate amplitudes (not fit). However they are present in the estimated functional
75  Bool_t fFixSigma; //logical value of sigma parameter, which allows to fix the parameter (not to fit).
76  Bool_t fFixT; //logical value of t parameter, which allows to fix the parameter (not to fit).
77  Bool_t fFixB; //logical value of b parameter, which allows to fix the parameter (not to fit).
78  Bool_t fFixS; //logical value of s parameter, which allows to fix the parameter (not to fit).
79  Bool_t fFixA0; //logical value of a0 parameter, which allows to fix the parameter (not to fit).
80  Bool_t fFixA1; //logical value of a1 parameter, which allows to fix the parameter (not to fit).
81  Bool_t fFixA2; //logical value of a2 parameter, which allows to fix the parameter (not to fit).
82 
83 public:
84  enum {
99  };
100  TSpectrumFit(void); //default constructor
101  TSpectrumFit(Int_t numberPeaks);
102  virtual ~TSpectrumFit();
103 
104  //auxiliary functions for 1. parameter fit functions
105 protected:
110  Double_t Derb(Int_t num_of_fitted_peaks,Double_t i,const Double_t* parameter,Double_t sigma,Double_t t,Double_t b);
112  Double_t Derdersigma(Int_t num_of_fitted_peaks,Double_t i,const Double_t* parameter,Double_t sigma);
119  Double_t Ders(Int_t num_of_fitted_peaks,Double_t i,const Double_t* parameter,Double_t sigma);
120  Double_t Dersigma(Int_t num_of_fitted_peaks,Double_t i,const Double_t* parameter,Double_t sigma,Double_t t,Double_t s,Double_t b);
121  Double_t Dert(Int_t num_of_fitted_peaks,Double_t i,const Double_t* parameter,Double_t sigma,Double_t b);
124  Double_t Shape(Int_t num_of_fitted_peaks,Double_t i,const Double_t *parameter,Double_t sigma,Double_t t,Double_t s,Double_t b,Double_t a0,Double_t a1,Double_t a2);
125  void StiefelInversion(Double_t **a,Int_t rozmer);
126 
127 public:
128  void FitAwmi(Double_t *source);
129  void FitStiefel(Double_t *source);
130  Double_t *GetAmplitudes() const {return fAmpCalc;}
132  Double_t *GetAreas() const {return fArea;}
133  Double_t *GetAreasErrors() const {return fAreaErr;}
134  void GetBackgroundParameters(Double_t &a0, Double_t &a0Err, Double_t &a1, Double_t &a1Err, Double_t &a2, Double_t &a2Err);
135  Double_t GetChi() const {return fChi;}
138  void GetSigma(Double_t &sigma, Double_t &sigmaErr);
139  void GetTailParameters(Double_t &t, Double_t &tErr, Double_t &b, Double_t &bErr, Double_t &s, Double_t &sErr);
140  void SetBackgroundParameters(Double_t a0Init, Bool_t fixA0, Double_t a1Init, Bool_t fixA1, Double_t a2Init, Bool_t fixA2);
141  void SetFitParameters(Int_t xmin,Int_t xmax, Int_t numberIterations, Double_t alpha, Int_t statisticType, Int_t alphaOptim, Int_t power, Int_t fitTaylor);
142  void SetPeakParameters(Double_t sigma, Bool_t fixSigma, const Double_t *positionInit, const Bool_t *fixPosition, const Double_t *ampInit, const Bool_t *fixAmp);
143  void SetTailParameters(Double_t tInit, Bool_t fixT, Double_t bInit, Bool_t fixB, Double_t sInit, Bool_t fixS);
144 
145  ClassDef(TSpectrumFit,1) //Spectrum Fitter using algorithm without matrix inversion and conjugate gradient method for symmetrical matrices (Stiefel-Hestens method)
146 };
147 
148 #endif
149 
Double_t Shape(Int_t num_of_fitted_peaks, Double_t i, const Double_t *parameter, Double_t sigma, Double_t t, Double_t s, Double_t b, Double_t a0, Double_t a1, Double_t a2)
AUXILIARY FUNCTION // // This function calculates peaks shape function (see manual) // Function param...
Double_t fAlpha
Definition: TSpectrumFit.h:42
float xmin
Definition: THbookFile.cxx:93
Double_t * fAreaErr
Definition: TSpectrumFit.h:51
virtual ~TSpectrumFit()
destructor
Bool_t fFixSigma
Definition: TSpectrumFit.h:75
Double_t Derpb(Double_t a, Double_t sigma, Double_t t, Double_t b)
AUXILIARY FUNCTION // // This function calculates derivative of the area of peak // according to b pa...
void SetPeakParameters(Double_t sigma, Bool_t fixSigma, const Double_t *positionInit, const Bool_t *fixPosition, const Double_t *ampInit, const Bool_t *fixAmp)
SETTER FUNCTION.
Double_t * fAmpInit
Definition: TSpectrumFit.h:47
Double_t * fPositionCalc
Definition: TSpectrumFit.h:45
Double_t Area(Double_t a, Double_t sigma, Double_t t, Double_t b)
AUXILIARY FUNCTION // // This function calculates area of a peak // Function parameters: // -a-amplit...
Double_t * GetPositionsErrors() const
Definition: TSpectrumFit.h:137
void StiefelInversion(Double_t **a, Int_t rozmer)
AUXILIARY FUNCTION // // This function calculates solution of the system of linear equations...
Int_t fFitTaylor
Definition: TSpectrumFit.h:41
Int_t fStatisticType
Definition: TSpectrumFit.h:38
int Int_t
Definition: RtypesCore.h:41
Double_t Dert(Int_t num_of_fitted_peaks, Double_t i, const Double_t *parameter, Double_t sigma, Double_t b)
AUXILIARY FUNCTION // // This function calculates derivative of peaks shape function (see manual) // ...
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
Double_t * GetAmplitudesErrors() const
Definition: TSpectrumFit.h:131
Double_t Derpsigma(Double_t a, Double_t t, Double_t b)
Double_t * fPositionErr
Definition: TSpectrumFit.h:46
Double_t fTErr
Definition: TSpectrumFit.h:57
void GetTailParameters(Double_t &t, Double_t &tErr, Double_t &b, Double_t &bErr, Double_t &s, Double_t &sErr)
GETTER FUNCTION.
Double_t fA0Init
Definition: TSpectrumFit.h:64
Bool_t fFixS
Definition: TSpectrumFit.h:78
Int_t fAlphaOptim
Definition: TSpectrumFit.h:39
Double_t * GetAmplitudes() const
Definition: TSpectrumFit.h:130
Double_t x[n]
Definition: legend1.C:17
#define ClassDef(name, id)
Definition: Rtypes.h:254
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:33
Double_t fA1Calc
Definition: TSpectrumFit.h:68
Double_t fA2Init
Definition: TSpectrumFit.h:70
void SetBackgroundParameters(Double_t a0Init, Bool_t fixA0, Double_t a1Init, Bool_t fixA1, Double_t a2Init, Bool_t fixA2)
SETTER FUNCTION.
Double_t * fPositionInit
Definition: TSpectrumFit.h:44
Double_t Dera2(Double_t i)
derivative of background according to a2
Double_t fSigmaCalc
Definition: TSpectrumFit.h:53
Int_t fNumberIterations
Definition: TSpectrumFit.h:35
Double_t Derfc(Double_t x)
AUXILIARY FUNCTION // // This function calculates derivative of error function of x...
Double_t fSCalc
Definition: TSpectrumFit.h:62
Double_t fChi
Definition: TSpectrumFit.h:43
Double_t fA2Err
Definition: TSpectrumFit.h:72
Bool_t fFixB
Definition: TSpectrumFit.h:77
Double_t fA1Init
Definition: TSpectrumFit.h:67
Double_t fTInit
Definition: TSpectrumFit.h:55
Double_t fSErr
Definition: TSpectrumFit.h:63
Double_t fBCalc
Definition: TSpectrumFit.h:59
Double_t * fAmpErr
Definition: TSpectrumFit.h:49
Double_t Dersigma(Int_t num_of_fitted_peaks, Double_t i, const Double_t *parameter, Double_t sigma, Double_t t, Double_t s, Double_t b)
AUXILIARY FUNCTION // // This function calculates derivative of peaks shape function (see manual) // ...
Double_t fBInit
Definition: TSpectrumFit.h:58
Double_t Erfc(Double_t x)
Double_t fBErr
Definition: TSpectrumFit.h:60
Double_t Deramp(Double_t i, Double_t i0, Double_t sigma, Double_t t, Double_t s, Double_t b)
AUXILIARY FUNCTION // // This function calculates derivative of peak shape function (see manual) // a...
float xmax
Definition: THbookFile.cxx:93
Double_t fA1Err
Definition: TSpectrumFit.h:69
Double_t * GetAreasErrors() const
Definition: TSpectrumFit.h:133
void SetFitParameters(Int_t xmin, Int_t xmax, Int_t numberIterations, Double_t alpha, Int_t statisticType, Int_t alphaOptim, Int_t power, Int_t fitTaylor)
SETTER FUNCTION.
Double_t Deri0(Double_t i, Double_t amp, Double_t i0, Double_t sigma, Double_t t, Double_t s, Double_t b)
AUXILIARY FUNCTION // // This function calculates derivative of peak shape function (see manual) // a...
Double_t * GetPositions() const
Definition: TSpectrumFit.h:136
Bool_t fFixA0
Definition: TSpectrumFit.h:79
Double_t Derpa(Double_t sigma, Double_t t, Double_t b)
AUXILIARY FUNCTION // // This function calculates derivative of the area of peak // according to its ...
Double_t Ders(Int_t num_of_fitted_peaks, Double_t i, const Double_t *parameter, Double_t sigma)
AUXILIARY FUNCTION // // This function calculates derivative of peaks shape function (see manual) // ...
Double_t fA0Err
Definition: TSpectrumFit.h:66
double Double_t
Definition: RtypesCore.h:55
Double_t Dera1(Double_t i)
derivative of background according to a1
Double_t Derderi0(Double_t i, Double_t amp, Double_t i0, Double_t sigma)
AUXILIARY FUNCTION // // This function calculates second derivative of peak shape function // (see ma...
The TH1 histogram class.
Definition: TH1.h:80
Double_t GetChi() const
Definition: TSpectrumFit.h:135
Double_t Derdersigma(Int_t num_of_fitted_peaks, Double_t i, const Double_t *parameter, Double_t sigma)
AUXILIARY FUNCTION // // This function calculates second derivative of peaks shape function // (see m...
Double_t Derpt(Double_t a, Double_t sigma, Double_t b)
AUXILIARY FUNCTION // // This function calculates derivative of the area of peak // according to t pa...
Bool_t * fFixPosition
Definition: TSpectrumFit.h:73
Bool_t * fFixAmp
Definition: TSpectrumFit.h:74
Double_t fSInit
Definition: TSpectrumFit.h:61
TSpectrumFit(void)
void SetTailParameters(Double_t tInit, Bool_t fixT, Double_t bInit, Bool_t fixB, Double_t sInit, Bool_t fixS)
SETTER FUNCTION.
Double_t Ourpowl(Double_t a, Int_t pw)
power function
Advanced 1-dimentional spectra fitting functions.
Definition: TSpectrumFit.h:32
Double_t Derb(Int_t num_of_fitted_peaks, Double_t i, const Double_t *parameter, Double_t sigma, Double_t t, Double_t b)
AUXILIARY FUNCTION // // This function calculates derivative of peaks shape function (see manual) // ...
Double_t * GetAreas() const
Definition: TSpectrumFit.h:132
void GetBackgroundParameters(Double_t &a0, Double_t &a0Err, Double_t &a1, Double_t &a1Err, Double_t &a2, Double_t &a2Err)
GETTER FUNCTION.
Double_t * fArea
Definition: TSpectrumFit.h:50
Double_t fSigmaInit
Definition: TSpectrumFit.h:52
Double_t fA2Calc
Definition: TSpectrumFit.h:71
Double_t fA0Calc
Definition: TSpectrumFit.h:65
Double_t fTCalc
Definition: TSpectrumFit.h:56
void GetSigma(Double_t &sigma, Double_t &sigmaErr)
GETTER FUNCTION.
void FitAwmi(Double_t *source)
ONE-DIMENSIONAL FIT FUNCTION ALGORITHM WITHOUT MATRIX INVERSION This function fits the source spectru...
void FitStiefel(Double_t *source)
ONE-DIMENSIONAL FIT FUNCTION ALGORITHM WITH MATRIX INVERSION (STIEFEL-HESTENS METHOD) This function f...
Double_t * fAmpCalc
Definition: TSpectrumFit.h:48
Double_t fSigmaErr
Definition: TSpectrumFit.h:54
Bool_t fFixT
Definition: TSpectrumFit.h:76
Bool_t fFixA1
Definition: TSpectrumFit.h:80
Bool_t fFixA2
Definition: TSpectrumFit.h:81