ROOT logo
ROOT » HIST » SPECTRUM » TSpectrum2Fit

class TSpectrum2Fit: public TNamed

THIS CLASS CONTAINS ADVANCED TWO-DIMENSIONAL SPECTRA
FITTING FUNCTIONS

These functions were written by:
Miroslav Morhac
Institute of Physics
Slovak Academy of Sciences
Dubravska cesta 9, 842 28 BRATISLAVA
SLOVAKIA

email:fyzimiro@savba.sk,    fax:+421 7 54772479

The original code in C has been repackaged as a C++ class by R.Brun

The algorithms in this class have been published in the following
references:
[1] M. Morhac et al.: Efficient fitting algorithms applied to
analysis of coincidence gamma-ray spectra. Computer Physics
Communications, Vol 172/1 (2005) pp. 19-41.

[2]  M. Morhac et al.: Study of fitting algorithms applied to
simultaneous analysis of large number of peaks in gamma-ray spectra.
Applied Spectroscopy, Vol. 57, No. 7, pp. 753-760, 2003.



Function Members (Methods)

public:
TSpectrum2Fit()
TSpectrum2Fit(Int_t numberPeaks)
TSpectrum2Fit(const TSpectrum2Fit&)
virtual~TSpectrum2Fit()
voidTObject::AbstractMethod(const char* method) const
virtual voidTObject::AppendPad(Option_t* option = "")
virtual voidTObject::Browse(TBrowser* b)
static TClass*Class()
virtual const char*TObject::ClassName() const
virtual voidTNamed::Clear(Option_t* option = "")
virtual TObject*TNamed::Clone(const char* newname = "") const
virtual Int_tTNamed::Compare(const TObject* obj) const
virtual voidTNamed::Copy(TObject& named) const
virtual voidTObject::Delete(Option_t* option = "")MENU
virtual Int_tTObject::DistancetoPrimitive(Int_t px, Int_t py)
virtual voidTObject::Draw(Option_t* option = "")
virtual voidTObject::DrawClass() constMENU
virtual TObject*TObject::DrawClone(Option_t* option = "") constMENU
virtual voidTObject::Dump() constMENU
virtual voidTObject::Error(const char* method, const char* msgfmt) const
virtual voidTObject::Execute(const char* method, const char* params, Int_t* error = 0)
virtual voidTObject::Execute(TMethod* method, TObjArray* params, Int_t* error = 0)
virtual voidTObject::ExecuteEvent(Int_t event, Int_t px, Int_t py)
virtual voidTObject::Fatal(const char* method, const char* msgfmt) const
virtual voidTNamed::FillBuffer(char*& buffer)
virtual TObject*TObject::FindObject(const char* name) const
virtual TObject*TObject::FindObject(const TObject* obj) const
voidFitAwmi(Float_t** source)
voidFitStiefel(Float_t** source)
voidGetAmplitudeErrors(Float_t* amplitudeErrors, Float_t* amplitudeErrorsX1, Float_t* amplitudeErrorsY1)
voidGetAmplitudes(Float_t* amplitudes, Float_t* amplitudesX1, Float_t* amplitudesY1)
voidGetBackgroundParameters(Double_t& a0, Double_t& a0Err, Double_t& ax, Double_t& axErr, Double_t& ay, Double_t& ayErr)
Double_tGetChi() const
virtual Option_t*TObject::GetDrawOption() const
static Long_tTObject::GetDtorOnly()
virtual const char*TObject::GetIconName() const
virtual const char*TNamed::GetName() const
virtual char*TObject::GetObjectInfo(Int_t px, Int_t py) const
static Bool_tTObject::GetObjectStat()
virtual Option_t*TObject::GetOption() const
voidGetPositionErrors(Float_t* positionErrorsX, Float_t* positionErrorsY, Float_t* positionErrorsX1, Float_t* positionErrorsY1)
voidGetPositions(Float_t* positionsX, Float_t* positionsY, Float_t* positionsX1, Float_t* positionsY1)
voidGetRo(Double_t& ro, Double_t& roErr)
voidGetSigmaX(Double_t& sigmaX, Double_t& sigmaErrX)
voidGetSigmaY(Double_t& sigmaY, Double_t& sigmaErrY)
voidGetTailParameters(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)
virtual const char*TNamed::GetTitle() const
virtual UInt_tTObject::GetUniqueID() const
voidGetVolumeErrors(Float_t* volumeErrors)
voidGetVolumes(Float_t* volumes)
virtual Bool_tTObject::HandleTimer(TTimer* timer)
virtual ULong_tTNamed::Hash() const
virtual voidTObject::Info(const char* method, const char* msgfmt) const
virtual Bool_tTObject::InheritsFrom(const char* classname) const
virtual Bool_tTObject::InheritsFrom(const TClass* cl) const
virtual voidTObject::Inspect() constMENU
voidTObject::InvertBit(UInt_t f)
virtual TClass*IsA() const
virtual Bool_tTObject::IsEqual(const TObject* obj) const
virtual Bool_tTObject::IsFolder() const
Bool_tTObject::IsOnHeap() const
virtual Bool_tTNamed::IsSortable() const
Bool_tTObject::IsZombie() const
virtual voidTNamed::ls(Option_t* option = "") const
voidTObject::MayNotUse(const char* method) const
virtual Bool_tTObject::Notify()
static voidTObject::operator delete(void* ptr)
static voidTObject::operator delete(void* ptr, void* vp)
static voidTObject::operator delete[](void* ptr)
static voidTObject::operator delete[](void* ptr, void* vp)
void*TObject::operator new(size_t sz)
void*TObject::operator new(size_t sz, void* vp)
void*TObject::operator new[](size_t sz)
void*TObject::operator new[](size_t sz, void* vp)
TSpectrum2Fit&operator=(const TSpectrum2Fit&)
virtual voidTObject::Paint(Option_t* option = "")
virtual voidTObject::Pop()
virtual voidTNamed::Print(Option_t* option = "") const
virtual Int_tTObject::Read(const char* name)
virtual voidTObject::RecursiveRemove(TObject* obj)
voidTObject::ResetBit(UInt_t f)
virtual voidTObject::SaveAs(const char* filename = "", Option_t* option = "") constMENU
virtual voidTObject::SavePrimitive(basic_ostream<char,char_traits<char> >& out, Option_t* option = "")
voidSetBackgroundParameters(Double_t a0Init, Bool_t fixA0, Double_t axInit, Bool_t fixAx, Double_t ayInit, Bool_t fixAy)
voidTObject::SetBit(UInt_t f)
voidTObject::SetBit(UInt_t f, Bool_t set)
virtual voidTObject::SetDrawOption(Option_t* option = "")MENU
static voidTObject::SetDtorOnly(void* obj)
voidSetFitParameters(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)
virtual voidTNamed::SetName(const char* name)MENU
virtual voidTNamed::SetNameTitle(const char* name, const char* title)
static voidTObject::SetObjectStat(Bool_t stat)
voidSetPeakParameters(Double_t sigmaX, Bool_t fixSigmaX, Double_t sigmaY, Bool_t fixSigmaY, Double_t ro, Bool_t fixRo, const Float_t* positionInitX, const Bool_t* fixPositionX, const Float_t* positionInitY, const Bool_t* fixPositionY, const Float_t* positionInitX1, const Bool_t* fixPositionX1, const Float_t* positionInitY1, const Bool_t* fixPositionY1, const Float_t* ampInit, const Bool_t* fixAmp, const Float_t* ampInitX1, const Bool_t* fixAmpX1, const Float_t* ampInitY1, const Bool_t* fixAmpY1)
voidSetTailParameters(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)
virtual voidTNamed::SetTitle(const char* title = "")MENU
virtual voidTObject::SetUniqueID(UInt_t uid)
virtual voidShowMembers(TMemberInspector& insp, char* parent)
virtual Int_tTNamed::Sizeof() const
virtual voidStreamer(TBuffer& b)
voidStreamerNVirtual(TBuffer& b)
virtual voidTObject::SysError(const char* method, const char* msgfmt) const
Bool_tTObject::TestBit(UInt_t f) const
Int_tTObject::TestBits(UInt_t f) const
virtual voidTObject::UseCurrentStyle()
virtual voidTObject::Warning(const char* method, const char* msgfmt) const
virtual Int_tTObject::Write(const char* name = 0, Int_t option = 0, Int_t bufsize = 0)
virtual Int_tTObject::Write(const char* name = 0, Int_t option = 0, Int_t bufsize = 0) const
protected:
Double_tDeramp2(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)
Double_tDerampx(Double_t x, Double_t x0, Double_t sigmax, Double_t tx, Double_t sx, Double_t bx)
Double_tDerbx(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)
Double_tDerby(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)
Double_tDerderi01(Double_t x, Double_t ax, Double_t x0, Double_t sigmax)
Double_tDerderi02(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_tDerderj02(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_tDerdersigmax(Int_t numOfFittedPeaks, Double_t x, Double_t y, const Double_t* parameter, Double_t sigmax, Double_t sigmay, Double_t ro)
Double_tDerdersigmay(Int_t numOfFittedPeaks, Double_t x, Double_t y, const Double_t* parameter, Double_t sigmax, Double_t sigmay, Double_t ro)
Double_tDerfc(Double_t x)
Double_tDeri01(Double_t x, Double_t ax, Double_t x0, Double_t sigmax, Double_t tx, Double_t sx, Double_t bx)
Double_tDeri02(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)
Double_tDerj02(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)
Double_tDerpa2(Double_t sx, Double_t sy, Double_t ro)
Double_tDerpro(Double_t a, Double_t sx, Double_t sy, Double_t ro)
Double_tDerpsigmax(Double_t a, Double_t sy, Double_t ro)
Double_tDerpsigmay(Double_t a, Double_t sx, Double_t ro)
Double_tDerro(Int_t numOfFittedPeaks, Double_t x, Double_t y, const Double_t* parameter, Double_t sx, Double_t sy, Double_t r)
Double_tDersigmax(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)
Double_tDersigmay(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)
Double_tDersx(Int_t numOfFittedPeaks, Double_t x, const Double_t* parameter, Double_t sigmax)
Double_tDersxy(Int_t numOfFittedPeaks, Double_t x, Double_t y, const Double_t* parameter, Double_t sigmax, Double_t sigmay)
Double_tDersy(Int_t numOfFittedPeaks, Double_t x, const Double_t* parameter, Double_t sigmax)
Double_tDertx(Int_t numOfFittedPeaks, Double_t x, const Double_t* parameter, Double_t sigmax, Double_t bx)
Double_tDertxy(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)
Double_tDerty(Int_t numOfFittedPeaks, Double_t x, const Double_t* parameter, Double_t sigmax, Double_t bx)
virtual voidTObject::DoError(int level, const char* location, const char* fmt, va_list va) const
Double_tErfc(Double_t x)
voidTObject::MakeZombie()
Double_tOurpowl(Double_t a, Int_t pw)
Double_tShape2(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)
voidStiefelInversion(Double_t** a, Int_t size)
Double_tVolume(Double_t a, Double_t sx, Double_t sy, Double_t ro)

Data Members

public:
enum { kFitOptimChiCounts
kFitOptimChiFuncValues
kFitOptimMaxLikelihood
kFitAlphaHalving
kFitAlphaOptimal
kFitPower2
kFitPower4
kFitPower6
kFitPower8
kFitPower10
kFitPower12
kFitTaylorOrderFirst
kFitTaylorOrderSecond
kFitNumRegulCycles
};
enum TObject::EStatusBits { kCanDelete
kMustCleanup
kObjInCanvas
kIsReferenced
kHasUUID
kCannotPick
kNoContextMenu
kInvalidObject
};
enum TObject::[unnamed] { kIsOnHeap
kNotDeleted
kZombie
kBitMask
kSingleKey
kOverwrite
kWriteDelete
};
protected:
Double_tfA0Calccalculated value of background a0 parameter
Double_tfA0Errerror value of background a0 parameter
Double_tfA0Initinitial value of background a0 parameter(backgroud is estimated as a0+ax*x+ay*y)
Double_tfAlphaconvergence coefficient, input parameter, it should be positive number and <=1, for details see references
Int_tfAlphaOptimoptimization of convergence algorithm, possible values kFitAlphaHalving, kFitAlphaOptimal
Double_t*fAmpCalc[fNPeaks] array of calculated values of amplitudes of 2D peaks, output parameters
Double_t*fAmpCalcX1[fNPeaks] array of calculated values of amplitudes of 1D ridges in x direction, output parameters
Double_t*fAmpCalcY1[fNPeaks] array of calculated values of amplitudes of 1D ridges in y direction, output parameters
Double_t*fAmpErr[fNPeaks] array of amplitudes errors of 2D peaks, output parameters
Double_t*fAmpErrX1[fNPeaks] array of amplitudes errors of 1D ridges in x direction, output parameters
Double_t*fAmpErrY1[fNPeaks] array of amplitudes errors of 1D ridges in y direction, output parameters
Double_t*fAmpInit[fNPeaks] array of initial values of amplitudes of 2D peaks, input parameters
Double_t*fAmpInitX1[fNPeaks] array of initial values of amplitudes of 1D ridges in x direction, input parameters
Double_t*fAmpInitY1[fNPeaks] array of initial values of amplitudes of 1D ridges in y direction, input parameters
Double_tfAxCalccalculated value of background ax parameter
Double_tfAxErrerror value of background ax parameter
Double_tfAxInitinitial value of background ax parameter(backgroud is estimated as a0+ax*x+ay*y)
Double_tfAyCalccalculated value of background ay parameter
Double_tfAyErrerror value of background ay parameter
Double_tfAyInitinitial value of background ay parameter(backgroud is estimated as a0+ax*x+ay*y)
Double_tfBxCalccalculated value of b parameter for 1D ridges in x direction
Double_tfBxErrerror value of b parameter for 1D ridges in x direction
Double_tfBxInitinitial value of b parameter for 1D ridges in x direction (slope), for details see html manual and references
Double_tfByCalccalculated value of b parameter for 1D ridges in y direction
Double_tfByErrerror value of b parameter for 1D ridges in y direction
Double_tfByInitinitial value of b parameter for 1D ridges in y direction (slope), for details see html manual and references
Double_tfChihere the fitting functions return resulting chi square
Int_tfFitTaylororder of Taylor expansion, possible values kFitTaylorOrderFirst, kFitTaylorOrderSecond. It applies only for Awmi fitting function.
Bool_tfFixA0logical value of a0 parameter, which allows to fix the parameter (not to fit).
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
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
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
Bool_tfFixAxlogical value of ax parameter, which allows to fix the parameter (not to fit).
Bool_tfFixAylogical value of ay parameter, which allows to fix the parameter (not to fit).
Bool_tfFixBxlogical value of b parameter for 1D ridges in x direction, which allows to fix the parameter (not to fit).
Bool_tfFixBylogical value of b parameter for 1D ridges in y direction, which allows to fix the parameter (not to fit).
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
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
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
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
Bool_tfFixRological value of correlation coefficient, which allows to fix the parameter (not to fit).
Bool_tfFixSigmaXlogical value of sigma x parameter, which allows to fix the parameter (not to fit).
Bool_tfFixSigmaYlogical value of sigma y parameter, which allows to fix the parameter (not to fit).
Bool_tfFixSxlogical value of s parameter for 1D ridges in x direction, which allows to fix the parameter (not to fit).
Bool_tfFixSxylogical value of s parameter for 2D peaks, which allows to fix the parameter (not to fit).
Bool_tfFixSylogical value of s parameter for 1D ridges in y direction, which allows to fix the parameter (not to fit).
Bool_tfFixTxlogical value of t parameter for 1D ridges in x direction, which allows to fix the parameter (not to fit).
Bool_tfFixTxylogical value of t parameter for 2D peaks, which allows to fix the parameter (not to fit).
Bool_tfFixTylogical value of t parameter for 1D ridges in y direction, which allows to fix the parameter (not to fit).
Int_tfNPeaksnumber of peaks present in fit, input parameter, it should be > 0
TStringTNamed::fNameobject identifier
Int_tfNumberIterationsnumber of iterations in fitting procedure, input parameter, it should be > 0
Double_t*fPositionCalcX[fNPeaks] array of calculated values of x positions of 2D peaks, output parameters
Double_t*fPositionCalcX1[fNPeaks] array of calculated x positions of 1D ridges, output parameters
Double_t*fPositionCalcY[fNPeaks] array of calculated values of y positions of 2D peaks, output parameters
Double_t*fPositionCalcY1[fNPeaks] array of calculated y positions of 1D ridges, output parameters
Double_t*fPositionErrX[fNPeaks] array of error values of x positions of 2D peaks, output parameters
Double_t*fPositionErrX1[fNPeaks] array of x positions errors of 1D ridges, output parameters
Double_t*fPositionErrY[fNPeaks] array of error values of y positions of 2D peaks, output parameters
Double_t*fPositionErrY1[fNPeaks] array of y positions errors of 1D ridges, output parameters
Double_t*fPositionInitX[fNPeaks] array of initial values of x positions of 2D peaks, input parameters
Double_t*fPositionInitX1[fNPeaks] array of initial x positions of 1D ridges, input parameters
Double_t*fPositionInitY[fNPeaks] array of initial values of y positions of 2D peaks, input parameters
Double_t*fPositionInitY1[fNPeaks] array of initial y positions of 1D ridges, input parameters
Int_tfPowerpossible values kFitPower2,4,6,8,10,12, for details see references. It applies only for Awmi fitting function.
Double_tfRoCalccalculated value of correlation coefficient
Double_tfRoErrerror value of correlation coefficient
Double_tfRoInitinitial value of correlation coefficient
Double_tfSigmaCalcXcalculated value of sigma x parameter
Double_tfSigmaCalcYcalculated value of sigma y parameter
Double_tfSigmaErrXerror value of sigma x parameter
Double_tfSigmaErrYerror value of sigma y parameter
Double_tfSigmaInitXinitial value of sigma x parameter
Double_tfSigmaInitYinitial value of sigma y parameter
Int_tfStatisticTypetype of statistics, possible values kFitOptimChiCounts (chi square statistics with counts as weighting coefficients), kFitOptimChiFuncValues (chi square statistics with function values as weighting coefficients),kFitOptimMaxLikelihood
Double_tfSxCalccalculated value of s parameter for 1D ridges in x direction
Double_tfSxErrerror value of s parameter for 1D ridges in x direction
Double_tfSxInitinitial value of s parameter for 1D ridges in x direction (relative amplitude of step), for details see html manual and references
Double_tfSxyCalccalculated value of s parameter for 2D peaks
Double_tfSxyErrerror value of s parameter for 2D peaks
Double_tfSxyInitinitial value of s parameter for 2D peaks (relative amplitude of step), for details see html manual and references
Double_tfSyCalccalculated value of s parameter for 1D ridges in y direction
Double_tfSyErrerror value of s parameter for 1D ridges in y direction
Double_tfSyInitinitial value of s parameter for 1D ridges in y direction (relative amplitude of step), for details see html manual and references
TStringTNamed::fTitleobject title
Double_tfTxCalccalculated value of t parameter for 1D ridges in x direction
Double_tfTxErrerror value of t parameter for 1D ridges in x direction
Double_tfTxInitinitial value of t parameter for 1D ridges in x direction (relative amplitude of tail), for details see html manual and references
Double_tfTxyCalccalculated value of t parameter for 2D peaks
Double_tfTxyErrerror value of t parameter for 2D peaks
Double_tfTxyInitinitial value of t parameter for 2D peaks (relative amplitude of tail), for details see html manual and references
Double_tfTyCalccalculated value of t parameter for 1D ridges in y direction
Double_tfTyErrerror value of t parameter for 1D ridges in y direction
Double_tfTyInitinitial value of t parameter for 1D ridges in y direction (relative amplitude of tail), for details see html manual and references
Double_t*fVolume[fNPeaks] array of calculated volumes of 2D peaks, output parameters
Double_t*fVolumeErr[fNPeaks] array of volumes errors of 2D peaks, output parameters
Int_tfXmaxlast fitted channel in x direction
Int_tfXminfirst fitted channel in x direction
Int_tfYmaxlast fitted channel in y direction
Int_tfYminfirst fitted channel in y direction

Class Charts

Inheritance Inherited Members Includes Libraries
Class Charts

Function documentation

TSpectrum2Fit()
default constructor
TSpectrum2Fit(Int_t numberPeaks)
numberPeaks: number of fitted peaks (must be greater than zero)
the constructor allocates arrays for all fitted parameters (peak positions, amplitudes etc) and sets the member
variables to their default values. One can change these variables by member functions (setters) of TSpectrumFit class.

Shape function of the fitted peaks contains the two-dimensional symmetrical Gaussian two one-dimensional symmetrical Gaussian ridges as well as nonsymmetrical terms and background.

~TSpectrum2Fit()
 destructor
Double_t Erfc(Double_t x)
AUXILIARY FUNCTION

This function calculates error function of x.


Double_t Derfc(Double_t x)
AUXILIARY FUNCTION

This function calculates derivative of error function of x.


Double_t Ourpowl(Double_t a, Int_t pw)
power function
void StiefelInversion(Double_t** a, Int_t size)
AUXILIARY FUNCTION

This function calculates solution of the system of linear equations.
The matrix a should have a dimension size*(size+4)
The calling function should fill in the matrix, the column size should
contain vector y (right side of the system of equations). The result is
placed into size+1 column of the matrix.
according to sigma of peaks.
Function parameters:
-a-matrix with dimension size*(size+4)                          //
-size-number of rows of the matrix


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)
AUXILIARY FUNCTION

This function calculates 2D peaks shape function (see manual)
Function parameters:
-numOfFittedPeaks-number of fitted peaks
-x-channel in x-dimension
-y-channel in y-dimension
-parameter-array of peaks parameters (amplitudes and positions)
-sigmax-sigmax of peaks
-sigmay-sigmay of peaks
-ro-correlation coefficient
-a0,ax,ay-bac kground coefficients
-txy,tx,ty, sxy,sy,sx-relative amplitudes
-bx, by-slopes


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)
AUXILIARY FUNCTION

This function calculates derivative of 2D peaks shape function (see manual)
according to amplitude of 2D peak
Function parameters:
-x-channel in x-dimension
-y-channel in y-dimension
-x0-position of peak in x-dimension
-y0-position of peak in y-dimension
-sigmax-sigmax of peaks
-sigmay-sigmay of peaks
-ro-correlation coefficient
-txy, sxy-relative amplitudes
-bx, by-slopes


Double_t Derampx(Double_t x, Double_t x0, Double_t sigmax, Double_t tx, Double_t sx, Double_t bx)
AUXILIARY FUNCTION

This function calculates derivative of 2D peaks shape function (see manual)
according to amplitude of the ridge
Function parameters:
-x-channel in x-dimension
-x0-position of peak in x-dimension
-y0-position of peak in y-dimension
-sigmax-sigmax of peaks
-ro-correlation coefficient
-tx, sx-relative amplitudes
-bx-slope


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)
AUXILIARY FUNCTION

This function calculates derivative of 2D peaks shape function (see manual)
according to x position of 2D peak
Function parameters:
-x-channel in x-dimension
-y-channel in y-dimension
-a-amplitude
-x0-position of peak in x-dimension
-y0-position of peak in y-dimension
-sigmax-sigmax of peaks
-sigmay-sigmay of peaks
-ro-correlation coefficient
-txy, sxy-relative amplitudes
-bx, by-slopes


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)
AUXILIARY FUNCTION

This function calculates second derivative of 2D peaks shape function
(see manual) according to x position of 2D peak
Function parameters:
-x-channel in x-dimension
-y-channel in y-dimension
-a-amplitude
-x0-position of peak in x-dimension
-y0-position of peak in y-dimension
-sigmax-sigmax of peaks
-sigmay-sigmay of peaks
-ro-correlation coefficient


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)
AUXILIARY FUNCTION

This function calculates derivative of 2D peaks shape function (see manual)
according to y position of 2D peak
Function parameters:
-x-channel in x-dimension
-y-channel in y-dimension
-a-amplitude
-x0-position of peak in x-dimension
-y0-position of peak in y-dimension
-sigmax-sigmax of peaks
-sigmay-sigmay of peaks
-ro-correlation coefficient
-txy, sxy-relative amplitudes
-bx, by-slopes


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)
AUXILIARY FUNCTION

This function calculates second derivative of 2D peaks shape function
(see manual) according to y position of 2D peak
Function parameters:
-x-channel in x-dimension
-y-channel in y-dimension
-a-amplitude
-x0-position of peak in x-dimension
-y0-position of peak in y-dimension
-sigmax-sigmax of peaks
-sigmay-sigmay of peaks
-ro-correlation coefficient


Double_t Deri01(Double_t x, Double_t ax, Double_t x0, Double_t sigmax, Double_t tx, Double_t sx, Double_t bx)
AUXILIARY FUNCTION

This function calculates derivative of 2D peaks shape function (see manual)
according to x position of 1D ridge
Function parameters:
-x-channel in x-dimension
-ax-amplitude of ridge
-x0-position of peak in x-dimension
-sigmax-sigmax of peaks
-ro-correlation coefficient
-tx, sx-relative amplitudes
-bx-slope


Double_t Derderi01(Double_t x, Double_t ax, Double_t x0, Double_t sigmax)
AUXILIARY FUNCTION

This function calculates second derivative of 2D peaks shape function
(see manual) according to x position of 1D ridge
Function parameters:
-x-channel in x-dimension
-ax-amplitude of ridge
-x0-position of peak in x-dimension
-sigmax-sigmax of peaks


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)
AUXILIARY FUNCTION

This function calculates derivative of peaks shape function (see manual)
according to sigmax of peaks.
Function parameters:
-numOfFittedPeaks-number of fitted peaks
-x,y-position of channel
-parameter-array of peaks parameters (amplitudes and positions)
-sigmax-sigmax of peaks
-sigmay-sigmay of peaks
-ro-correlation coefficient
-txy, sxy, tx, sx-relative amplitudes
-bx, by-slopes


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)
AUXILIARY FUNCTION

This function calculates second derivative of peaks shape function
(see manual) according to sigmax of peaks.
Function parameters:
-numOfFittedPeaks-number of fitted peaks
-x,y-position of channel
-parameter-array of peaks parameters (amplitudes and positions)
-sigmax-sigmax of peaks
-sigmay-sigmay of peaks
-ro-correlation coefficient


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)
AUXILIARY FUNCTION

This function calculates derivative of peaks shape function (see manual)
according to sigmax of peaks.
Function parameters:
-numOfFittedPeaks-number of fitted peaks
-x,y-position of channel
-parameter-array of peaks parameters (amplitudes and positions)
-sigmax-sigmax of peaks
-sigmay-sigmay of peaks
-ro-correlation coefficient
-txy, sxy, ty, sy-relative amplitudes
-bx, by-slopes


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)
AUXILIARY FUNCTION

This function calculates second derivative of peaks shape function
(see manual) according to sigmay of peaks.
Function parameters:
-numOfFittedPeaks-number of fitted peaks
-x,y-position of channel
-parameter-array of peaks parameters (amplitudes and positions)
-sigmax-sigmax of peaks
-sigmay-sigmay of peaks
-ro-correlation coefficient


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)
AUXILIARY FUNCTION

This function calculates derivative of peaks shape function (see manual)
according to correlation coefficient ro.
Function parameters:
-numOfFittedPeaks-number of fitted peaks
-x,y-position of channel
-parameter-array of peaks parameters (amplitudes and positions)
-sx-sigmax of peaks
-sy-sigmay of peaks
-r-correlation coefficient ro


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)
AUXILIARY FUNCTION

This function calculates derivative of peaks shape function (see manual)
according to relative amplitude txy.
Function parameters:
-numOfFittedPeaks-number of fitted peaks
-x,y-position of channel
-parameter-array of peaks parameters (amplitudes and positions)
-sigmax-sigmax of peaks
-sigmay-sigmay of peaks
-bx, by-slopes


Double_t Dersxy(Int_t numOfFittedPeaks, Double_t x, Double_t y, const Double_t* parameter, Double_t sigmax, Double_t sigmay)
AUXILIARY FUNCTION

This function calculates derivative of peaks shape function (see manual)
according to relative amplitude sxy.
Function parameters:
-numOfFittedPeaks-number of fitted peaks
-x,y-position of channel
-parameter-array of peaks parameters (amplitudes and positions)
-sigmax-sigmax of peaks
-sigmay-sigmay of peaks


Double_t Dertx(Int_t numOfFittedPeaks, Double_t x, const Double_t* parameter, Double_t sigmax, Double_t bx)
AUXILIARY FUNCTION

This function calculates derivative of peaks shape function (see manual)
according to relative amplitude tx.
Function parameters:
-numOfFittedPeaks-number of fitted peaks
-x-position of channel
-parameter-array of peaks parameters (amplitudes and positions)
-sigmax-sigma of 1D ridge
-bx-slope


Double_t Derty(Int_t numOfFittedPeaks, Double_t x, const Double_t* parameter, Double_t sigmax, Double_t bx)
AUXILIARY FUNCTION

This function calculates derivative of peaks shape function (see manual)
according to relative amplitude ty.
Function parameters:
-numOfFittedPeaks-number of fitted peaks
-x-position of channel
-parameter-array of peaks parameters (amplitudes and positions)
-sigmax-sigma of 1D ridge
-bx-slope


Double_t Dersx(Int_t numOfFittedPeaks, Double_t x, const Double_t* parameter, Double_t sigmax)
AUXILIARY FUNCTION

This function calculates derivative of peaks shape function (see manual)
according to relative amplitude sx.
Function parameters:
-numOfFittedPeaks-number of fitted peaks
-x-position of channel
-parameter-array of peaks parameters (amplitudes and positions)
-sigmax-sigma of 1D ridge


Double_t Dersy(Int_t numOfFittedPeaks, Double_t x, const Double_t* parameter, Double_t sigmax)
AUXILIARY FUNCTION

This function calculates derivative of peaks shape function (see manual)
according to relative amplitude sy.
Function parameters:
-numOfFittedPeaks-number of fitted peaks
-x-position of channel
-parameter-array of peaks parameters (amplitudes and positions)
-sigmax-sigma of 1D ridge


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)
AUXILIARY FUNCTION

This function calculates derivative of peaks shape function (see manual)
according to slope bx.
Function parameters:
-numOfFittedPeaks-number of fitted peaks
-x,y-position of channel
-parameter-array of peaks parameters (amplitudes and positions)
-sigmax-sigmax of peaks
-sigmay-sigmay of peaks
-txy, tx-relative amplitudes
-bx, by-slopes


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)
AUXILIARY FUNCTION

This function calculates derivative of peaks shape function (see manual)
according to slope by.
Function parameters:
-numOfFittedPeaks-number of fitted peaks
-x,y-position of channel
-parameter-array of peaks parameters (amplitudes and positions)
-sigmax-sigmax of peaks
-sigmay-sigmay of peaks
-txy, ty-relative amplitudes
-bx, by-slopes


Double_t Volume(Double_t a, Double_t sx, Double_t sy, Double_t ro)
AUXILIARY FUNCTION

This function calculates volume of a peak
Function parameters:
-a-amplitude of the peak
-sx,sy-sigmas of peak
-ro-correlation coefficient


Double_t Derpa2(Double_t sx, Double_t sy, Double_t ro)
AUXILIARY FUNCTION

This function calculates derivative of the volume of a peak
according to amplitute
Function parameters:
-sx,sy-sigmas of peak
-ro-correlation coefficient


Double_t Derpsigmax(Double_t a, Double_t sy, Double_t ro)
AUXILIARY FUNCTION

This function calculates derivative of the volume of a peak
according to sigmax
Function parameters:
-a-amplitude of peak
-sy-sigma of peak
-ro-correlation coefficient


Double_t Derpsigmay(Double_t a, Double_t sx, Double_t ro)
AUXILIARY FUNCTION

This function calculates derivative of the volume of a peak
according to sigmay
Function parameters:
-a-amplitude of peak
-sx-sigma of peak
-ro-correlation coefficient


Double_t Derpro(Double_t a, Double_t sx, Double_t sy, Double_t ro)
AUXILIARY FUNCTION

This function calculates derivative of the volume of a peak
according to ro
Function parameters:
-a-amplitude of peak
-sx,sy-sigmas of peak
-ro-correlation coefficient


void FitAwmi(Float_t** source)
 TWO-DIMENSIONAL FIT FUNCTION
  ALGORITHM WITHOUT MATRIX INVERSION
  This function fits the source spectrum. The calling program should
  fill in input parameters of the TSpectrum2Fit class.
  The fitted parameters are written into
  TSpectrum2Fit class output parameters and fitted data are written into
  source spectrum.

        Function parameters:
        source-pointer to the matrix of source spectrum




Fitting

 

Goal: to estimate simultaneously peak shape parameters in spectra with large number of peaks

 

         peaks can be fitted separately, each peak (or multiplets) in a region or together all peaks in a spectrum. To fit separately each peak one needs to determine the fitted region. However it can happen that the regions of neighboring peaks are overlapping. Then the results of fitting are very poor. On the other hand, when fitting together all peaks found in a  spectrum, one needs to have a method that is  stable (converges) and fast enough to carry out fitting in reasonable time

         we have implemented the nonsymmetrical semiempirical peak shape function

         it contains the two-dimensional symmetrical Gaussian two one-dimensional symmetrical Gaussian ridges as well as nonsymmetrical terms and background.

where Txy, Tx, Ty, Sxy, Sx, Sy are relative amplitudes and Bx, By are slopes.

 

         algorithm without matrix inversion (AWMI) allows fitting tens, hundreds of peaks simultaneously that represent sometimes thousands of parameters [2], [5].

 

 

Function:

void TSpectrumFit2::FitAwmi(float **fSource)

 

This function fits the source spectrum using AWMI algorithm. The calling program should fill in input parameters of the TSpectrumFit2 class using a set of TSpectrumFit2 setters. The fitted parameters are written into the class and fitted data are written into source spectrum.

 

 

Parameter:

        fSource-pointer to the matrix of source spectrum                 

 

 

Member variables of  TSpectrumFit2 class:

   Int_t     fNPeaks;                        //number of peaks present in fit, input parameter, it should be > 0

   Int_t     fNumberIterations;              //number of iterations in fitting procedure, input parameter, it should be > 0

   Int_t     fXmin;                          //first fitted channel in x direction

   Int_t     fXmax;                          //last fitted channel in x direction

   Int_t     fYmin;                          //first fitted channel in y direction

   Int_t     fYmax;                          //last fitted channel in y direction

   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

   Int_t     fAlphaOptim;                    optimization of convergence algorithm, possible values kFitAlphaHalving, kFitAlphaOptimal

   Int_t     fPower;                         //possible values kFitPower2,4,6,8,10,12, for details see references. It applies only for Awmi fitting function.

   Int_t     fFitTaylor;                     //order of Taylor expansion, possible values kFitTaylorOrderFirst, kFitTaylorOrderSecond. It applies only for Awmi fitting function.

   Double_t  fAlpha;                         //convergence coefficient, input parameter, it should be positive number and <=1, for details see references

   Double_t  fChi;                           //here the fitting functions return resulting chi square  

   Double_t *fPositionInitX;                 //[fNPeaks] array of initial values of x positions of 2D peaks, input parameters

   Double_t *fPositionCalcX;                 //[fNPeaks] array of calculated values of x positions of 2D peaks, output parameters

   Double_t *fPositionErrX;                  //[fNPeaks] array of error values of x positions of 2D peaks, output parameters

   Double_t *fPositionInitY;                 //[fNPeaks] array of initial values of y positions of 2D peaks, input parameters

   Double_t *fPositionCalcY;                 //[fNPeaks] array of calculated values of y positions of 2D peaks, output parameters

   Double_t *fPositionErrY;                  //[fNPeaks] array of error values of y positions of 2D peaks, output parameters

   Double_t *fPositionInitX1;                //[fNPeaks] array of initial x positions of 1D ridges, input parameters

   Double_t *fPositionCalcX1;                //[fNPeaks] array of calculated x positions of 1D ridges, output parameters

   Double_t *fPositionErrX1;                 //[fNPeaks] array of x positions errors of 1D ridges, output parameters

   Double_t *fPositionInitY1;                //[fNPeaks] array of initial y positions of 1D ridges, input parameters

   Double_t *fPositionCalcY1;                //[fNPeaks] array of calculated y positions of 1D ridges, output parameters

   Double_t *fPositionErrY1;                 //[fNPeaks] array of y positions errors of 1D ridges, output parameters

   Double_t *fAmpInit;                       //[fNPeaks] array of initial values of amplitudes of 2D peaks, input parameters

   Double_t *fAmpCalc;                       //[fNPeaks] array of calculated values of amplitudes of 2D peaks, output parameters

   Double_t *fAmpErr;                        //[fNPeaks] array of amplitudes errors of 2D peaks, output parameters

   Double_t *fAmpInitX1;                     //[fNPeaks] array of initial values of amplitudes of 1D ridges in x direction, input parameters

   Double_t *fAmpCalcX1;                     //[fNPeaks] array of calculated values of amplitudes of 1D ridges in x direction, output parameters

   Double_t *fAmpErrX1;                      //[fNPeaks] array of amplitudes errors of 1D ridges in x direction, output parameters

   Double_t *fAmpInitY1;                     //[fNPeaks] array of initial values of amplitudes of 1D ridges in y direction, input parameters

   Double_t *fAmpCalcY1;                     //[fNPeaks] array of calculated values of amplitudes of 1D ridges in y direction, output parameters

   Double_t *fAmpErrY1;                      //[fNPeaks] array of amplitudes errors of 1D ridges in y direction, output parameters

   Double_t *fVolume;                        //[fNPeaks] array of calculated volumes of 2D peaks, output parameters

   Double_t *fVolumeErr;                     //[fNPeaks] array of volumes errors of 2D peaks, output parameters

   Double_t  fSigmaInitX;                    //initial value of sigma x parameter

   Double_t  fSigmaCalcX;                    //calculated value of sigma x parameter

   Double_t  fSigmaErrX;                     //error value of sigma x parameter

   Double_t  fSigmaInitY;                    //initial value of sigma y parameter

   Double_t  fSigmaCalcY;                    //calculated value of sigma y parameter

   Double_t  fSigmaErrY;                     //error value of sigma y parameter

   Double_t  fRoInit;                        //initial value of correlation coefficient

   Double_t  fRoCalc;                        //calculated value of correlation coefficient

   Double_t  fRoErr;                         //error value of correlation coefficient

   Double_t  fTxyInit;                       //initial value of t parameter for 2D peaks (relative amplitude of tail), for details see html manual and references

   Double_t  fTxyCalc;                       //calculated value of t parameter for 2D peaks

   Double_t  fTxyErr;                        //error value of t parameter for 2D peaks

   Double_t  fSxyInit;                       //initial value of s parameter for 2D peaks (relative amplitude of step), for details see html manual and references

   Double_t  fSxyCalc;                       //calculated value of s parameter for 2D peaks

   Double_t  fSxyErr;                        //error value of s parameter for 2D peaks

   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

   Double_t  fTxCalc;                        //calculated value of t parameter for 1D ridges in x direction

   Double_t  fTxErr;                         //error value of t parameter for 1D ridges in x direction

   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

   Double_t  fTyCalc;                        calculated value of t parameter for 1D ridges in y direction

   Double_t  fTyErr;                         //error value of t parameter for 1D ridges in y direction

   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

   Double_t  fSxCalc;                        //calculated value of s parameter for 1D ridges in x direction

   Double_t  fSxErr;                         //error value of s parameter for 1D ridges in x direction

   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

   Double_t  fSyCalc;                        //calculated value of s parameter for 1D ridges in y direction

   Double_t  fSyErr;                         //error value of s parameter for 1D ridges in y direction

   Double_t  fBxInit;                        //initial value of b parameter for 1D ridges in x direction (slope), for details see html manual and references

   Double_t  fBxCalc;                        //calculated value of b parameter for 1D ridges in x direction

   Double_t  fBxErr;                         //error value of b parameter for 1D ridges in x direction

   Double_t  fByInit;                        //initial value of b parameter for 1D ridges in y direction (slope), for details see html manual and references

   Double_t  fByCalc;                        //calculated value of b parameter for 1D ridges in y direction

   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  fA0Calc;                        //calculated value of background a0 parameter

   Double_t  fA0Err;                         //error value of background a0 parameter

   Double_t  fAxInit;                        //initial value of background ax parameter(backgroud is estimated as a0+ax*x+ay*y)

   Double_t  fAxCalc;                        //calculated value of background ax parameter

   Double_t  fAxErr;                         //error value of background ax parameter

   Double_t  fAyInit;                        //initial value of background ay parameter(backgroud is estimated as a0+ax*x+ay*y)

   Double_t  fAyCalc;                        //calculated value of background ay parameter

   Double_t  fAyErr;                         error value of background ay parameter  

   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

   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

   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

   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

   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

   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

   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

   Bool_t    fFixSigmaX;                     //logical value of sigma x parameter, which allows to fix the parameter (not to fit).

   Bool_t    fFixSigmaY;                     //logical value of sigma y parameter, which allows to fix the parameter (not to fit).

   Bool_t    fFixRo;                         //logical value of correlation coefficient, which allows to fix the parameter (not to fit).

   Bool_t    fFixTxy;                        //logical value of t parameter for 2D peaks, which allows to fix the parameter (not to fit).

   Bool_t    fFixSxy;                        //logical value of s parameter for 2D peaks, which allows to fix the parameter (not to fit).

   Bool_t    fFixTx;                         //logical value of t parameter for 1D ridges in x direction, which allows to fix the parameter (not to fit).

   Bool_t    fFixTy;                         //logical value of t parameter for 1D ridges in y direction, which allows to fix the parameter (not to fit).

   Bool_t    fFixSx;                         //logical value of s parameter for 1D ridges in x direction, which allows to fix the parameter (not to fit).

   Bool_t    fFixSy;                         //logical value of s parameter for 1D ridges in y direction, which allows to fix the parameter (not to fit).

   Bool_t    fFixBx;                         //logical value of b parameter for 1D ridges in x direction, which allows to fix the parameter (not to fit).

   Bool_t    fFixBy;                         //logical value of b parameter for 1D ridges in y direction, which allows to fix the parameter (not to fit).

   Bool_t    fFixA0;                         //logical value of a0 parameter, which allows to fix the parameter (not to fit).

   Bool_t    fFixAx;                         //logical value of ax parameter, which allows to fix the parameter (not to fit).

   Bool_t    fFixAy;                         //logical value of ay parameter, which allows to fix the parameter (not to fit).

 

References:

[1] Phillps G.W., Marlow K.W., NIM 137 (1976) 525.

[2] I. A. Slavic: Nonlinear least-squares fitting without matrix inversion applied to complex Gaussian spectra analysis. NIM 134 (1976) 285-289.

[3] T. Awaya: A new method for curve fitting to the data with low statistics not using chi-square method. NIM 165 (1979) 317-323.

[4] T. Hauschild, M. Jentschel: Comparison of maximum likelihood estimation and chi-square statistics applied to counting experiments. NIM A 457 (2001) 384-401.

 [5]  M. Morháč,  J. Kliman,  M. Jandel,  Ľ. Krupa, V. Matoušek: Study of fitting algorithms applied to simultaneous analysis of large number of peaks in -ray spectra. Applied Spectroscopy, Vol. 57, No. 7, pp. 753-760, 2003

 

Example  – script FitAwmi2.c:

Fig. 1 Original two-dimensional spectrum with found peaks (using TSpectrum2 peak searching function). The positions of peaks were used as initial estimates in fitting procedure.

Fig. 2 Fitted (generated from fitted parameters) spectrum of the data from Fig. 1 using Algorithm Without Matrix Inversion. Each peak was represented by 7 parameters, which together with Sigmax, Sigmay and a0 resulted in 38 fitted parameters. The chi-square after 1000 iterations was 0.642342.

 

Script:

 

// Example to illustrate fitting function, algorithm without matrix inversion (AWMI) (class TSpectrumFit2).

// To execute this example, do

// root > .x FitAwmi2.C

 

void FitAwmi2() {

   Int_t i, j, nfound;

   Int_t nbinsx = 64;

   Int_t nbinsy = 64;  

   Int_t xmin  = 0;

   Int_t xmax  = nbinsx;

   Int_t ymin  = 0;

   Int_t ymax  = nbinsy;

   Float_t ** source = new float *[nbinsx];  

   Float_t ** dest = new float *[nbinsx];     

   for (i=0;i<nbinsx;i++)

                                                source[i]=new float[nbinsy];

   for (i=0;i<nbinsx;i++)

                                                dest[i]=new float[nbinsy];

   TH2F *search = new TH2F("search","High resolution peak searching",nbinsx,xmin,xmax,nbinsy,ymin,ymax);

   TFile *f = new TFile("TSpectrum2.root");

   search=(TH2F*) f->Get("search4;1");

   TCanvas *Searching = new TCanvas("Searching","Two-dimensional fitting using Algorithm Without Matrix Inversion",10,10,1000,700);

   TSpectrum2 *s = new TSpectrum2();

   for (i = 0; i < nbinsx; i++){

     for (j = 0; j < nbinsy; j++){

                    source[i][j] = search->GetBinContent(i + 1,j + 1);

                 }

   }

   //searching for candidate peaks positions    

   nfound = s->SearchHighRes(source, dest, nbinsx, nbinsy, 2, 5, kTRUE, 3, kFALSE, 3);  

   Bool_t *FixPosX = new Bool_t[nfound];

   Bool_t *FixPosY = new Bool_t[nfound];  

   Bool_t *FixAmp = new Bool_t[nfound];     

   Float_t *PosX = new Float_t[nfound];        

   Float_t *PosY = new Float_t[nfound];

   Float_t *Amp = new Float_t[nfound];     

   Float_t *AmpXY = new Float_t[nfound];        

   PosX = s->GetPositionX();

   PosY = s->GetPositionY();      

   printf("Found %d candidate peaks\n",nfound);  

   for(i = 0; i< nfound ; i++){

      FixPosX[i] = kFALSE;

      FixPosY[i] = kFALSE;     

      FixAmp[i] = kFALSE;   

      Amp[i] = source[(int)(PosX[i]+0.5)][(int)(PosY[i]+0.5)];      //initial values of peaks amplitudes, input parameters         

      AmpXY[i] = 0;

   }

   //filling in the initial estimates of the input parameters

   TSpectrumFit2 *pfit=new TSpectrumFit2(nfound);

   pfit->SetFitParameters(xmin, xmax-1, ymin, ymax-1, 1000, 0.1, pfit->kFitOptimChiCounts, pfit->kFitAlphaHalving, pfit->kFitPower2, pfit->kFitTaylorOrderFirst);  

   pfit->SetPeakParameters(2, kFALSE, 2, kFALSE, 0, kTRUE, PosX, (Bool_t *) FixPosX, PosY, (Bool_t *) FixPosY, PosX, (Bool_t *) FixPosX, PosY, (Bool_t *) FixPosY, Amp, (Bool_t *) FixAmp, AmpXY, (Bool_t *) FixAmp, AmpXY, (Bool_t *) FixAmp);     

   pfit->SetBackgroundParameters(0, kFALSE, 0, kTRUE, 0, kTRUE);  

   pfit->FitAwmi(source);

    for (i = 0; i < nbinsx; i++){

     for (j = 0; j < nbinsy; j++){

                  search->SetBinContent(i + 1, j + 1,source[i][j]);

                 }

   }  

   search->Draw("SURF");

}

 

Example 2 – script FitA2.c:

Fig. 3 Original two-dimensional gamma-gamma-ray spectrum with found peaks (using TSpectrum2 peak searching function).

Fig. 4 Fitted (generated from fitted parameters) spectrum of the data from Fig. 3 using Algorithm Without Matrix Inversion. 152 peaks were identified. Each peak was represented by 7 parameters, which together with Sigmax, Sigmay and a0 resulted in 1067 fitted parameters. The chi-square after 1000 iterations was 0.728675. One can observe good correspondence with the original data.

 

Script:

// Example to illustrate fitting function, algorithm without matrix inversion (AWMI) (class TSpectrumFit2).

// To execute this example, do

// root > .x FitA2.C

void FitA2() {

   Int_t i, j, nfound;

   Int_t nbinsx = 256;

   Int_t nbinsy = 256;  

   Int_t xmin  = 0;

   Int_t xmax  = nbinsx;

   Int_t ymin  = 0;

   Int_t ymax  = nbinsy;

   Float_t ** source = new float *[nbinsx];  

   Float_t ** dest = new float *[nbinsx];     

   for (i=0;i<nbinsx;i++)

                                                source[i]=new float[nbinsy];

   for (i=0;i<nbinsx;i++)

                                                dest[i]=new float[nbinsy];  

   TH2F *search = new TH2F("search","High resolution peak searching",nbinsx,xmin,xmax,nbinsy,ymin,ymax);

   TFile *f = new TFile("TSpectrum2.root");

   search=(TH2F*) f->Get("fit1;1");

   TCanvas *Searching = new TCanvas("Searching","Two-dimensional fitting using Algorithm Without Matrix Inversion",10,10,1000,700);

   TSpectrum2 *s = new TSpectrum2(1000,1);

   for (i = 0; i < nbinsx; i++){

     for (j = 0; j < nbinsy; j++){

                    source[i][j] = search->GetBinContent(i + 1,j + 1);

                 }

   }  

   nfound = s->SearchHighRes(source, dest, nbinsx, nbinsy, 2, 2, kTRUE, 100, kFALSE, 3);  

   printf("Found %d candidate peaks\n",nfound);

   Bool_t *FixPosX = new Bool_t[nfound];

   Bool_t *FixPosY = new Bool_t[nfound];  

   Bool_t *FixAmp = new Bool_t[nfound];     

   Float_t *PosX = new Float_t[nfound];        

   Float_t *PosY = new Float_t[nfound];

   Float_t *Amp = new Float_t[nfound];     

   Float_t *AmpXY = new Float_t[nfound];        

   PosX = s->GetPositionX();

   PosY = s->GetPositionY();     

   for(i = 0; i< nfound ; i++){

      FixPosX[i] = kFALSE;

      FixPosY[i] = kFALSE;     

      FixAmp[i] = kFALSE;   

      Amp[i] = source[(int)(PosX[i]+0.5)][(int)(PosY[i]+0.5)];      //initial values of peaks amplitudes, input parameters         

      AmpXY[i] = 0;

   }

   filling in the initial estimates of the input parameters

   TSpectrumFit2 *pfit=new TSpectrumFit2(nfound);

   pfit->SetFitParameters(xmin, xmax-1, ymin, ymax-1, 1000, 0.1, pfit->kFitOptimChiCounts, pfit->kFitAlphaHalving, pfit->kFitPower2, pfit->kFitTaylorOrderFirst);  

   pfit->SetPeakParameters(2, kFALSE, 2, kFALSE, 0, kTRUE, PosX, (Bool_t *) FixPosX, PosY, (Bool_t *) FixPosY, PosX, (Bool_t *) FixPosX, PosY, (Bool_t *) FixPosY, Amp, (Bool_t *) FixAmp, AmpXY, (Bool_t *) FixAmp, AmpXY, (Bool_t *) FixAmp);     

   pfit->SetBackgroundParameters(0, kFALSE, 0, kTRUE, 0, kTRUE);  

   pfit->FitAwmi(source);

   for (i = 0; i < nbinsx; i++){

     for (j = 0; j < nbinsy; j++){

                  search->SetBinContent(i + 1, j + 1,source[i][j]);

                 }

   }  

   search->Draw("SURF");  

}

 

void FitStiefel(Float_t** source)
 TWO-DIMENSIONAL FIT FUNCTION USING STIEFEL-HESTENS
  ALGORITHM.
  This function fits the source spectrum. The calling program should
  fill in input parameters of the TSpectrum2Fit class.
  The fitted parameters are written into
  TSpectrum2Fit class output parameters and fitted data are written into
  source spectrum.

        Function parameters:
        source-pointer to the matrix of source spectrum




Stiefel fitting algorithm

 

Function:

void TSpectrumFit2::FitStiefel(float **fSource)

This function fits the source spectrum using Stiefel-Hestens method [1].  The calling program should fill in input fitting parameters of the TSpectrumFit2 class using a set of TSpectrumFit2 setters. The fitted parameters are written into the class and the fitted data are written into source spectrum. It converges faster than Awmi method.

 

Parameter:

        fSource-pointer to the matrix of source spectrum                 

       

Reference:

[1] B. Mihaila: Analysis of complex gamma spectra, Rom. Jorn. Phys., Vol. 39, No. 2, (1994), 139-148.

 

Example 1 – script FitS.c:

Fig. 1 Original two-dimensional spectrum with found peaks (using TSpectrum2 peak searching function). The positions of peaks were used as initial estimates in fitting procedure.

Fig. 2 Fitted (generated from fitted parameters) spectrum of the data from Fig. 1 using Stiefel-Hestens method. Each peak was represented by 7 parameters, which together with Sigmax, Sigmay and a0 resulted in 38 fitted parameters. The chi-square after 1000 iterations was 0.642157.

 

Script:

// Example to illustrate fitting function, algorithm without matrix inversion (AWMI) (class TSpectrumFit2).

// To execute this example, do

// root > .x FitStiefel2.C

 

void FitStiefel2() {

   Int_t i, j, nfound;

   Int_t nbinsx = 64;

   Int_t nbinsy = 64;  

   Int_t xmin  = 0;

   Int_t xmax  = nbinsx;

   Int_t ymin  = 0;

   Int_t ymax  = nbinsy;

   Float_t ** source = new float *[nbinsx];  

   Float_t ** dest = new float *[nbinsx];     

   for (i=0;i<nbinsx;i++)

                                                source[i]=new float[nbinsy];

   for (i=0;i<nbinsx;i++)

                                                dest[i]=new float[nbinsy];

   TH2F *search = new TH2F("search","High resolution peak searching",nbinsx,xmin,xmax,nbinsy,ymin,ymax);

   TFile *f = new TFile("TSpectrum2.root");

   search=(TH2F*) f->Get("search4;1");

   TCanvas *Searching = new TCanvas("Searching","Two-dimensional fitting using Stiefel-Hestens method",10,10,1000,700);

   TSpectrum2 *s = new TSpectrum2();

   for (i = 0; i < nbinsx; i++){

     for (j = 0; j < nbinsy; j++){

                    source[i][j] = search->GetBinContent(i + 1,j + 1);

                 }

   }  

   nfound = s->SearchHighRes(source, dest, nbinsx, nbinsy, 2, 5, kTRUE, 3, kFALSE, 3);  

   printf("Found %d candidate peaks\n",nfound);

   Bool_t *FixPosX = new Bool_t[nfound];

   Bool_t *FixPosY = new Bool_t[nfound];   

   Bool_t *FixAmp = new Bool_t[nfound];     

   Float_t *PosX = new Float_t[nfound];        

   Float_t *PosY = new Float_t[nfound];

   Float_t *Amp = new Float_t[nfound];     

   Float_t *AmpXY = new Float_t[nfound];        

   PosX = s->GetPositionX();

   PosY = s->GetPositionY();     

   for(i = 0; i< nfound ; i++){

      FixPosX[i] = kFALSE;

      FixPosY[i] = kFALSE;     

      FixAmp[i] = kFALSE;   

      Amp[i] = source[(int)(PosX[i]+0.5)][(int)(PosY[i]+0.5)];      //initial values of peaks amplitudes, input parameters         

      AmpXY[i] = 0;

   }

   //filling in the initial estimates of the input parameters

   TSpectrumFit2 *pfit=new TSpectrumFit2(nfound);

   pfit->SetFitParameters(xmin, xmax-1, ymin, ymax-1, 1000, 0.1, pfit->kFitOptimChiCounts, pfit->kFitAlphaHalving, pfit->kFitPower2, pfit->kFitTaylorOrderFirst);  

   pfit->SetPeakParameters(2, kFALSE, 2, kFALSE, 0, kTRUE, PosX, (Bool_t *) FixPosX, PosY, (Bool_t *) FixPosY, PosX, (Bool_t *) FixPosX, PosY, (Bool_t *) FixPosY, Amp, (Bool_t *) FixAmp, AmpXY, (Bool_t *) FixAmp, AmpXY, (Bool_t *) FixAmp);     

   pfit->SetBackgroundParameters(0, kFALSE, 0, kTRUE, 0, kTRUE);  

   pfit->FitStiefel(source);

    for (i = 0; i < nbinsx; i++){

     for (j = 0; j < nbinsy; j++){

                  search->SetBinContent(i + 1, j + 1,source[i][j]);

                 }

   }  

   search->Draw("SURF");

}

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)
   SETTER FUNCTION

   This function sets the following fitting parameters:
         -xmin, xmax, ymin, ymax - fitting region
         -numberIterations - # of desired iterations in the fit
         -alpha - convergence coefficient, it should be positive number and <=1, for details see references
         -statisticType - 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
         -alphaOptim - optimization of convergence algorithm, possible values kFitAlphaHalving, kFitAlphaOptimal
         -power - possible values kFitPower2,4,6,8,10,12, for details see references. It applies only for Awmi fitting function.
         -fitTaylor - order of Taylor expansion, possible values kFitTaylorOrderFirst, kFitTaylorOrderSecond. It applies only for Awmi fitting function.

void SetPeakParameters(Double_t sigmaX, Bool_t fixSigmaX, Double_t sigmaY, Bool_t fixSigmaY, Double_t ro, Bool_t fixRo, const Float_t* positionInitX, const Bool_t* fixPositionX, const Float_t* positionInitY, const Bool_t* fixPositionY, const Float_t* positionInitX1, const Bool_t* fixPositionX1, const Float_t* positionInitY1, const Bool_t* fixPositionY1, const Float_t* ampInit, const Bool_t* fixAmp, const Float_t* ampInitX1, const Bool_t* fixAmpX1, const Float_t* ampInitY1, const Bool_t* fixAmpY1)
   SETTER FUNCTION

   This function sets the following fitting parameters of peaks:
         -sigmaX - initial value of sigma x parameter
         -fixSigmaX - logical value of sigma x parameter, which allows to fix the parameter (not to fit)
         -sigmaY - initial value of sigma y parameter
         -fixSigmaY - logical value of sigma y parameter, which allows to fix the parameter (not to fit)
         -ro - initial value of ro parameter (correlation coefficient)
         -fixRo - logical value of ro parameter, which allows to fix the parameter (not to fit)
         -positionInitX - aray of initial values of peaks x positions
         -fixPositionX - array of logical values which allow to fix appropriate x positions (not fit). However they are present in the estimated functional.
         -positionInitY - aray of initial values of peaks y positions
         -fixPositionY - array of logical values which allow to fix appropriate y positions (not fit). However they are present in the estimated functional.
         -ampInit - aray of initial values of  2D peaks amplitudes
         -fixAmp - aray of logical values which allow to fix appropriate amplitudes of 2D peaks (not fit). However they are present in the estimated functional
         -ampInitX1 - aray of initial values of amplitudes of  1D ridges in x direction
         -fixAmpX1 - aray 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
         -ampInitY1 - aray of initial values of amplitudes of  1D ridges in y direction
         -fixAmpY1 - aray 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

void SetBackgroundParameters(Double_t a0Init, Bool_t fixA0, Double_t axInit, Bool_t fixAx, Double_t ayInit, Bool_t fixAy)
   SETTER FUNCTION

   This function sets the following fitting parameters of background:
         -a0Init - initial value of a0 parameter (backgroud is estimated as a0+ax*x+ay*y)
         -fixA0 - logical value of a0 parameter, which allows to fix the parameter (not to fit)
         -axInit - initial value of ax parameter
         -fixAx - logical value of ax parameter, which allows to fix the parameter (not to fit)
         -ayInit - initial value of ay parameter
         -fixAy - logical value of ay parameter, which allows to fix the parameter (not to fit)

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)
   SETTER FUNCTION

   This function sets the following fitting parameters of tails of peaks
         -tInitXY - initial value of txy parameter
         -fixTxy - logical value of txy parameter, which allows to fix the parameter (not to fit)
         -tInitX - initial value of tx parameter
         -fixTx - logical value of tx parameter, which allows to fix the parameter (not to fit)
         -tInitY - initial value of ty parameter
         -fixTy - logical value of ty parameter, which allows to fix the parameter (not to fit)
         -bInitX - initial value of bx parameter
         -fixBx - logical value of bx parameter, which allows to fix the parameter (not to fit)
         -bInitY - initial value of by parameter
         -fixBy - logical value of by parameter, which allows to fix the parameter (not to fit)
         -sInitXY - initial value of sxy parameter
         -fixSxy - logical value of sxy parameter, which allows to fix the parameter (not to fit)
         -sInitX - initial value of sx parameter
         -fixSx - logical value of sx parameter, which allows to fix the parameter (not to fit)
         -sInitY - initial value of sy parameter
         -fixSy - logical value of sy parameter, which allows to fix the parameter (not to fit)

void GetPositions(Float_t* positionsX, Float_t* positionsY, Float_t* positionsX1, Float_t* positionsY1)
   GETTER FUNCTION

   This function gets the positions of fitted 2D peaks and 1D ridges
         -positionX - gets vector of x positions of 2D peaks
         -positionY - gets vector of y positions of 2D peaks
         -positionX1 - gets vector of x positions of 1D ridges
         -positionY1 - gets vector of y positions of 1D ridges

void GetPositionErrors(Float_t* positionErrorsX, Float_t* positionErrorsY, Float_t* positionErrorsX1, Float_t* positionErrorsY1)
   GETTER FUNCTION

   This function gets the errors of positions of fitted 2D peaks and 1D ridges
         -positionErrorsX - gets vector of errors of x positions of 2D peaks
         -positionErrorsY - gets vector of errors of y positions of 2D peaks
         -positionErrorsX1 - gets vector of errors of x positions of 1D ridges
         -positionErrorsY1 - gets vector of errors of y positions of 1D ridges

void GetAmplitudes(Float_t* amplitudes, Float_t* amplitudesX1, Float_t* amplitudesY1)
   GETTER FUNCTION

   This function gets the amplitudes of fitted 2D peaks and 1D ridges
         -amplitudes - gets vector of amplitudes of 2D peaks
         -amplitudesX1 - gets vector of amplitudes of 1D ridges in x direction
         -amplitudesY1 - gets vector of amplitudes of 1D ridges in y direction

void GetAmplitudeErrors(Float_t* amplitudeErrors, Float_t* amplitudeErrorsX1, Float_t* amplitudeErrorsY1)
   GETTER FUNCTION

   This function gets the amplitudes of fitted 2D peaks and 1D ridges
         -amplitudeErrors - gets vector of amplitudes errors of 2D peaks
         -amplitudeErrorsX1 - gets vector of amplitudes errors of 1D ridges in x direction
         -amplitudesErrorY1 - gets vector of amplitudes errors of 1D ridges in y direction

void GetVolumes(Float_t* volumes)
   GETTER FUNCTION

   This function gets the volumes of fitted 2D peaks
         -volumes - gets vector of volumes of 2D peaks

void GetVolumeErrors(Float_t* volumeErrors)
   GETTER FUNCTION

   This function gets errors of the volumes of fitted 2D peaks
         -volumeErrors - gets vector of volumes errors of 2D peaks

void GetSigmaX(Double_t& sigmaX, Double_t& sigmaErrX)
   GETTER FUNCTION

   This function gets the sigma x parameter and its error
         -sigmaX - gets the fitted value of sigma x parameter
         -sigmaErrX - gets error value of sigma x parameter

void GetSigmaY(Double_t& sigmaY, Double_t& sigmaErrY)
   GETTER FUNCTION

   This function gets the sigma y parameter and its error
         -sigmaY - gets the fitted value of sigma y parameter
         -sigmaErrY - gets error value of sigma y parameter

void GetRo(Double_t& ro, Double_t& roErr)
   GETTER FUNCTION

   This function gets the ro parameter and its error
         -ro - gets the fitted value of ro parameter
         -roErr - gets error value of ro parameter

void GetBackgroundParameters(Double_t& a0, Double_t& a0Err, Double_t& ax, Double_t& axErr, Double_t& ay, Double_t& ayErr)
   GETTER FUNCTION

   This function gets the background parameters and their errors
         -a0 - gets the fitted value of a0 parameter
         -a0Err - gets error value of a0 parameter
         -ax - gets the fitted value of ax parameter
         -axErr - gets error value of ax parameter
         -ay - gets the fitted value of ay parameter
         -ayErr - gets error value of ay parameter

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)
   GETTER FUNCTION

   This function gets the tail parameters and their errors
         -txy - gets the fitted value of txy parameter
         -txyErr - gets error value of txy parameter
         -tx - gets the fitted value of tx parameter
         -txErr - gets error value of tx parameter
         -ty - gets the fitted value of ty parameter
         -tyErr - gets error value of ty parameter
         -bx - gets the fitted value of bx parameter
         -bxErr - gets error value of bx parameter
         -by - gets the fitted value of by parameter
         -byErr - gets error value of by parameter
         -sxy - gets the fitted value of sxy parameter
         -sxyErr - gets error value of sxy parameter
         -sx - gets the fitted value of sx parameter
         -sxErr - gets error value of sx parameter
         -sy - gets the fitted value of sy parameter
         -syErr - gets error value of sy parameter

TSpectrum2Fit(void)
Double_t GetChi() const
{return fChi;}