Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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
16class TH1;
17
18class TSpectrumFit : public TNamed {
19protected:
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
69public:
70 enum {
85 };
86 TSpectrumFit(void); //default constructor
87 TSpectrumFit(Int_t numberPeaks);
88 ~TSpectrumFit() override;
89
90 //auxiliary functions for 1. parameter fit functions
91protected:
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
113public:
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;}
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 ClassDefOverride(TSpectrumFit,1) //Spectrum Fitter using algorithm without matrix inversion and conjugate gradient method for symmetrical matrices (Stiefel-Hestens method)
132};
133
134#endif
135
#define b(i)
Definition RSha256.hxx:100
#define a(i)
Definition RSha256.hxx:99
bool Bool_t
Definition RtypesCore.h:63
int Int_t
Definition RtypesCore.h:45
double Double_t
Definition RtypesCore.h:59
#define ClassDefOverride(name, id)
Definition Rtypes.h:341
float xmin
float xmax
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:59
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
Advanced 1-dimensional spectra fitting functions.
Double_t GetChi() const
Bool_t fFixSigma
logical value of sigma parameter, which allows to fix the parameter (not to fit).
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...
Double_t Area(Double_t a, Double_t sigma, Double_t t, Double_t b)
This function calculates area of a peak Function parameters:
Double_t * fAmpErr
[fNPeaks] array of amplitude errors
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:
Double_t fAlpha
convergence coefficient, input parameter, it should be positive number and <=1, for details see refer...
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.
Double_t Ourpowl(Double_t a, Int_t pw)
Power function.
Int_t fAlphaOptim
optimization of convergence algorithm, possible values kFitAlphaHalving, kFitAlphaOptimal
Double_t fA1Init
initial value of background a1 parameter(backgroud is estimated as a0+a1*x+a2*x*x)
void FitStiefel(Double_t *source)
This function fits the source spectrum.
Double_t * GetAmplitudesErrors() const
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...
Double_t * fPositionCalc
[fNPeaks] array of calculated values of fitted positions, output parameters
Double_t fSInit
initial value of s parameter (relative amplitude of step), for details see html manual and references
Double_t fA1Calc
calculated value of background a1 parameter
Double_t fSigmaErr
error value of sigma parameter
Bool_t * fFixAmp
[fNPeaks] array of logical values which allow to fix appropriate amplitudes (not fit)....
Double_t * fAmpCalc
[fNPeaks] array of calculated values of fitted amplitudes, output parameters
Bool_t fFixB
logical value of b parameter, which allows to fix the parameter (not to fit).
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.
Double_t fBErr
error value of b parameter
void FitAwmi(Double_t *source)
This function fits the source spectrum.
Int_t fNumberIterations
number of iterations in fitting procedure, input parameter, it should be > 0
Double_t * fAmpInit
[fNPeaks] array of initial values of peaks amplitudes, input parameters
Double_t fA1Err
error value of background a1 parameter
Double_t fA2Err
error value of background a2 parameter
Bool_t * fFixPosition
[fNPeaks] array of logical values which allow to fix appropriate positions (not fit)....
Double_t fA0Calc
calculated value of background a0 parameter
Int_t fNPeaks
number of peaks present in fit, input parameter, it should be > 0
Double_t fBCalc
calculated value of b parameter
Bool_t fFixA2
logical value of a2 parameter, which allows to fix the parameter (not to fit).
Double_t Erfc(Double_t x)
TSpectrumFit(void)
Default constructor.
Double_t * GetAreasErrors() const
~TSpectrumFit() override
Destructor.
Double_t fA2Init
initial value of background a2 parameter(backgroud is estimated as a0+a1*x+a2*x*x)
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.
Bool_t fFixA1
logical value of a1 parameter, which allows to fix the parameter (not to fit).
Double_t fA2Calc
calculated value of background a2 parameter
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.
Double_t * fPositionErr
[fNPeaks] array of position errors
Double_t fTInit
initial value of t parameter (relative amplitude of tail), for details see html manual and references
Bool_t fFixT
logical value of t parameter, which allows to fix the parameter (not to fit).
void GetSigma(Double_t &sigma, Double_t &sigmaErr)
This function gets the sigma parameter and its error.
Int_t fFitTaylor
order of Taylor expansion, possible values kFitTaylorOrderFirst, kFitTaylorOrderSecond....
Double_t Dera1(Double_t i)
Derivative of background according to a1.
Double_t fSigmaCalc
calculated value of sigma parameter
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.
Double_t fSErr
error value of s parameter
Bool_t fFixA0
logical value of a0 parameter, which allows to fix the parameter (not to fit).
Int_t fPower
possible values kFitPower2,4,6,8,10,12, for details see references. It applies only for Awmi fitting ...
Double_t fChi
here the fitting functions return resulting chi square
Double_t * GetAreas() const
Double_t * fArea
[fNPeaks] array of calculated areas of peaks
Double_t Dera2(Double_t i)
Derivative of background according to a2.
Double_t fBInit
initial value of b parameter (slope), for details see html manual and references
Double_t * fPositionInit
[fNPeaks] array of initial values of peaks positions, input parameters
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.
Double_t * GetAmplitudes() const
Bool_t fFixS
logical value of s parameter, which allows to fix the parameter (not to fit).
Double_t fA0Init
initial value of background a0 parameter(backgroud is estimated as a0+a1*x+a2*x*x)
Int_t fXmin
first fitted channel
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:
Double_t fSCalc
calculated value of s parameter
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.
Int_t fXmax
last fitted channel
Int_t fStatisticType
type of statistics, possible values kFitOptimChiCounts (chi square statistics with counts as weightin...
Double_t fTErr
error value of t parameter
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...
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.
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.
Double_t fSigmaInit
initial value of sigma parameter
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...
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...
Double_t fA0Err
error value of background a0 parameter
Double_t fTCalc
calculated value of t parameter
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:
Double_t * fAreaErr
[fNPeaks] array of errors of peak areas
Double_t * GetPositionsErrors() const
void StiefelInversion(Double_t **a, Int_t rozmer)
This function calculates solution of the system of linear equations.
Double_t Derfc(Double_t x)
This function calculates derivative of error function of x.
Double_t * GetPositions() const
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:
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.
const Double_t sigma
Double_t x[n]
Definition legend1.C:17