Logo ROOT  
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 #include "TNamed.h"
15 
16 class TH1;
17 
18 class TSpectrumFit : public TNamed {
19 protected:
20  Int_t fNPeaks; ///< number of peaks present in fit, input parameter, it should be > 0
21  Int_t fNumberIterations; ///< number of iterations in fitting procedure, input parameter, it should be > 0
22  Int_t fXmin; ///< first fitted channel
23  Int_t fXmax; ///< last fitted channel
24  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
25  Int_t fAlphaOptim; ///< optimization of convergence algorithm, possible values kFitAlphaHalving, kFitAlphaOptimal
26  Int_t fPower; ///< possible values kFitPower2,4,6,8,10,12, for details see references. It applies only for Awmi fitting function.
27  Int_t fFitTaylor; ///< order of Taylor expansion, possible values kFitTaylorOrderFirst, kFitTaylorOrderSecond. It applies only for Awmi fitting function.
28  Double_t fAlpha; ///< convergence coefficient, input parameter, it should be positive number and <=1, for details see references
29  Double_t fChi; ///< here the fitting functions return resulting chi square
30  Double_t *fPositionInit; ///< [fNPeaks] array of initial values of peaks positions, input parameters
31  Double_t *fPositionCalc; ///< [fNPeaks] array of calculated values of fitted positions, output parameters
32  Double_t *fPositionErr; ///< [fNPeaks] array of position errors
33  Double_t *fAmpInit; ///< [fNPeaks] array of initial values of peaks amplitudes, input parameters
34  Double_t *fAmpCalc; ///< [fNPeaks] array of calculated values of fitted amplitudes, output parameters
35  Double_t *fAmpErr; ///< [fNPeaks] array of amplitude errors
36  Double_t *fArea; ///< [fNPeaks] array of calculated areas of peaks
37  Double_t *fAreaErr; ///< [fNPeaks] array of errors of peak areas
38  Double_t fSigmaInit; ///< initial value of sigma parameter
39  Double_t fSigmaCalc; ///< calculated value of sigma parameter
40  Double_t fSigmaErr; ///< error value of sigma parameter
41  Double_t fTInit; ///< initial value of t parameter (relative amplitude of tail), for details see html manual and references
42  Double_t fTCalc; ///< calculated value of t parameter
43  Double_t fTErr; ///< error value of t parameter
44  Double_t fBInit; ///< initial value of b parameter (slope), for details see html manual and references
45  Double_t fBCalc; ///< calculated value of b parameter
46  Double_t fBErr; ///< error value of b parameter
47  Double_t fSInit; ///< initial value of s parameter (relative amplitude of step), for details see html manual and references
48  Double_t fSCalc; ///< calculated value of s parameter
49  Double_t fSErr; ///< error value of s parameter
50  Double_t fA0Init; ///< initial value of background a0 parameter(backgroud is estimated as a0+a1*x+a2*x*x)
51  Double_t fA0Calc; ///< calculated value of background a0 parameter
52  Double_t fA0Err; ///< error value of background a0 parameter
53  Double_t fA1Init; ///< initial value of background a1 parameter(backgroud is estimated as a0+a1*x+a2*x*x)
54  Double_t fA1Calc; ///< calculated value of background a1 parameter
55  Double_t fA1Err; ///< error value of background a1 parameter
56  Double_t fA2Init; ///< initial value of background a2 parameter(backgroud is estimated as a0+a1*x+a2*x*x)
57  Double_t fA2Calc; ///< calculated value of background a2 parameter
58  Double_t fA2Err; ///< error value of background a2 parameter
59  Bool_t *fFixPosition; ///< [fNPeaks] array of logical values which allow to fix appropriate positions (not fit). However they are present in the estimated functional
60  Bool_t *fFixAmp; ///< [fNPeaks] array of logical values which allow to fix appropriate amplitudes (not fit). However they are present in the estimated functional
61  Bool_t fFixSigma; ///< logical value of sigma parameter, which allows to fix the parameter (not to fit).
62  Bool_t fFixT; ///< logical value of t parameter, which allows to fix the parameter (not to fit).
63  Bool_t fFixB; ///< logical value of b parameter, which allows to fix the parameter (not to fit).
64  Bool_t fFixS; ///< logical value of s parameter, which allows to fix the parameter (not to fit).
65  Bool_t fFixA0; ///< logical value of a0 parameter, which allows to fix the parameter (not to fit).
66  Bool_t fFixA1; ///< logical value of a1 parameter, which allows to fix the parameter (not to fit).
67  Bool_t fFixA2; ///< logical value of a2 parameter, which allows to fix the parameter (not to fit).
68 
69 public:
70  enum {
85  };
86  TSpectrumFit(void); //default constructor
87  TSpectrumFit(Int_t numberPeaks);
88  virtual ~TSpectrumFit();
89 
90  //auxiliary functions for 1. parameter fit functions
91 protected:
96  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);
98  Double_t Derdersigma(Int_t num_of_fitted_peaks,Double_t i,const Double_t* parameter,Double_t sigma);
105  Double_t Ders(Int_t num_of_fitted_peaks,Double_t i,const Double_t* parameter,Double_t sigma);
106  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);
107  Double_t Dert(Int_t num_of_fitted_peaks,Double_t i,const Double_t* parameter,Double_t sigma,Double_t b);
110  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);
111  void StiefelInversion(Double_t **a,Int_t rozmer);
112 
113 public:
114  void FitAwmi(Double_t *source);
115  void FitStiefel(Double_t *source);
116  Double_t *GetAmplitudes() const {return fAmpCalc;}
118  Double_t *GetAreas() const {return fArea;}
119  Double_t *GetAreasErrors() const {return fAreaErr;}
120  void GetBackgroundParameters(Double_t &a0, Double_t &a0Err, Double_t &a1, Double_t &a1Err, Double_t &a2, Double_t &a2Err);
121  Double_t GetChi() const {return fChi;}
124  void GetSigma(Double_t &sigma, Double_t &sigmaErr);
125  void GetTailParameters(Double_t &t, Double_t &tErr, Double_t &b, Double_t &bErr, Double_t &s, Double_t &sErr);
126  void SetBackgroundParameters(Double_t a0Init, Bool_t fixA0, Double_t a1Init, Bool_t fixA1, Double_t a2Init, Bool_t fixA2);
127  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);
128  void SetPeakParameters(Double_t sigma, Bool_t fixSigma, const Double_t *positionInit, const Bool_t *fixPosition, const Double_t *ampInit, const Bool_t *fixAmp);
129  void SetTailParameters(Double_t tInit, Bool_t fixT, Double_t bInit, Bool_t fixB, Double_t sInit, Bool_t fixS);
130 
131  ClassDef(TSpectrumFit,1) //Spectrum Fitter using algorithm without matrix inversion and conjugate gradient method for symmetrical matrices (Stiefel-Hestens method)
132 };
133 
134 #endif
135 
TSpectrumFit::fFixB
Bool_t fFixB
logical value of b parameter, which allows to fix the parameter (not to fit).
Definition: TSpectrumFit.h:63
TSpectrumFit::fBInit
Double_t fBInit
initial value of b parameter (slope), for details see html manual and references
Definition: TSpectrumFit.h:44
TSpectrumFit::fBErr
Double_t fBErr
error value of b parameter
Definition: TSpectrumFit.h:46
TSpectrumFit::GetPositions
Double_t * GetPositions() const
Definition: TSpectrumFit.h:122
TSpectrumFit::fA2Err
Double_t fA2Err
error value of background a2 parameter
Definition: TSpectrumFit.h:58
TSpectrumFit::fBCalc
Double_t fBCalc
calculated value of b parameter
Definition: TSpectrumFit.h:45
TSpectrumFit::fTErr
Double_t fTErr
error value of t parameter
Definition: TSpectrumFit.h:43
TSpectrumFit::kFitPower10
@ kFitPower10
Definition: TSpectrumFit.h:80
TSpectrumFit::fA0Calc
Double_t fA0Calc
calculated value of background a0 parameter
Definition: TSpectrumFit.h:51
TSpectrumFit::Erfc
Double_t Erfc(Double_t x)
Definition: TSpectrumFit.cxx:177
TSpectrumFit::fSigmaErr
Double_t fSigmaErr
error value of sigma parameter
Definition: TSpectrumFit.h:40
TSpectrumFit::kFitTaylorOrderFirst
@ kFitTaylorOrderFirst
Definition: TSpectrumFit.h:82
TSpectrumFit::Shape
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)
This function calculates peaks shape function (see manual) Function parameters:
Definition: TSpectrumFit.cxx:516
TSpectrumFit::SetFitParameters
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)
This function sets the following fitting parameters:
Definition: TSpectrumFit.cxx:2608
TSpectrumFit::fA2Init
Double_t fA2Init
initial value of background a2 parameter(backgroud is estimated as a0+a1*x+a2*x*x)
Definition: TSpectrumFit.h:56
TSpectrumFit::SetBackgroundParameters
void SetBackgroundParameters(Double_t a0Init, Bool_t fixA0, Double_t a1Init, Bool_t fixA1, Double_t a2Init, Bool_t fixA2)
This function sets the following fitting parameters of background:
Definition: TSpectrumFit.cxx:2691
TSpectrumFit::Deri0
Double_t Deri0(Double_t i, Double_t amp, Double_t i0, Double_t sigma, Double_t t, Double_t s, Double_t b)
This function calculates derivative of peak shape function (see manual) according to peak position.
Definition: TSpectrumFit.cxx:266
xmax
float xmax
Definition: THbookFile.cxx:95
TNamed.h
TSpectrumFit::fA0Err
Double_t fA0Err
error value of background a0 parameter
Definition: TSpectrumFit.h:52
TSpectrumFit::kFitPower4
@ kFitPower4
Definition: TSpectrumFit.h:77
TSpectrumFit::fA1Calc
Double_t fA1Calc
calculated value of background a1 parameter
Definition: TSpectrumFit.h:54
TSpectrumFit::fNumberIterations
Int_t fNumberIterations
number of iterations in fitting procedure, input parameter, it should be > 0
Definition: TSpectrumFit.h:21
TGeant4Unit::s
static constexpr double s
Definition: TGeant4SystemOfUnits.h:168
TSpectrumFit::Derfc
Double_t Derfc(Double_t x)
This function calculates derivative of error function of x.
Definition: TSpectrumFit.cxx:200
TSpectrumFit::Derderi0
Double_t Derderi0(Double_t i, Double_t amp, Double_t i0, Double_t sigma)
This function calculates second derivative of peak shape function (see manual) according to peak posi...
Definition: TSpectrumFit.cxx:303
Int_t
int Int_t
Definition: RtypesCore.h:45
TSpectrumFit::Derpa
Double_t Derpa(Double_t sigma, Double_t t, Double_t b)
This function calculates derivative of the area of peak according to its amplitude.
Definition: TSpectrumFit.cxx:592
TSpectrumFit::SetTailParameters
void SetTailParameters(Double_t tInit, Bool_t fixT, Double_t bInit, Bool_t fixB, Double_t sInit, Bool_t fixS)
This function sets the following fitting parameters of tails of peaks.
Definition: TSpectrumFit.cxx:2710
TSpectrumFit::Ourpowl
Double_t Ourpowl(Double_t a, Int_t pw)
Power function.
Definition: TSpectrumFit.cxx:674
x
Double_t x[n]
Definition: legend1.C:17
TSpectrumFit::GetAreas
Double_t * GetAreas() const
Definition: TSpectrumFit.h:118
TSpectrumFit::fSigmaCalc
Double_t fSigmaCalc
calculated value of sigma parameter
Definition: TSpectrumFit.h:39
TSpectrumFit::Derdersigma
Double_t Derdersigma(Int_t num_of_fitted_peaks, Double_t i, const Double_t *parameter, Double_t sigma)
This function calculates second derivative of peaks shape function (see manual) according to sigma of...
Definition: TSpectrumFit.cxx:376
TSpectrumFit::Area
Double_t Area(Double_t a, Double_t sigma, Double_t t, Double_t b)
This function calculates area of a peak Function parameters:
Definition: TSpectrumFit.cxx:569
TSpectrumFit::Derb
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)
This function calculates derivative of peaks shape function (see manual) according to slope b.
Definition: TSpectrumFit.cxx:463
TSpectrumFit::kFitOptimChiCounts
@ kFitOptimChiCounts
Definition: TSpectrumFit.h:71
Bool_t
bool Bool_t
Definition: RtypesCore.h:63
TSpectrumFit::fSErr
Double_t fSErr
error value of s parameter
Definition: TSpectrumFit.h:49
b
#define b(i)
Definition: RSha256.hxx:118
TSpectrumFit::fFixAmp
Bool_t * fFixAmp
[fNPeaks] array of logical values which allow to fix appropriate amplitudes (not fit)....
Definition: TSpectrumFit.h:60
TSpectrumFit::kFitOptimChiFuncValues
@ kFitOptimChiFuncValues
Definition: TSpectrumFit.h:72
bool
TSpectrumFit::fAmpErr
Double_t * fAmpErr
[fNPeaks] array of amplitude errors
Definition: TSpectrumFit.h:35
TSpectrumFit::fA2Calc
Double_t fA2Calc
calculated value of background a2 parameter
Definition: TSpectrumFit.h:57
TSpectrumFit::kFitOptimMaxLikelihood
@ kFitOptimMaxLikelihood
Definition: TSpectrumFit.h:73
TSpectrumFit::GetAreasErrors
Double_t * GetAreasErrors() const
Definition: TSpectrumFit.h:119
TSpectrumFit::fArea
Double_t * fArea
[fNPeaks] array of calculated areas of peaks
Definition: TSpectrumFit.h:36
TSpectrumFit::fSigmaInit
Double_t fSigmaInit
initial value of sigma parameter
Definition: TSpectrumFit.h:38
TSpectrumFit::fAmpCalc
Double_t * fAmpCalc
[fNPeaks] array of calculated values of fitted amplitudes, output parameters
Definition: TSpectrumFit.h:34
TSpectrumFit::kFitAlphaHalving
@ kFitAlphaHalving
Definition: TSpectrumFit.h:74
TSpectrumFit::fFitTaylor
Int_t fFitTaylor
order of Taylor expansion, possible values kFitTaylorOrderFirst, kFitTaylorOrderSecond....
Definition: TSpectrumFit.h:27
TSpectrumFit::GetTailParameters
void GetTailParameters(Double_t &t, Double_t &tErr, Double_t &b, Double_t &bErr, Double_t &s, Double_t &sErr)
This function gets the tail parameters and their errors.
Definition: TSpectrumFit.cxx:2760
TSpectrumFit::fPower
Int_t fPower
possible values kFitPower2,4,6,8,10,12, for details see references. It applies only for Awmi fitting ...
Definition: TSpectrumFit.h:26
TSpectrumFit::GetAmplitudesErrors
Double_t * GetAmplitudesErrors() const
Definition: TSpectrumFit.h:117
TSpectrumFit::fFixA0
Bool_t fFixA0
logical value of a0 parameter, which allows to fix the parameter (not to fit).
Definition: TSpectrumFit.h:65
TSpectrumFit::Deramp
Double_t Deramp(Double_t i, Double_t i0, Double_t sigma, Double_t t, Double_t s, Double_t b)
This function calculates derivative of peak shape function (see manual) according to amplitude of pea...
Definition: TSpectrumFit.cxx:229
TSpectrumFit::fFixA1
Bool_t fFixA1
logical value of a1 parameter, which allows to fix the parameter (not to fit).
Definition: TSpectrumFit.h:66
TSpectrumFit::fFixA2
Bool_t fFixA2
logical value of a2 parameter, which allows to fix the parameter (not to fit).
Definition: TSpectrumFit.h:67
TSpectrumFit::GetAmplitudes
Double_t * GetAmplitudes() const
Definition: TSpectrumFit.h:116
TSpectrumFit::fTCalc
Double_t fTCalc
calculated value of t parameter
Definition: TSpectrumFit.h:42
TSpectrumFit::GetSigma
void GetSigma(Double_t &sigma, Double_t &sigmaErr)
This function gets the sigma parameter and its error.
Definition: TSpectrumFit.cxx:2725
xmin
float xmin
Definition: THbookFile.cxx:95
TSpectrumFit::fFixS
Bool_t fFixS
logical value of s parameter, which allows to fix the parameter (not to fit).
Definition: TSpectrumFit.h:64
TSpectrumFit::fNPeaks
Int_t fNPeaks
number of peaks present in fit, input parameter, it should be > 0
Definition: TSpectrumFit.h:20
TSpectrumFit::kFitPower6
@ kFitPower6
Definition: TSpectrumFit.h:78
TSpectrumFit::fPositionCalc
Double_t * fPositionCalc
[fNPeaks] array of calculated values of fitted positions, output parameters
Definition: TSpectrumFit.h:31
a
auto * a
Definition: textangle.C:12
TSpectrumFit::Ders
Double_t Ders(Int_t num_of_fitted_peaks, Double_t i, const Double_t *parameter, Double_t sigma)
This function calculates derivative of peaks shape function (see manual) according to relative amplit...
Definition: TSpectrumFit.cxx:437
TNamed
Definition: TNamed.h:29
TSpectrumFit::fAmpInit
Double_t * fAmpInit
[fNPeaks] array of initial values of peaks amplitudes, input parameters
Definition: TSpectrumFit.h:33
TSpectrumFit::fAlpha
Double_t fAlpha
convergence coefficient, input parameter, it should be positive number and <=1, for details see refer...
Definition: TSpectrumFit.h:28
TSpectrumFit::fAreaErr
Double_t * fAreaErr
[fNPeaks] array of errors of peak areas
Definition: TSpectrumFit.h:37
TSpectrumFit::fSInit
Double_t fSInit
initial value of s parameter (relative amplitude of step), for details see html manual and references
Definition: TSpectrumFit.h:47
TSpectrumFit::Derpb
Double_t Derpb(Double_t a, Double_t sigma, Double_t t, Double_t b)
This function calculates derivative of the area of peak according to b parameter.
Definition: TSpectrumFit.cxx:658
TSpectrumFit::fFixPosition
Bool_t * fFixPosition
[fNPeaks] array of logical values which allow to fix appropriate positions (not fit)....
Definition: TSpectrumFit.h:59
TSpectrumFit::GetChi
Double_t GetChi() const
Definition: TSpectrumFit.h:121
TSpectrumFit::fFixSigma
Bool_t fFixSigma
logical value of sigma parameter, which allows to fix the parameter (not to fit).
Definition: TSpectrumFit.h:61
TSpectrumFit::SetPeakParameters
void SetPeakParameters(Double_t sigma, Bool_t fixSigma, const Double_t *positionInit, const Bool_t *fixPosition, const Double_t *ampInit, const Bool_t *fixAmp)
This function sets the following fitting parameters of peaks:
Definition: TSpectrumFit.cxx:2656
TSpectrumFit::fXmax
Int_t fXmax
last fitted channel
Definition: TSpectrumFit.h:23
TSpectrumFit::fTInit
Double_t fTInit
initial value of t parameter (relative amplitude of tail), for details see html manual and references
Definition: TSpectrumFit.h:41
TSpectrumFit::fPositionErr
Double_t * fPositionErr
[fNPeaks] array of position errors
Definition: TSpectrumFit.h:32
TSpectrumFit::GetBackgroundParameters
void GetBackgroundParameters(Double_t &a0, Double_t &a0Err, Double_t &a1, Double_t &a1Err, Double_t &a2, Double_t &a2Err)
This function gets the background parameters and their errors.
Definition: TSpectrumFit.cxx:2740
TSpectrumFit::kFitPower8
@ kFitPower8
Definition: TSpectrumFit.h:79
TSpectrumFit::Dert
Double_t Dert(Int_t num_of_fitted_peaks, Double_t i, const Double_t *parameter, Double_t sigma, Double_t b)
This function calculates derivative of peaks shape function (see manual) according to relative amplit...
Definition: TSpectrumFit.cxx:409
sigma
const Double_t sigma
Definition: h1analysisProxy.h:11
TSpectrumFit::Dersigma
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)
This function calculates derivative of peaks shape function (see manual) according to sigma of peaks.
Definition: TSpectrumFit.cxx:331
TSpectrumFit::Derpsigma
Double_t Derpsigma(Double_t a, Double_t t, Double_t b)
This function calculates derivative of the area of peak according to sigma of peaks.
Definition: TSpectrumFit.cxx:614
TSpectrumFit::FitAwmi
void FitAwmi(Double_t *source)
This function fits the source spectrum.
Definition: TSpectrumFit.cxx:820
TSpectrumFit::FitStiefel
void FitStiefel(Double_t *source)
This function fits the source spectrum.
Definition: TSpectrumFit.cxx:1857
TSpectrumFit::fXmin
Int_t fXmin
first fitted channel
Definition: TSpectrumFit.h:22
Double_t
double Double_t
Definition: RtypesCore.h:59
TSpectrumFit::kFitAlphaOptimal
@ kFitAlphaOptimal
Definition: TSpectrumFit.h:75
TSpectrumFit::fFixT
Bool_t fFixT
logical value of t parameter, which allows to fix the parameter (not to fit).
Definition: TSpectrumFit.h:62
TSpectrumFit::StiefelInversion
void StiefelInversion(Double_t **a, Int_t rozmer)
This function calculates solution of the system of linear equations.
Definition: TSpectrumFit.cxx:1721
TSpectrumFit::fA1Err
Double_t fA1Err
error value of background a1 parameter
Definition: TSpectrumFit.h:55
TSpectrumFit::fStatisticType
Int_t fStatisticType
type of statistics, possible values kFitOptimChiCounts (chi square statistics with counts as weightin...
Definition: TSpectrumFit.h:24
ClassDef
#define ClassDef(name, id)
Definition: Rtypes.h:325
TSpectrumFit::fA0Init
Double_t fA0Init
initial value of background a0 parameter(backgroud is estimated as a0+a1*x+a2*x*x)
Definition: TSpectrumFit.h:50
TH1
Definition: TH1.h:57
TSpectrumFit::~TSpectrumFit
virtual ~TSpectrumFit()
Destructor.
Definition: TSpectrumFit.cxx:160
TSpectrumFit::kFitPower12
@ kFitPower12
Definition: TSpectrumFit.h:81
TSpectrumFit
Advanced 1-dimensional spectra fitting functions.
Definition: TSpectrumFit.h:18
TSpectrumFit::GetPositionsErrors
Double_t * GetPositionsErrors() const
Definition: TSpectrumFit.h:123
TSpectrumFit::fAlphaOptim
Int_t fAlphaOptim
optimization of convergence algorithm, possible values kFitAlphaHalving, kFitAlphaOptimal
Definition: TSpectrumFit.h:25
TSpectrumFit::Dera1
Double_t Dera1(Double_t i)
Derivative of background according to a1.
Definition: TSpectrumFit.cxx:492
TSpectrumFit::kFitTaylorOrderSecond
@ kFitTaylorOrderSecond
Definition: TSpectrumFit.h:83
TSpectrumFit::kFitPower2
@ kFitPower2
Definition: TSpectrumFit.h:76
TSpectrumFit::Dera2
Double_t Dera2(Double_t i)
Derivative of background according to a2.
Definition: TSpectrumFit.cxx:500
TSpectrumFit::Derpt
Double_t Derpt(Double_t a, Double_t sigma, Double_t b)
This function calculates derivative of the area of peak according to t parameter.
Definition: TSpectrumFit.cxx:636
TSpectrumFit::fA1Init
Double_t fA1Init
initial value of background a1 parameter(backgroud is estimated as a0+a1*x+a2*x*x)
Definition: TSpectrumFit.h:53
TSpectrumFit::fPositionInit
Double_t * fPositionInit
[fNPeaks] array of initial values of peaks positions, input parameters
Definition: TSpectrumFit.h:30
TSpectrumFit::TSpectrumFit
TSpectrumFit(void)
Default constructor.
Definition: TSpectrumFit.cxx:35
TSpectrumFit::kFitNumRegulCycles
@ kFitNumRegulCycles
Definition: TSpectrumFit.h:84
TSpectrumFit::fSCalc
Double_t fSCalc
calculated value of s parameter
Definition: TSpectrumFit.h:48
int
TSpectrumFit::fChi
Double_t fChi
here the fitting functions return resulting chi square
Definition: TSpectrumFit.h:29