Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TSpectrum2Fit.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_TSpectrum2Fit
12#define ROOT_TSpectrum2Fit
13
14#include "TNamed.h"
15
16class TSpectrum2Fit : public TNamed {
17protected:
18 Int_t fNPeaks; ///< number of peaks present in fit, input parameter, it should be > 0
19 Int_t fNumberIterations; ///< number of iterations in fitting procedure, input parameter, it should be > 0
20 Int_t fXmin; ///< first fitted channel in x direction
21 Int_t fXmax; ///< last fitted channel in x direction
22 Int_t fYmin; ///< first fitted channel in y direction
23 Int_t fYmax; ///< last fitted channel in y direction
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 *fPositionInitX; ///< [fNPeaks] array of initial values of x positions of 2D peaks, input parameters
31 Double_t *fPositionCalcX; ///< [fNPeaks] array of calculated values of x positions of 2D peaks, output parameters
32 Double_t *fPositionErrX; ///< [fNPeaks] array of error values of x positions of 2D peaks, output parameters
33 Double_t *fPositionInitY; ///< [fNPeaks] array of initial values of y positions of 2D peaks, input parameters
34 Double_t *fPositionCalcY; ///< [fNPeaks] array of calculated values of y positions of 2D peaks, output parameters
35 Double_t *fPositionErrY; ///< [fNPeaks] array of error values of y positions of 2D peaks, output parameters
36 Double_t *fPositionInitX1; ///< [fNPeaks] array of initial x positions of 1D ridges, input parameters
37 Double_t *fPositionCalcX1; ///< [fNPeaks] array of calculated x positions of 1D ridges, output parameters
38 Double_t *fPositionErrX1; ///< [fNPeaks] array of x positions errors of 1D ridges, output parameters
39 Double_t *fPositionInitY1; ///< [fNPeaks] array of initial y positions of 1D ridges, input parameters
40 Double_t *fPositionCalcY1; ///< [fNPeaks] array of calculated y positions of 1D ridges, output parameters
41 Double_t *fPositionErrY1; ///< [fNPeaks] array of y positions errors of 1D ridges, output parameters
42 Double_t *fAmpInit; ///< [fNPeaks] array of initial values of amplitudes of 2D peaks, input parameters
43 Double_t *fAmpCalc; ///< [fNPeaks] array of calculated values of amplitudes of 2D peaks, output parameters
44 Double_t *fAmpErr; ///< [fNPeaks] array of amplitudes errors of 2D peaks, output parameters
45 Double_t *fAmpInitX1; ///< [fNPeaks] array of initial values of amplitudes of 1D ridges in x direction, input parameters
46 Double_t *fAmpCalcX1; ///< [fNPeaks] array of calculated values of amplitudes of 1D ridges in x direction, output parameters
47 Double_t *fAmpErrX1; ///< [fNPeaks] array of amplitudes errors of 1D ridges in x direction, output parameters
48 Double_t *fAmpInitY1; ///< [fNPeaks] array of initial values of amplitudes of 1D ridges in y direction, input parameters
49 Double_t *fAmpCalcY1; ///< [fNPeaks] array of calculated values of amplitudes of 1D ridges in y direction, output parameters
50 Double_t *fAmpErrY1; ///< [fNPeaks] array of amplitudes errors of 1D ridges in y direction, output parameters
51 Double_t *fVolume; ///< [fNPeaks] array of calculated volumes of 2D peaks, output parameters
52 Double_t *fVolumeErr; ///< [fNPeaks] array of volumes errors of 2D peaks, output parameters
53 Double_t fSigmaInitX; ///< initial value of sigma x parameter
54 Double_t fSigmaCalcX; ///< calculated value of sigma x parameter
55 Double_t fSigmaErrX; ///< error value of sigma x parameter
56 Double_t fSigmaInitY; ///< initial value of sigma y parameter
57 Double_t fSigmaCalcY; ///< calculated value of sigma y parameter
58 Double_t fSigmaErrY; ///< error value of sigma y parameter
59 Double_t fRoInit; ///< initial value of correlation coefficient
60 Double_t fRoCalc; ///< calculated value of correlation coefficient
61 Double_t fRoErr; ///< error value of correlation coefficient
62 Double_t fTxyInit; ///< initial value of t parameter for 2D peaks (relative amplitude of tail), for details see html manual and references
63 Double_t fTxyCalc; ///< calculated value of t parameter for 2D peaks
64 Double_t fTxyErr; ///< error value of t parameter for 2D peaks
65 Double_t fSxyInit; ///< initial value of s parameter for 2D peaks (relative amplitude of step), for details see html manual and references
66 Double_t fSxyCalc; ///< calculated value of s parameter for 2D peaks
67 Double_t fSxyErr; ///< error value of s parameter for 2D peaks
68 Double_t fTxInit; ///< initial value of t parameter for 1D ridges in x direction (relative amplitude of tail), for details see html manual and references
69 Double_t fTxCalc; ///< calculated value of t parameter for 1D ridges in x direction
70 Double_t fTxErr; ///< error value of t parameter for 1D ridges in x direction
71 Double_t fTyInit; ///< initial value of t parameter for 1D ridges in y direction (relative amplitude of tail), for details see html manual and references
72 Double_t fTyCalc; ///< calculated value of t parameter for 1D ridges in y direction
73 Double_t fTyErr; ///< error value of t parameter for 1D ridges in y direction
74 Double_t fSxInit; ///< initial value of s parameter for 1D ridges in x direction (relative amplitude of step), for details see html manual and references
75 Double_t fSxCalc; ///< calculated value of s parameter for 1D ridges in x direction
76 Double_t fSxErr; ///< error value of s parameter for 1D ridges in x direction
77 Double_t fSyInit; ///< initial value of s parameter for 1D ridges in y direction (relative amplitude of step), for details see html manual and references
78 Double_t fSyCalc; ///< calculated value of s parameter for 1D ridges in y direction
79 Double_t fSyErr; ///< error value of s parameter for 1D ridges in y direction
80 Double_t fBxInit; ///< initial value of b parameter for 1D ridges in x direction (slope), for details see html manual and references
81 Double_t fBxCalc; ///< calculated value of b parameter for 1D ridges in x direction
82 Double_t fBxErr; ///< error value of b parameter for 1D ridges in x direction
83 Double_t fByInit; ///< initial value of b parameter for 1D ridges in y direction (slope), for details see html manual and references
84 Double_t fByCalc; ///< calculated value of b parameter for 1D ridges in y direction
85 Double_t fByErr; ///< error value of b parameter for 1D ridges in y direction
86 Double_t fA0Init; ///< initial value of background a0 parameter(backgroud is estimated as a0+ax*x+ay*y)
87 Double_t fA0Calc; ///< calculated value of background a0 parameter
88 Double_t fA0Err; ///< error value of background a0 parameter
89 Double_t fAxInit; ///< initial value of background ax parameter(backgroud is estimated as a0+ax*x+ay*y)
90 Double_t fAxCalc; ///< calculated value of background ax parameter
91 Double_t fAxErr; ///< error value of background ax parameter
92 Double_t fAyInit; ///< initial value of background ay parameter(backgroud is estimated as a0+ax*x+ay*y)
93 Double_t fAyCalc; ///< calculated value of background ay parameter
94 Double_t fAyErr; ///< error value of background ay parameter
95 Bool_t *fFixPositionX; ///< [fNPeaks] array of logical values which allow to fix appropriate x positions of 2D peaks (not fit). However they are present in the estimated functional
96 Bool_t *fFixPositionY; ///< [fNPeaks] array of logical values which allow to fix appropriate y positions of 2D peaks (not fit). However they are present in the estimated functional
97 Bool_t *fFixPositionX1; ///< [fNPeaks] array of logical values which allow to fix appropriate x positions of 1D ridges (not fit). However they are present in the estimated functional
98 Bool_t *fFixPositionY1; ///< [fNPeaks] array of logical values which allow to fix appropriate y positions of 1D ridges (not fit). However they are present in the estimated functional
99 Bool_t *fFixAmp; ///< [fNPeaks] array of logical values which allow to fix appropriate amplitudes of 2D peaks (not fit). However they are present in the estimated functional
100 Bool_t *fFixAmpX1; ///< [fNPeaks] array of logical values which allow to fix appropriate amplitudes of 1D ridges in x direction (not fit). However they are present in the estimated functional
101 Bool_t *fFixAmpY1; ///< [fNPeaks] array of logical values which allow to fix appropriate amplitudes of 1D ridges in y direction (not fit). However they are present in the estimated functional
102 Bool_t fFixSigmaX; ///< logical value of sigma x parameter, which allows to fix the parameter (not to fit).
103 Bool_t fFixSigmaY; ///< logical value of sigma y parameter, which allows to fix the parameter (not to fit).
104 Bool_t fFixRo; ///< logical value of correlation coefficient, which allows to fix the parameter (not to fit).
105 Bool_t fFixTxy; ///< logical value of t parameter for 2D peaks, which allows to fix the parameter (not to fit).
106 Bool_t fFixSxy; ///< logical value of s parameter for 2D peaks, which allows to fix the parameter (not to fit).
107 Bool_t fFixTx; ///< logical value of t parameter for 1D ridges in x direction, which allows to fix the parameter (not to fit).
108 Bool_t fFixTy; ///< logical value of t parameter for 1D ridges in y direction, which allows to fix the parameter (not to fit).
109 Bool_t fFixSx; ///< logical value of s parameter for 1D ridges in x direction, which allows to fix the parameter (not to fit).
110 Bool_t fFixSy; ///< logical value of s parameter for 1D ridges in y direction, which allows to fix the parameter (not to fit).
111 Bool_t fFixBx; ///< logical value of b parameter for 1D ridges in x direction, which allows to fix the parameter (not to fit).
112 Bool_t fFixBy; ///< logical value of b parameter for 1D ridges in y direction, which allows to fix the parameter (not to fit).
113 Bool_t fFixA0; ///< logical value of a0 parameter, which allows to fix the parameter (not to fit).
114 Bool_t fFixAx; ///< logical value of ax parameter, which allows to fix the parameter (not to fit).
115 Bool_t fFixAy; ///< logical value of ay parameter, which allows to fix the parameter (not to fit).
116public:
117 enum {
132 };
133 TSpectrum2Fit(void); //default constructor
134 TSpectrum2Fit(Int_t numberPeaks);
135 ~TSpectrum2Fit() override;
136 //auxiliary functions for 2. parameter fit functions
137protected:
140 Double_t Derbx(Int_t numOfFittedPeaks,Double_t x,Double_t y,const Double_t *parameter,Double_t sigmax,Double_t sigmay,Double_t txy,Double_t tx,Double_t bx,Double_t by);
141 Double_t Derby(Int_t numOfFittedPeaks,Double_t x,Double_t y,const Double_t *parameter,Double_t sigmax,Double_t sigmay,Double_t txy,Double_t ty,Double_t bx,Double_t by);
145 Double_t Derdersigmax(Int_t numOfFittedPeaks,Double_t x,Double_t y,const Double_t *parameter,Double_t sigmax,Double_t sigmay,Double_t ro);
146 Double_t Derdersigmay(Int_t numOfFittedPeaks,Double_t x,Double_t y,const Double_t *parameter,Double_t sigmax,Double_t sigmay,Double_t ro);
155 Double_t Derro(Int_t numOfFittedPeaks,Double_t x,Double_t y,const Double_t *parameter,Double_t sx,Double_t sy,Double_t r);
156 Double_t Dersigmax(Int_t numOfFittedPeaks,Double_t x,Double_t y,const Double_t *parameter,Double_t sigmax,Double_t sigmay,Double_t ro,Double_t txy,Double_t sxy,Double_t tx,Double_t sx,Double_t bx,Double_t by);
157 Double_t Dersigmay(Int_t numOfFittedPeaks,Double_t x,Double_t y,const Double_t *parameter,Double_t sigmax,Double_t sigmay,Double_t ro,Double_t txy,Double_t sxy,Double_t ty,Double_t sy,Double_t bx,Double_t by);
158 Double_t Dersx(Int_t numOfFittedPeaks,Double_t x,const Double_t *parameter,Double_t sigmax);
159 Double_t Dersxy(Int_t numOfFittedPeaks,Double_t x,Double_t y,const Double_t *parameter,Double_t sigmax,Double_t sigmay);
160 Double_t Dersy(Int_t numOfFittedPeaks,Double_t x,const Double_t *parameter,Double_t sigmax);
161 Double_t Dertx(Int_t numOfFittedPeaks,Double_t x,const Double_t *parameter,Double_t sigmax,Double_t bx);
162 Double_t Dertxy(Int_t numOfFittedPeaks,Double_t x,Double_t y,const Double_t *parameter,Double_t sigmax,Double_t sigmay,Double_t bx,Double_t by);
163 Double_t Derty(Int_t numOfFittedPeaks,Double_t x,const Double_t *parameter,Double_t sigmax,Double_t bx);
166 Double_t Shape2(Int_t numOfFittedPeaks,Double_t x,Double_t y,const Double_t *parameter,Double_t sigmax,Double_t sigmay,Double_t ro,Double_t a0,Double_t ax,Double_t ay,Double_t txy,Double_t sxy,Double_t tx,Double_t ty,Double_t sx,Double_t sy,Double_t bx,Double_t by);
169
170public:
171 void FitAwmi(Double_t **source);
172 void FitStiefel(Double_t **source);
173 void GetAmplitudes(Double_t *amplitudes, Double_t *amplitudesX1, Double_t *amplitudesY1);
174 void GetAmplitudeErrors(Double_t *amplitudeErrors, Double_t *amplitudeErrorsX1, Double_t *amplitudeErrorsY1);
175 void GetBackgroundParameters(Double_t &a0, Double_t &a0Err, Double_t &ax, Double_t &axErr, Double_t &ay, Double_t &ayErr);
176 Double_t GetChi() const {return fChi;}
177 void GetPositions(Double_t *positionsX, Double_t *positionsY, Double_t *positionsX1, Double_t *positionsY1);
178 void GetPositionErrors(Double_t *positionErrorsX, Double_t *positionErrorsY, Double_t *positionErrorsX1, Double_t *positionErrorsY1);
179 void GetRo(Double_t &ro, Double_t &roErr);
180 void GetSigmaX(Double_t &sigmaX, Double_t &sigmaErrX);
181 void GetSigmaY(Double_t &sigmaY, Double_t &sigmaErrY);
182 void GetTailParameters(Double_t &txy, Double_t &txyErr, Double_t &tx, Double_t &txErr, Double_t &ty, Double_t &tyErr, Double_t &bx, Double_t &bxErr, Double_t &by, Double_t &byErr, Double_t &sxy, Double_t &sxyErr, Double_t &sx, Double_t &sxErr, Double_t &sy, Double_t &syErr);
183 void GetVolumes(Double_t *volumes);
184 void GetVolumeErrors(Double_t *volumeErrors);
185 void SetBackgroundParameters(Double_t a0Init, Bool_t fixA0, Double_t axInit, Bool_t fixAx, Double_t ayInit, Bool_t fixAy);
186 void SetFitParameters(Int_t xmin,Int_t xmax,Int_t ymin,Int_t ymax, Int_t numberIterations, Double_t alpha, Int_t statisticType, Int_t alphaOptim, Int_t power, Int_t fitTaylor);
187 void SetPeakParameters(Double_t sigmaX, Bool_t fixSigmaX, Double_t sigmaY, Bool_t fixSigmaY, Double_t ro, Bool_t fixRo, const Double_t *positionInitX, const Bool_t *fixPositionX, const Double_t *positionInitY, const Bool_t *fixPositionY, const Double_t *positionInitX1, const Bool_t *fixPositionX1, const Double_t *positionInitY1, const Bool_t *fixPositionY1, const Double_t *ampInit, const Bool_t *fixAmp, const Double_t *ampInitX1, const Bool_t *fixAmpX1, const Double_t *ampInitY1, const Bool_t *fixAmpY1);
188 void SetTailParameters(Double_t tInitXY, Bool_t fixTxy, Double_t tInitX, Bool_t fixTx, Double_t tInitY, Bool_t fixTy, Double_t bInitX, Bool_t fixBx, Double_t bInitY, Bool_t fixBy, Double_t sInitXY, Bool_t fixSxy, Double_t sInitX, Bool_t fixSx, Double_t sInitY, Bool_t fixSy);
189
190 ClassDefOverride(TSpectrum2Fit,1) //Spectrum2 Fitter using algorithm without matrix inversion and conjugate gradient method for symmetrical matrices (Stiefel-Hestens method)
191};
192
193#endif
194
#define a(i)
Definition RSha256.hxx:99
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
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
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
float xmin
float ymin
float xmax
float ymax
The TNamed class is the base class for all named ROOT classes.
Definition TNamed.h:29
Advanced 2-dimensional spectra fitting functions.
Double_t Derbx(Int_t numOfFittedPeaks, Double_t x, Double_t y, const Double_t *parameter, Double_t sigmax, Double_t sigmay, Double_t txy, Double_t tx, Double_t bx, Double_t by)
This function calculates derivative of peaks shape function (see manual) according to slope bx.
Double_t Dersxy(Int_t numOfFittedPeaks, Double_t x, Double_t y, const Double_t *parameter, Double_t sigmax, Double_t sigmay)
This function calculates derivative of peaks shape function (see manual) according to relative amplit...
Bool_t fFixBx
logical value of b parameter for 1D ridges in x direction, which allows to fix the parameter (not to ...
Bool_t * fFixPositionY
[fNPeaks] array of logical values which allow to fix appropriate y positions of 2D peaks (not fit)....
Double_t fSigmaErrY
error value of sigma y parameter
Bool_t * fFixAmpY1
[fNPeaks] array of logical values which allow to fix appropriate amplitudes of 1D ridges in y directi...
Double_t * fAmpCalcX1
[fNPeaks] array of calculated values of amplitudes of 1D ridges in x direction, output parameters
Double_t Volume(Double_t a, Double_t sx, Double_t sy, Double_t ro)
This function calculates volume of a peak.
Bool_t fFixSxy
logical value of s parameter for 2D peaks, which allows to fix the parameter (not to fit).
Double_t fTxErr
error value of t parameter for 1D ridges in x direction
Double_t Deri01(Double_t x, Double_t ax, Double_t x0, Double_t sigmax, Double_t tx, Double_t sx, Double_t bx)
This function calculates derivative of 2D peaks shape function (see manual) according to x position o...
Double_t * fAmpErrY1
[fNPeaks] array of amplitudes errors of 1D ridges in y direction, output parameters
Double_t fBxInit
initial value of b parameter for 1D ridges in x direction (slope), for details see html manual and re...
Double_t Shape2(Int_t numOfFittedPeaks, Double_t x, Double_t y, const Double_t *parameter, Double_t sigmax, Double_t sigmay, Double_t ro, Double_t a0, Double_t ax, Double_t ay, Double_t txy, Double_t sxy, Double_t tx, Double_t ty, Double_t sx, Double_t sy, Double_t bx, Double_t by)
This function calculates 2D peaks shape function (see manual)
Double_t fA0Err
error value of background a0 parameter
Double_t GetChi() const
void SetPeakParameters(Double_t sigmaX, Bool_t fixSigmaX, Double_t sigmaY, Bool_t fixSigmaY, Double_t ro, Bool_t fixRo, const Double_t *positionInitX, const Bool_t *fixPositionX, const Double_t *positionInitY, const Bool_t *fixPositionY, const Double_t *positionInitX1, const Bool_t *fixPositionX1, const Double_t *positionInitY1, const Bool_t *fixPositionY1, const Double_t *ampInit, const Bool_t *fixAmp, const Double_t *ampInitX1, const Bool_t *fixAmpX1, const Double_t *ampInitY1, const Bool_t *fixAmpY1)
This function sets the following fitting parameters of peaks:
Double_t fSigmaCalcX
calculated value of sigma x parameter
Bool_t * fFixAmpX1
[fNPeaks] array of logical values which allow to fix appropriate amplitudes of 1D ridges in x directi...
Double_t Derderi02(Double_t x, Double_t y, Double_t a, Double_t x0, Double_t y0, Double_t sigmax, Double_t sigmay, Double_t ro)
This function calculates second derivative of 2D peaks shape function (see manual) according to x pos...
Bool_t fFixRo
logical value of correlation coefficient, which allows to fix the parameter (not to fit).
Double_t * fAmpErr
[fNPeaks] array of amplitudes errors of 2D peaks, output parameters
Bool_t fFixSx
logical value of s parameter for 1D ridges in x direction, which allows to fix the parameter (not to ...
Double_t * fAmpInit
[fNPeaks] array of initial values of amplitudes of 2D peaks, input parameters
void GetTailParameters(Double_t &txy, Double_t &txyErr, Double_t &tx, Double_t &txErr, Double_t &ty, Double_t &tyErr, Double_t &bx, Double_t &bxErr, Double_t &by, Double_t &byErr, Double_t &sxy, Double_t &sxyErr, Double_t &sx, Double_t &sxErr, Double_t &sy, Double_t &syErr)
This function gets the tail parameters and their errors.
Double_t * fVolume
[fNPeaks] array of calculated volumes of 2D peaks, output parameters
Double_t Derdersigmax(Int_t numOfFittedPeaks, Double_t x, Double_t y, const Double_t *parameter, Double_t sigmax, Double_t sigmay, Double_t ro)
This function calculates second derivative of peaks shape function (see manual) according to sigmax o...
Double_t * fPositionCalcY
[fNPeaks] array of calculated values of y positions of 2D peaks, output parameters
Int_t fFitTaylor
order of Taylor expansion, possible values kFitTaylorOrderFirst, kFitTaylorOrderSecond....
Double_t fSxyErr
error value of s parameter for 2D peaks
Double_t * fPositionErrY1
[fNPeaks] array of y positions errors of 1D ridges, output parameters
Double_t fBxCalc
calculated value of b parameter for 1D ridges in x direction
Double_t fRoErr
error value of correlation coefficient
Double_t fTxCalc
calculated value of t parameter for 1D ridges in x direction
Double_t fRoCalc
calculated value of correlation coefficient
Bool_t fFixBy
logical value of b parameter for 1D ridges in y direction, which allows to fix the parameter (not to ...
Double_t Derby(Int_t numOfFittedPeaks, Double_t x, Double_t y, const Double_t *parameter, Double_t sigmax, Double_t sigmay, Double_t txy, Double_t ty, Double_t bx, Double_t by)
This function calculates derivative of peaks shape function (see manual) according to slope by.
Double_t fTxyInit
initial value of t parameter for 2D peaks (relative amplitude of tail), for details see html manual a...
Bool_t fFixSigmaY
logical value of sigma y parameter, which allows to fix the parameter (not to fit).
Double_t Deramp2(Double_t x, Double_t y, Double_t x0, Double_t y0, Double_t sigmax, Double_t sigmay, Double_t ro, Double_t txy, Double_t sxy, Double_t bx, Double_t by)
This function calculates derivative of 2D peaks shape function (see manual) according to amplitude of...
void GetRo(Double_t &ro, Double_t &roErr)
This function gets the ro parameter and its error.
Int_t fXmax
last fitted channel in x direction
Double_t fChi
here the fitting functions return resulting chi square
void GetVolumes(Double_t *volumes)
This function gets the volumes of fitted 2D peaks.
Double_t fByErr
error value of b parameter for 1D ridges in y direction
Double_t fA0Init
initial value of background a0 parameter(backgroud is estimated as a0+ax*x+ay*y)
Double_t fRoInit
initial value of correlation coefficient
Double_t Deri02(Double_t x, Double_t y, Double_t a, Double_t x0, Double_t y0, Double_t sigmax, Double_t sigmay, Double_t ro, Double_t txy, Double_t sxy, Double_t bx, Double_t by)
This function calculates derivative of 2D peaks shape function (see manual) according to x position o...
Double_t fSigmaInitX
initial value of sigma x parameter
void GetBackgroundParameters(Double_t &a0, Double_t &a0Err, Double_t &ax, Double_t &axErr, Double_t &ay, Double_t &ayErr)
This function gets the background parameters and their errors.
Double_t Derampx(Double_t x, Double_t x0, Double_t sigmax, Double_t tx, Double_t sx, Double_t bx)
This function calculates derivative of 2D peaks shape function (see manual) according to amplitude of...
Double_t * fPositionInitY1
[fNPeaks] array of initial y positions of 1D ridges, input parameters
Bool_t fFixSy
logical value of s parameter for 1D ridges in y direction, which allows to fix the parameter (not to ...
Bool_t * fFixAmp
[fNPeaks] array of logical values which allow to fix appropriate amplitudes of 2D peaks (not fit)....
Double_t Derderj02(Double_t x, Double_t y, Double_t a, Double_t x0, Double_t y0, Double_t sigmax, Double_t sigmay, Double_t ro)
This function calculates second derivative of 2D peaks shape function (see manual) according to y pos...
Double_t fAyErr
error value of background ay parameter
Double_t fSxyCalc
calculated value of s parameter for 2D peaks
Double_t * fPositionInitX1
[fNPeaks] array of initial x positions of 1D ridges, input parameters
Double_t fTyCalc
calculated value of t parameter for 1D ridges in y direction
Double_t * fPositionErrX1
[fNPeaks] array of x positions errors of 1D ridges, output parameters
Double_t Dertxy(Int_t numOfFittedPeaks, Double_t x, Double_t y, const Double_t *parameter, Double_t sigmax, Double_t sigmay, Double_t bx, Double_t by)
This function calculates derivative of peaks shape function (see manual) according to relative amplit...
Double_t Derj02(Double_t x, Double_t y, Double_t a, Double_t x0, Double_t y0, Double_t sigmax, Double_t sigmay, Double_t ro, Double_t txy, Double_t sxy, Double_t bx, Double_t by)
This function calculates derivative of 2D peaks shape function (see manual) according to y position o...
Double_t * fPositionInitY
[fNPeaks] array of initial values of y positions of 2D peaks, input parameters
Double_t fSxCalc
calculated value of s parameter for 1D ridges in x direction
void GetAmplitudeErrors(Double_t *amplitudeErrors, Double_t *amplitudeErrorsX1, Double_t *amplitudeErrorsY1)
This function gets the amplitudes of fitted 2D peaks and 1D ridges.
Double_t fAxCalc
calculated value of background ax parameter
Double_t fSyErr
error value of s parameter for 1D ridges in y direction
Int_t fStatisticType
type of statistics, possible values kFitOptimChiCounts (chi square statistics with counts as weightin...
Int_t fYmax
last fitted channel in y direction
Double_t * fPositionErrY
[fNPeaks] array of error values of y positions of 2D peaks, output parameters
Double_t Dersx(Int_t numOfFittedPeaks, Double_t x, const Double_t *parameter, Double_t sigmax)
This function calculates derivative of peaks shape function (see manual) according to relative amplit...
Double_t fByCalc
calculated value of b parameter for 1D ridges in y direction
Double_t fTyInit
initial value of t parameter for 1D ridges in y direction (relative amplitude of tail),...
Double_t fTxyCalc
calculated value of t parameter for 2D peaks
Double_t Derpa2(Double_t sx, Double_t sy, Double_t ro)
This function calculates derivative of the volume of a peak according to amplitude.
Double_t fByInit
initial value of b parameter for 1D ridges in y direction (slope), for details see html manual and re...
Double_t fSigmaCalcY
calculated value of sigma y parameter
Double_t * fAmpErrX1
[fNPeaks] array of amplitudes errors of 1D ridges in x direction, output parameters
Double_t Dertx(Int_t numOfFittedPeaks, Double_t x, const Double_t *parameter, Double_t sigmax, Double_t bx)
This function calculates derivative of peaks shape function (see manual) according to relative amplit...
Int_t fNumberIterations
number of iterations in fitting procedure, input parameter, it should be > 0
Bool_t fFixTx
logical value of t parameter for 1D ridges in x direction, which allows to fix the parameter (not to ...
Double_t * fPositionCalcY1
[fNPeaks] array of calculated y positions of 1D ridges, output parameters
Bool_t fFixTxy
logical value of t parameter for 2D peaks, which allows to fix the parameter (not to fit).
void FitAwmi(Double_t **source)
This function fits the source spectrum.
Bool_t fFixTy
logical value of t parameter for 1D ridges in y direction, which allows to fix the parameter (not to ...
Double_t fAlpha
convergence coefficient, input parameter, it should be positive number and <=1, for details see refer...
void GetPositions(Double_t *positionsX, Double_t *positionsY, Double_t *positionsX1, Double_t *positionsY1)
This function gets the positions of fitted 2D peaks and 1D ridges.
Bool_t fFixAx
logical value of ax parameter, which allows to fix the parameter (not to fit).
Double_t fTyErr
error value of t parameter for 1D ridges in y direction
Double_t fSyCalc
calculated value of s parameter for 1D ridges in y direction
Int_t fXmin
first fitted channel in x direction
Double_t fTxyErr
error value of t parameter for 2D peaks
void SetFitParameters(Int_t xmin, Int_t xmax, Int_t ymin, Int_t ymax, 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:
void GetSigmaX(Double_t &sigmaX, Double_t &sigmaErrX)
This function gets the sigma x parameter and its error.
Double_t Derderi01(Double_t x, Double_t ax, Double_t x0, Double_t sigmax)
This function calculates second derivative of 2D peaks shape function (see manual) according to x pos...
Double_t Erfc(Double_t x)
This function calculates error function of x.
Double_t Derro(Int_t numOfFittedPeaks, Double_t x, Double_t y, const Double_t *parameter, Double_t sx, Double_t sy, Double_t r)
This function calculates derivative of peaks shape function (see manual) according to correlation coe...
Bool_t * fFixPositionX1
[fNPeaks] array of logical values which allow to fix appropriate x positions of 1D ridges (not fit)....
Double_t Dersigmay(Int_t numOfFittedPeaks, Double_t x, Double_t y, const Double_t *parameter, Double_t sigmax, Double_t sigmay, Double_t ro, Double_t txy, Double_t sxy, Double_t ty, Double_t sy, Double_t bx, Double_t by)
This function calculates derivative of peaks shape function (see manual) according to sigmax of peaks...
Double_t Ourpowl(Double_t a, Int_t pw)
power function
Double_t * fAmpInitY1
[fNPeaks] array of initial values of amplitudes of 1D ridges in y direction, input parameters
Double_t fAxErr
error value of background ax parameter
Double_t Derdersigmay(Int_t numOfFittedPeaks, Double_t x, Double_t y, const Double_t *parameter, Double_t sigmax, Double_t sigmay, Double_t ro)
This function calculates second derivative of peaks shape function (see manual) according to sigmay o...
Double_t fAyCalc
calculated value of background ay parameter
Double_t fSigmaInitY
initial value of sigma y parameter
void SetTailParameters(Double_t tInitXY, Bool_t fixTxy, Double_t tInitX, Bool_t fixTx, Double_t tInitY, Bool_t fixTy, Double_t bInitX, Bool_t fixBx, Double_t bInitY, Bool_t fixBy, Double_t sInitXY, Bool_t fixSxy, Double_t sInitX, Bool_t fixSx, Double_t sInitY, Bool_t fixSy)
This function sets the following fitting parameters of tails of peaks.
Double_t * fPositionErrX
[fNPeaks] array of error values of x positions of 2D peaks, output parameters
Double_t * fAmpInitX1
[fNPeaks] array of initial values of amplitudes of 1D ridges in x direction, input parameters
Bool_t * fFixPositionY1
[fNPeaks] array of logical values which allow to fix appropriate y positions of 1D ridges (not fit)....
void GetVolumeErrors(Double_t *volumeErrors)
This function gets errors of the volumes of fitted 2D peaks.
Double_t Dersy(Int_t numOfFittedPeaks, Double_t x, const Double_t *parameter, Double_t sigmax)
This function calculates derivative of peaks shape function (see manual) according to relative amplit...
~TSpectrum2Fit() override
Destructor.
Double_t Derpsigmay(Double_t a, Double_t sx, Double_t ro)
This function calculates derivative of the volume of a peak according to sigmay.
Double_t Derty(Int_t numOfFittedPeaks, Double_t x, const Double_t *parameter, Double_t sigmax, Double_t bx)
This function calculates derivative of peaks shape function (see manual) according to relative amplit...
Double_t fTxInit
initial value of t parameter for 1D ridges in x direction (relative amplitude of tail),...
Double_t fSyInit
initial value of s parameter for 1D ridges in y direction (relative amplitude of step),...
Double_t Dersigmax(Int_t numOfFittedPeaks, Double_t x, Double_t y, const Double_t *parameter, Double_t sigmax, Double_t sigmay, Double_t ro, Double_t txy, Double_t sxy, Double_t tx, Double_t sx, Double_t bx, Double_t by)
This function calculates derivative of peaks shape function (see manual) according to sigmax of peaks...
Double_t * fAmpCalc
[fNPeaks] array of calculated values of amplitudes of 2D peaks, output parameters
void GetSigmaY(Double_t &sigmaY, Double_t &sigmaErrY)
This function gets the sigma y parameter and its error.
TSpectrum2Fit(void)
Default constructor.
Double_t fA0Calc
calculated value of background a0 parameter
Double_t * fPositionInitX
[fNPeaks] array of initial values of x positions of 2D peaks, input parameters
Int_t fAlphaOptim
optimization of convergence algorithm, possible values kFitAlphaHalving, kFitAlphaOptimal
Bool_t fFixSigmaX
logical value of sigma x parameter, which allows to fix the parameter (not to fit).
Double_t fSxErr
error value of s parameter for 1D ridges in x direction
Int_t fNPeaks
number of peaks present in fit, input parameter, it should be > 0
void SetBackgroundParameters(Double_t a0Init, Bool_t fixA0, Double_t axInit, Bool_t fixAx, Double_t ayInit, Bool_t fixAy)
This function sets the following fitting parameters of background:
Double_t Derpro(Double_t a, Double_t sx, Double_t sy, Double_t ro)
This function calculates derivative of the volume of a peak according to ro.
Double_t * fAmpCalcY1
[fNPeaks] array of calculated values of amplitudes of 1D ridges in y direction, output parameters
Bool_t fFixAy
logical value of ay parameter, which allows to fix the parameter (not to fit).
void FitStiefel(Double_t **source)
This function fits the source spectrum.
Double_t Derfc(Double_t x)
This function calculates derivative of error function of x.
Double_t fSxInit
initial value of s parameter for 1D ridges in x direction (relative amplitude of step),...
Double_t fSigmaErrX
error value of sigma x parameter
Int_t fYmin
first fitted channel in y direction
Double_t * fPositionCalcX1
[fNPeaks] array of calculated x positions of 1D ridges, output parameters
void GetAmplitudes(Double_t *amplitudes, Double_t *amplitudesX1, Double_t *amplitudesY1)
This function gets the amplitudes of fitted 2D peaks and 1D ridges.
Double_t fBxErr
error value of b parameter for 1D ridges in x direction
Double_t Derpsigmax(Double_t a, Double_t sy, Double_t ro)
This function calculates derivative of the volume of a peak according to sigmax.
Double_t * fPositionCalcX
[fNPeaks] array of calculated values of x positions of 2D peaks, output parameters
Double_t fAxInit
initial value of background ax parameter(backgroud is estimated as a0+ax*x+ay*y)
Bool_t * fFixPositionX
[fNPeaks] array of logical values which allow to fix appropriate x positions of 2D peaks (not fit)....
void GetPositionErrors(Double_t *positionErrorsX, Double_t *positionErrorsY, Double_t *positionErrorsX1, Double_t *positionErrorsY1)
This function gets the errors of positions of fitted 2D peaks and 1D ridges.
void StiefelInversion(Double_t **a, Int_t size)
This function calculates solution of the system of linear equations.
Double_t fSxyInit
initial value of s parameter for 2D peaks (relative amplitude of step), for details see html manual a...
Int_t fPower
possible values kFitPower2,4,6,8,10,12, for details see references. It applies only for Awmi fitting ...
Bool_t fFixA0
logical value of a0 parameter, which allows to fix the parameter (not to fit).
Double_t fAyInit
initial value of background ay parameter(backgroud is estimated as a0+ax*x+ay*y)
Double_t * fVolumeErr
[fNPeaks] array of volumes errors of 2D peaks, output parameters
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17