Logo ROOT   6.12/07
Reference Guide
List of all members | Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members

Advanced 2-dimensional spectra fitting functions.

Author
Miroslav Morhac

Class for fitting 2D spectra using AWMI (algorithm without matrix inversion) and conjugate gradient algorithms for symmetrical matrices (Stiefel-Hestens method). AWMI method allows to fit simultaneously 100s up to 1000s peaks. Stiefel method is very stable, it converges faster, but is more time consuming.

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.

Definition at line 16 of file TSpectrum2Fit.h.

Public Types

enum  {
  kFitOptimChiCounts =0, kFitOptimChiFuncValues =1, kFitOptimMaxLikelihood =2, kFitAlphaHalving =0,
  kFitAlphaOptimal =1, kFitPower2 =2, kFitPower4 =4, kFitPower6 =6,
  kFitPower8 =8, kFitPower10 =10, kFitPower12 =12, kFitTaylorOrderFirst =0,
  kFitTaylorOrderSecond =1, kFitNumRegulCycles =100
}
 
- Public Types inherited from TObject
enum  {
  kIsOnHeap = 0x01000000, kNotDeleted = 0x02000000, kZombie = 0x04000000, kInconsistent = 0x08000000,
  kBitMask = 0x00ffffff
}
 
enum  { kSingleKey = BIT(0), kOverwrite = BIT(1), kWriteDelete = BIT(2) }
 
enum  EDeprecatedStatusBits { kObjInCanvas = BIT(3) }
 
enum  EStatusBits {
  kCanDelete = BIT(0), kMustCleanup = BIT(3), kIsReferenced = BIT(4), kHasUUID = BIT(5),
  kCannotPick = BIT(6), kNoContextMenu = BIT(8), kInvalidObject = BIT(13)
}
 

Public Member Functions

 TSpectrum2Fit (void)
 Default constructor. More...
 
 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. More...
 
virtual ~TSpectrum2Fit ()
 Destructor. More...
 
void FitAwmi (Double_t **source)
 This function fits the source spectrum. More...
 
void FitStiefel (Double_t **source)
 This function fits the source spectrum. More...
 
void GetAmplitudeErrors (Double_t *amplitudeErrors, Double_t *amplitudeErrorsX1, Double_t *amplitudeErrorsY1)
 This function gets the amplitudes of fitted 2D peaks and 1D ridges. More...
 
void GetAmplitudes (Double_t *amplitudes, Double_t *amplitudesX1, Double_t *amplitudesY1)
 This function gets the amplitudes of fitted 2D peaks and 1D ridges. More...
 
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. More...
 
Double_t GetChi () const
 
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. More...
 
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. More...
 
void GetRo (Double_t &ro, Double_t &roErr)
 This function gets the ro parameter and its error. More...
 
void GetSigmaX (Double_t &sigmaX, Double_t &sigmaErrX)
 This function gets the sigma x parameter and its error. More...
 
void GetSigmaY (Double_t &sigmaY, Double_t &sigmaErrY)
 This function gets the sigma y parameter and its error. More...
 
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. More...
 
void GetVolumeErrors (Double_t *volumeErrors)
 This function gets errors of the volumes of fitted 2D peaks. More...
 
void GetVolumes (Double_t *volumes)
 This function gets the volumes of fitted 2D peaks. More...
 
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: More...
 
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: More...
 
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: More...
 
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. More...
 
- Public Member Functions inherited from TNamed
 TNamed ()
 
 TNamed (const char *name, const char *title)
 
 TNamed (const TString &name, const TString &title)
 
 TNamed (const TNamed &named)
 TNamed copy ctor. More...
 
virtual ~TNamed ()
 TNamed destructor. More...
 
virtual void Clear (Option_t *option="")
 Set name and title to empty strings (""). More...
 
virtual TObjectClone (const char *newname="") const
 Make a clone of an object using the Streamer facility. More...
 
virtual Int_t Compare (const TObject *obj) const
 Compare two TNamed objects. More...
 
virtual void Copy (TObject &named) const
 Copy this to obj. More...
 
virtual void FillBuffer (char *&buffer)
 Encode TNamed into output buffer. More...
 
virtual const char * GetName () const
 Returns name of object. More...
 
virtual const char * GetTitle () const
 Returns title of object. More...
 
virtual ULong_t Hash () const
 Return hash value for this object. More...
 
virtual Bool_t IsSortable () const
 
virtual void ls (Option_t *option="") const
 List TNamed name and title. More...
 
TNamedoperator= (const TNamed &rhs)
 TNamed assignment operator. More...
 
virtual void Print (Option_t *option="") const
 Print TNamed name and title. More...
 
virtual void SetName (const char *name)
 Set the name of the TNamed. More...
 
virtual void SetNameTitle (const char *name, const char *title)
 Set all the TNamed parameters (name and title). More...
 
virtual void SetTitle (const char *title="")
 Set the title of the TNamed. More...
 
virtual Int_t Sizeof () const
 Return size of the TNamed part of the TObject. More...
 
- Public Member Functions inherited from TObject
 TObject ()
 TObject constructor. More...
 
 TObject (const TObject &object)
 TObject copy ctor. More...
 
virtual ~TObject ()
 TObject destructor. More...
 
void AbstractMethod (const char *method) const
 Use this method to implement an "abstract" method that you don't want to leave purely abstract. More...
 
virtual void AppendPad (Option_t *option="")
 Append graphics object to current pad. More...
 
virtual void Browse (TBrowser *b)
 Browse object. May be overridden for another default action. More...
 
ULong_t CheckedHash ()
 Checked and record whether for this class has a consistent Hash/RecursiveRemove setup (*) and then return the regular Hash value for this object. More...
 
virtual const char * ClassName () const
 Returns name of class to which the object belongs. More...
 
virtual void Delete (Option_t *option="")
 Delete this object. More...
 
virtual Int_t DistancetoPrimitive (Int_t px, Int_t py)
 Computes distance from point (px,py) to the object. More...
 
virtual void Draw (Option_t *option="")
 Default Draw method for all objects. More...
 
virtual void DrawClass () const
 Draw class inheritance tree of the class to which this object belongs. More...
 
virtual TObjectDrawClone (Option_t *option="") const
 Draw a clone of this object in the current selected pad for instance with: gROOT->SetSelectedPad(gPad). More...
 
virtual void Dump () const
 Dump contents of object on stdout. More...
 
virtual void Error (const char *method, const char *msgfmt,...) const
 Issue error message. More...
 
virtual void Execute (const char *method, const char *params, Int_t *error=0)
 Execute method on this object with the given parameter string, e.g. More...
 
virtual void Execute (TMethod *method, TObjArray *params, Int_t *error=0)
 Execute method on this object with parameters stored in the TObjArray. More...
 
virtual void ExecuteEvent (Int_t event, Int_t px, Int_t py)
 Execute action corresponding to an event at (px,py). More...
 
virtual void Fatal (const char *method, const char *msgfmt,...) const
 Issue fatal error message. More...
 
virtual TObjectFindObject (const char *name) const
 Must be redefined in derived classes. More...
 
virtual TObjectFindObject (const TObject *obj) const
 Must be redefined in derived classes. More...
 
virtual Option_tGetDrawOption () const
 Get option used by the graphics system to draw this object. More...
 
virtual const char * GetIconName () const
 Returns mime type name of object. More...
 
virtual char * GetObjectInfo (Int_t px, Int_t py) const
 Returns string containing info about the object at position (px,py). More...
 
virtual Option_tGetOption () const
 
virtual UInt_t GetUniqueID () const
 Return the unique object id. More...
 
virtual Bool_t HandleTimer (TTimer *timer)
 Execute action in response of a timer timing out. More...
 
Bool_t HasInconsistentHash () const
 Return true is the type of this object is known to have an inconsistent setup for Hash and RecursiveRemove (i.e. More...
 
virtual void Info (const char *method, const char *msgfmt,...) const
 Issue info message. More...
 
virtual Bool_t InheritsFrom (const char *classname) const
 Returns kTRUE if object inherits from class "classname". More...
 
virtual Bool_t InheritsFrom (const TClass *cl) const
 Returns kTRUE if object inherits from TClass cl. More...
 
virtual void Inspect () const
 Dump contents of this object in a graphics canvas. More...
 
void InvertBit (UInt_t f)
 
virtual Bool_t IsEqual (const TObject *obj) const
 Default equal comparison (objects are equal if they have the same address in memory). More...
 
virtual Bool_t IsFolder () const
 Returns kTRUE in case object contains browsable objects (like containers or lists of other objects). More...
 
R__ALWAYS_INLINE Bool_t IsOnHeap () const
 
R__ALWAYS_INLINE Bool_t IsZombie () const
 
void MayNotUse (const char *method) const
 Use this method to signal that a method (defined in a base class) may not be called in a derived class (in principle against good design since a child class should not provide less functionality than its parent, however, sometimes it is necessary). More...
 
virtual Bool_t Notify ()
 This method must be overridden to handle object notification. More...
 
void Obsolete (const char *method, const char *asOfVers, const char *removedFromVers) const
 Use this method to declare a method obsolete. More...
 
void operator delete (void *ptr)
 Operator delete. More...
 
void operator delete[] (void *ptr)
 Operator delete []. More...
 
voidoperator new (size_t sz)
 
voidoperator new (size_t sz, void *vp)
 
voidoperator new[] (size_t sz)
 
voidoperator new[] (size_t sz, void *vp)
 
TObjectoperator= (const TObject &rhs)
 TObject assignment operator. More...
 
virtual void Paint (Option_t *option="")
 This method must be overridden if a class wants to paint itself. More...
 
virtual void Pop ()
 Pop on object drawn in a pad to the top of the display list. More...
 
virtual Int_t Read (const char *name)
 Read contents of object with specified name from the current directory. More...
 
virtual void RecursiveRemove (TObject *obj)
 Recursively remove this object from a list. More...
 
void ResetBit (UInt_t f)
 
virtual void SaveAs (const char *filename="", Option_t *option="") const
 Save this object in the file specified by filename. More...
 
virtual void SavePrimitive (std::ostream &out, Option_t *option="")
 Save a primitive as a C++ statement(s) on output stream "out". More...
 
void SetBit (UInt_t f, Bool_t set)
 Set or unset the user status bits as specified in f. More...
 
void SetBit (UInt_t f)
 
virtual void SetDrawOption (Option_t *option="")
 Set drawing option for object. More...
 
virtual void SetUniqueID (UInt_t uid)
 Set the unique object id. More...
 
virtual void SysError (const char *method, const char *msgfmt,...) const
 Issue system error message. More...
 
R__ALWAYS_INLINE Bool_t TestBit (UInt_t f) const
 
Int_t TestBits (UInt_t f) const
 
virtual void UseCurrentStyle ()
 Set current style settings in this object This function is called when either TCanvas::UseCurrentStyle or TROOT::ForceStyle have been invoked. More...
 
virtual void Warning (const char *method, const char *msgfmt,...) const
 Issue warning message. More...
 
virtual Int_t Write (const char *name=0, Int_t option=0, Int_t bufsize=0)
 Write this object to the current directory. More...
 
virtual Int_t Write (const char *name=0, Int_t option=0, Int_t bufsize=0) const
 Write this object to the current directory. More...
 

Protected Member Functions

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 2D peak. More...
 
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 the ridge. More...
 
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. More...
 
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. More...
 
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 position of 1D ridge. More...
 
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 position of 2D peak. More...
 
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 position of 2D peak. More...
 
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 of peaks. More...
 
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 of peaks. More...
 
Double_t Derfc (Double_t x)
 This function calculates derivative of error function of x. More...
 
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 of 1D ridge. More...
 
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 of 2D peak. More...
 
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 of 2D peak Function parameters: More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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 coefficient ro. More...
 
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. More...
 
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. More...
 
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 amplitude sx. More...
 
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 amplitude sxy. More...
 
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 amplitude sy. More...
 
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 amplitude tx. More...
 
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 amplitude txy. More...
 
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 amplitude ty. More...
 
Double_t Erfc (Double_t x)
 This function calculates error function of x. More...
 
Double_t Ourpowl (Double_t a, Int_t pw)
 power function More...
 
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) More...
 
void StiefelInversion (Double_t **a, Int_t size)
 This function calculates solution of the system of linear equations. More...
 
Double_t Volume (Double_t a, Double_t sx, Double_t sy, Double_t ro)
 This function calculates volume of a peak. More...
 
- Protected Member Functions inherited from TObject
virtual void DoError (int level, const char *location, const char *fmt, va_list va) const
 Interface to ErrorHandler (protected). More...
 
void MakeZombie ()
 

Protected Attributes

Double_t fA0Calc
 calculated value of background a0 parameter More...
 
Double_t fA0Err
 error value of background a0 parameter More...
 
Double_t fA0Init
 initial value of background a0 parameter(backgroud is estimated as a0+ax*x+ay*y) More...
 
Double_t fAlpha
 convergence coefficient, input parameter, it should be positive number and <=1, for details see references More...
 
Int_t fAlphaOptim
 optimization of convergence algorithm, possible values kFitAlphaHalving, kFitAlphaOptimal More...
 
Double_tfAmpCalc
 [fNPeaks] array of calculated values of amplitudes of 2D peaks, output parameters More...
 
Double_tfAmpCalcX1
 [fNPeaks] array of calculated values of amplitudes of 1D ridges in x direction, output parameters More...
 
Double_tfAmpCalcY1
 [fNPeaks] array of calculated values of amplitudes of 1D ridges in y direction, output parameters More...
 
Double_tfAmpErr
 [fNPeaks] array of amplitudes errors of 2D peaks, output parameters More...
 
Double_tfAmpErrX1
 [fNPeaks] array of amplitudes errors of 1D ridges in x direction, output parameters More...
 
Double_tfAmpErrY1
 [fNPeaks] array of amplitudes errors of 1D ridges in y direction, output parameters More...
 
Double_tfAmpInit
 [fNPeaks] array of initial values of amplitudes of 2D peaks, input parameters More...
 
Double_tfAmpInitX1
 [fNPeaks] array of initial values of amplitudes of 1D ridges in x direction, input parameters More...
 
Double_tfAmpInitY1
 [fNPeaks] array of initial values of amplitudes of 1D ridges in y direction, input parameters More...
 
Double_t fAxCalc
 calculated value of background ax parameter More...
 
Double_t fAxErr
 error value of background ax parameter More...
 
Double_t fAxInit
 initial value of background ax parameter(backgroud is estimated as a0+ax*x+ay*y) More...
 
Double_t fAyCalc
 calculated value of background ay parameter More...
 
Double_t fAyErr
 error value of background ay parameter More...
 
Double_t fAyInit
 initial value of background ay parameter(backgroud is estimated as a0+ax*x+ay*y) More...
 
Double_t fBxCalc
 calculated value of b parameter for 1D ridges in x direction More...
 
Double_t fBxErr
 error value of b parameter for 1D ridges in x direction More...
 
Double_t fBxInit
 initial value of b parameter for 1D ridges in x direction (slope), for details see html manual and references More...
 
Double_t fByCalc
 calculated value of b parameter for 1D ridges in y direction More...
 
Double_t fByErr
 error value of b parameter for 1D ridges in y direction More...
 
Double_t fByInit
 initial value of b parameter for 1D ridges in y direction (slope), for details see html manual and references More...
 
Double_t fChi
 here the fitting functions return resulting chi square More...
 
Int_t fFitTaylor
 order of Taylor expansion, possible values kFitTaylorOrderFirst, kFitTaylorOrderSecond. It applies only for Awmi fitting function. More...
 
Bool_t fFixA0
 logical value of a0 parameter, which allows to fix the parameter (not to fit). More...
 
Bool_tfFixAmp
 [fNPeaks] array of logical values which allow to fix appropriate amplitudes of 2D peaks (not fit). However they are present in the estimated functional More...
 
Bool_tfFixAmpX1
 [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 More...
 
Bool_tfFixAmpY1
 [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 More...
 
Bool_t fFixAx
 logical value of ax parameter, which allows to fix the parameter (not to fit). More...
 
Bool_t fFixAy
 logical value of ay parameter, which allows to fix the parameter (not to fit). More...
 
Bool_t fFixBx
 logical value of b parameter for 1D ridges in x direction, which allows to fix the parameter (not to fit). More...
 
Bool_t fFixBy
 logical value of b parameter for 1D ridges in y direction, which allows to fix the parameter (not to fit). More...
 
Bool_tfFixPositionX
 [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 More...
 
Bool_tfFixPositionX1
 [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 More...
 
Bool_tfFixPositionY
 [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 More...
 
Bool_tfFixPositionY1
 [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 More...
 
Bool_t fFixRo
 logical value of correlation coefficient, which allows to fix the parameter (not to fit). More...
 
Bool_t fFixSigmaX
 logical value of sigma x parameter, which allows to fix the parameter (not to fit). More...
 
Bool_t fFixSigmaY
 logical value of sigma y parameter, which allows to fix the parameter (not to fit). More...
 
Bool_t fFixSx
 logical value of s parameter for 1D ridges in x direction, which allows to fix the parameter (not to fit). More...
 
Bool_t fFixSxy
 logical value of s parameter for 2D peaks, which allows to fix the parameter (not to fit). More...
 
Bool_t fFixSy
 logical value of s parameter for 1D ridges in y direction, which allows to fix the parameter (not to fit). More...
 
Bool_t fFixTx
 logical value of t parameter for 1D ridges in x direction, which allows to fix the parameter (not to fit). More...
 
Bool_t fFixTxy
 logical value of t parameter for 2D peaks, which allows to fix the parameter (not to fit). More...
 
Bool_t fFixTy
 logical value of t parameter for 1D ridges in y direction, which allows to fix the parameter (not to fit). More...
 
Int_t fNPeaks
 number of peaks present in fit, input parameter, it should be > 0 More...
 
Int_t fNumberIterations
 number of iterations in fitting procedure, input parameter, it should be > 0 More...
 
Double_tfPositionCalcX
 [fNPeaks] array of calculated values of x positions of 2D peaks, output parameters More...
 
Double_tfPositionCalcX1
 [fNPeaks] array of calculated x positions of 1D ridges, output parameters More...
 
Double_tfPositionCalcY
 [fNPeaks] array of calculated values of y positions of 2D peaks, output parameters More...
 
Double_tfPositionCalcY1
 [fNPeaks] array of calculated y positions of 1D ridges, output parameters More...
 
Double_tfPositionErrX
 [fNPeaks] array of error values of x positions of 2D peaks, output parameters More...
 
Double_tfPositionErrX1
 [fNPeaks] array of x positions errors of 1D ridges, output parameters More...
 
Double_tfPositionErrY
 [fNPeaks] array of error values of y positions of 2D peaks, output parameters More...
 
Double_tfPositionErrY1
 [fNPeaks] array of y positions errors of 1D ridges, output parameters More...
 
Double_tfPositionInitX
 [fNPeaks] array of initial values of x positions of 2D peaks, input parameters More...
 
Double_tfPositionInitX1
 [fNPeaks] array of initial x positions of 1D ridges, input parameters More...
 
Double_tfPositionInitY
 [fNPeaks] array of initial values of y positions of 2D peaks, input parameters More...
 
Double_tfPositionInitY1
 [fNPeaks] array of initial y positions of 1D ridges, input parameters More...
 
Int_t fPower
 possible values kFitPower2,4,6,8,10,12, for details see references. It applies only for Awmi fitting function. More...
 
Double_t fRoCalc
 calculated value of correlation coefficient More...
 
Double_t fRoErr
 error value of correlation coefficient More...
 
Double_t fRoInit
 initial value of correlation coefficient More...
 
Double_t fSigmaCalcX
 calculated value of sigma x parameter More...
 
Double_t fSigmaCalcY
 calculated value of sigma y parameter More...
 
Double_t fSigmaErrX
 error value of sigma x parameter More...
 
Double_t fSigmaErrY
 error value of sigma y parameter More...
 
Double_t fSigmaInitX
 initial value of sigma x parameter More...
 
Double_t fSigmaInitY
 initial value of sigma y parameter More...
 
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 More...
 
Double_t fSxCalc
 calculated value of s parameter for 1D ridges in x direction More...
 
Double_t fSxErr
 error value of s parameter for 1D ridges in x direction More...
 
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 More...
 
Double_t fSxyCalc
 calculated value of s parameter for 2D peaks More...
 
Double_t fSxyErr
 error value of s parameter for 2D peaks More...
 
Double_t fSxyInit
 initial value of s parameter for 2D peaks (relative amplitude of step), for details see html manual and references More...
 
Double_t fSyCalc
 calculated value of s parameter for 1D ridges in y direction More...
 
Double_t fSyErr
 error value of s parameter for 1D ridges in y direction More...
 
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 More...
 
Double_t fTxCalc
 calculated value of t parameter for 1D ridges in x direction More...
 
Double_t fTxErr
 error value of t parameter for 1D ridges in x direction More...
 
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 More...
 
Double_t fTxyCalc
 calculated value of t parameter for 2D peaks More...
 
Double_t fTxyErr
 error value of t parameter for 2D peaks More...
 
Double_t fTxyInit
 initial value of t parameter for 2D peaks (relative amplitude of tail), for details see html manual and references More...
 
Double_t fTyCalc
 calculated value of t parameter for 1D ridges in y direction More...
 
Double_t fTyErr
 error value of t parameter for 1D ridges in y direction More...
 
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 More...
 
Double_tfVolume
 [fNPeaks] array of calculated volumes of 2D peaks, output parameters More...
 
Double_tfVolumeErr
 [fNPeaks] array of volumes errors of 2D peaks, output parameters More...
 
Int_t fXmax
 last fitted channel in x direction More...
 
Int_t fXmin
 first fitted channel in x direction More...
 
Int_t fYmax
 last fitted channel in y direction More...
 
Int_t fYmin
 first fitted channel in y direction More...
 
- Protected Attributes inherited from TNamed
TString fName
 
TString fTitle
 

Additional Inherited Members

- Static Public Member Functions inherited from TObject
static Long_t GetDtorOnly ()
 Return destructor only flag. More...
 
static Bool_t GetObjectStat ()
 Get status of object stat flag. More...
 
static void SetDtorOnly (void *obj)
 Set destructor only flag. More...
 
static void SetObjectStat (Bool_t stat)
 Turn on/off tracking of objects in the TObjectTable. More...
 

#include <TSpectrum2Fit.h>

Inheritance diagram for TSpectrum2Fit:
[legend]

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
Enumerator
kFitOptimChiCounts 
kFitOptimChiFuncValues 
kFitOptimMaxLikelihood 
kFitAlphaHalving 
kFitAlphaOptimal 
kFitPower2 
kFitPower4 
kFitPower6 
kFitPower8 
kFitPower10 
kFitPower12 
kFitTaylorOrderFirst 
kFitTaylorOrderSecond 
kFitNumRegulCycles 

Definition at line 117 of file TSpectrum2Fit.h.

Constructor & Destructor Documentation

◆ TSpectrum2Fit() [1/2]

TSpectrum2Fit::TSpectrum2Fit ( void  )

Default constructor.

Definition at line 35 of file TSpectrum2Fit.cxx.

◆ TSpectrum2Fit() [2/2]

TSpectrum2Fit::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 non-symmetrical terms and background.

spectrum2fit_constructor_image001.gif

Definition at line 150 of file TSpectrum2Fit.cxx.

◆ ~TSpectrum2Fit()

TSpectrum2Fit::~TSpectrum2Fit ( )
virtual

Destructor.

Definition at line 259 of file TSpectrum2Fit.cxx.

Member Function Documentation

◆ Deramp2()

Double_t TSpectrum2Fit::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 
)
protected

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

Definition at line 541 of file TSpectrum2Fit.cxx.

◆ Derampx()

Double_t TSpectrum2Fit::Derampx ( Double_t  x,
Double_t  x0,
Double_t  sigmax,
Double_t  tx,
Double_t  sx,
Double_t  bx 
)
protected

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

Definition at line 588 of file TSpectrum2Fit.cxx.

◆ Derbx()

Double_t TSpectrum2Fit::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 
)
protected

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

Definition at line 1420 of file TSpectrum2Fit.cxx.

◆ Derby()

Double_t TSpectrum2Fit::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 
)
protected

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

Definition at line 1476 of file TSpectrum2Fit.cxx.

◆ Derderi01()

Double_t TSpectrum2Fit::Derderi01 ( Double_t  x,
Double_t  ax,
Double_t  x0,
Double_t  sigmax 
)
protected

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

Definition at line 865 of file TSpectrum2Fit.cxx.

◆ Derderi02()

Double_t TSpectrum2Fit::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 
)
protected

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

Definition at line 690 of file TSpectrum2Fit.cxx.

◆ Derderj02()

Double_t TSpectrum2Fit::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 
)
protected

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

Definition at line 785 of file TSpectrum2Fit.cxx.

◆ Derdersigmax()

Double_t TSpectrum2Fit::Derdersigmax ( Int_t  numOfFittedPeaks,
Double_t  x,
Double_t  y,
const Double_t parameter,
Double_t  sigmax,
Double_t  sigmay,
Double_t  ro 
)
protected

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

Definition at line 985 of file TSpectrum2Fit.cxx.

◆ Derdersigmay()

Double_t TSpectrum2Fit::Derdersigmay ( Int_t  numOfFittedPeaks,
Double_t  x,
Double_t  y,
const Double_t parameter,
Double_t  sigmax,
Double_t  sigmay,
Double_t  ro 
)
protected

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

Definition at line 1132 of file TSpectrum2Fit.cxx.

◆ Derfc()

Double_t TSpectrum2Fit::Derfc ( Double_t  x)
protected

This function calculates derivative of error function of x.

Definition at line 322 of file TSpectrum2Fit.cxx.

◆ Deri01()

Double_t TSpectrum2Fit::Deri01 ( Double_t  x,
Double_t  ax,
Double_t  x0,
Double_t  sigmax,
Double_t  tx,
Double_t  sx,
Double_t  bx 
)
protected

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

Definition at line 821 of file TSpectrum2Fit.cxx.

◆ Deri02()

Double_t TSpectrum2Fit::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 
)
protected

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

Definition at line 635 of file TSpectrum2Fit.cxx.

◆ Derj02()

Double_t TSpectrum2Fit::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 
)
protected

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

Definition at line 728 of file TSpectrum2Fit.cxx.

◆ Derpa2()

Double_t TSpectrum2Fit::Derpa2 ( Double_t  sx,
Double_t  sy,
Double_t  ro 
)
protected

This function calculates derivative of the volume of a peak according to amplitude.

Function parameters:

  • sx,sy-sigmas of peak
  • ro-correlation coefficient

Definition at line 1549 of file TSpectrum2Fit.cxx.

◆ Derpro()

Double_t TSpectrum2Fit::Derpro ( Double_t  a,
Double_t  sx,
Double_t  sy,
Double_t  ro 
)
protected

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

Definition at line 1618 of file TSpectrum2Fit.cxx.

◆ Derpsigmax()

Double_t TSpectrum2Fit::Derpsigmax ( Double_t  a,
Double_t  sy,
Double_t  ro 
)
protected

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

Definition at line 1572 of file TSpectrum2Fit.cxx.

◆ Derpsigmay()

Double_t TSpectrum2Fit::Derpsigmay ( Double_t  a,
Double_t  sx,
Double_t  ro 
)
protected

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

Definition at line 1595 of file TSpectrum2Fit.cxx.

◆ Derro()

Double_t TSpectrum2Fit::Derro ( Int_t  numOfFittedPeaks,
Double_t  x,
Double_t  y,
const Double_t parameter,
Double_t  sx,
Double_t  sy,
Double_t  r 
)
protected

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

Definition at line 1190 of file TSpectrum2Fit.cxx.

◆ Dersigmax()

Double_t TSpectrum2Fit::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 
)
protected

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

Definition at line 898 of file TSpectrum2Fit.cxx.

◆ Dersigmay()

Double_t TSpectrum2Fit::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 
)
protected

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

Definition at line 1045 of file TSpectrum2Fit.cxx.

◆ Dersx()

Double_t TSpectrum2Fit::Dersx ( Int_t  numOfFittedPeaks,
Double_t  x,
const Double_t parameter,
Double_t  sigmax 
)
protected

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

Definition at line 1363 of file TSpectrum2Fit.cxx.

◆ Dersxy()

Double_t TSpectrum2Fit::Dersxy ( Int_t  numOfFittedPeaks,
Double_t  x,
Double_t  y,
const Double_t parameter,
Double_t  sigmax,
Double_t  sigmay 
)
protected

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

Definition at line 1268 of file TSpectrum2Fit.cxx.

◆ Dersy()

Double_t TSpectrum2Fit::Dersy ( Int_t  numOfFittedPeaks,
Double_t  x,
const Double_t parameter,
Double_t  sigmax 
)
protected

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

Definition at line 1390 of file TSpectrum2Fit.cxx.

◆ Dertx()

Double_t TSpectrum2Fit::Dertx ( Int_t  numOfFittedPeaks,
Double_t  x,
const Double_t parameter,
Double_t  sigmax,
Double_t  bx 
)
protected

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

Definition at line 1298 of file TSpectrum2Fit.cxx.

◆ Dertxy()

Double_t TSpectrum2Fit::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 
)
protected

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

Definition at line 1232 of file TSpectrum2Fit.cxx.

◆ Derty()

Double_t TSpectrum2Fit::Derty ( Int_t  numOfFittedPeaks,
Double_t  x,
const Double_t parameter,
Double_t  sigmax,
Double_t  bx 
)
protected

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

Definition at line 1331 of file TSpectrum2Fit.cxx.

◆ Erfc()

Double_t TSpectrum2Fit::Erfc ( Double_t  x)
protected

This function calculates error function of x.

Definition at line 298 of file TSpectrum2Fit.cxx.

◆ FitAwmi()

void TSpectrum2Fit::FitAwmi ( Double_t **  source)

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 neighbouring 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 non-symmetrical semiempirical peak shape function
  • it contains the two-dimensional symmetrical Gaussian two one-dimensional symmetrical Gaussian ridges as well as non-symmetrical terms and background.

    spectrum2fit_awmi_image001.gif

    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].

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. Matouoek: 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 1 - script FitAwmi2.c:

spectrum2fit_awmi_image002.jpg
Original Fig. 1 two-dimensional spectrum with found peaks (using TSpectrum2 peak searching function). The positions of peaks were used as initial estimates in fitting procedure.
spectrum2fit_awmi_image003.jpg
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-squareafter 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;
Double_t ** source = new float *[nbinsx];
Double_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];
Double_t *PosX = new Double_t[nfound];
Double_t *PosY = new Double_t[nfound];
Double_t *Amp = new Double_t[nfound];
Double_t *AmpXY = new Double_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_t)(PosX[i]+0.5)][(Int_t)(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 FitAwmi2.c:

spectrum2fit_awmi_image004.jpg
Fig. 3 Original two-dimensional gamma-gamma-ray spectrum with found peaks (using TSpectrum2 peak searching function).
spectrum2fit_awmi_image005.jpg
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;
Double_t ** source = new float *[nbinsx];
Double_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];
Double_t *PosX = new Double_t[nfound];
Double_t *PosY = new Double_t[nfound];
Double_t *Amp = new Double_t[nfound];
Double_t *AmpXY = new Double_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_t)(PosX[i]+0.5)][(Int_t)(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");
}

Definition at line 1847 of file TSpectrum2Fit.cxx.

◆ FitStiefel()

void TSpectrum2Fit::FitStiefel ( Double_t **  source)

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

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.

Reference:

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

Example 1 - script FitS.c:

spectrum2fit_stiefel_image001.jpg
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.
spectrum2fit_stiefel_image002.jpg
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;
Double_t ** source = new float *[nbinsx];
Double_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];
Double_t *PosX = new Double_t[nfound];
Double_t *PosY = new Double_t[nfound];
Double_t *Amp = new Double_t[nfound];
Double_t *AmpXY = new Double_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_t)(PosX[i]+0.5)][(Int_t)(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");
}

Definition at line 3942 of file TSpectrum2Fit.cxx.

◆ GetAmplitudeErrors()

void TSpectrum2Fit::GetAmplitudeErrors ( Double_t amplitudeErrors,
Double_t amplitudeErrorsX1,
Double_t amplitudeErrorsY1 
)

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

Definition at line 5754 of file TSpectrum2Fit.cxx.

◆ GetAmplitudes()

void TSpectrum2Fit::GetAmplitudes ( Double_t amplitudes,
Double_t amplitudesX1,
Double_t amplitudesY1 
)

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

Definition at line 5739 of file TSpectrum2Fit.cxx.

◆ GetBackgroundParameters()

void TSpectrum2Fit::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.

  • 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

Definition at line 5827 of file TSpectrum2Fit.cxx.

◆ GetChi()

Double_t TSpectrum2Fit::GetChi ( ) const
inline

Definition at line 176 of file TSpectrum2Fit.h.

◆ GetPositionErrors()

void TSpectrum2Fit::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.

  • 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

Definition at line 5723 of file TSpectrum2Fit.cxx.

◆ GetPositions()

void TSpectrum2Fit::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.

  • 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

Definition at line 5706 of file TSpectrum2Fit.cxx.

◆ GetRo()

void TSpectrum2Fit::GetRo ( Double_t ro,
Double_t roErr 
)

This function gets the ro parameter and its error.

  • ro - gets the fitted value of ro parameter
  • roErr - gets error value of ro parameter

Definition at line 5812 of file TSpectrum2Fit.cxx.

◆ GetSigmaX()

void TSpectrum2Fit::GetSigmaX ( Double_t sigmaX,
Double_t sigmaErrX 
)

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

Definition at line 5790 of file TSpectrum2Fit.cxx.

◆ GetSigmaY()

void TSpectrum2Fit::GetSigmaY ( Double_t sigmaY,
Double_t sigmaErrY 
)

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

Definition at line 5801 of file TSpectrum2Fit.cxx.

◆ GetTailParameters()

void TSpectrum2Fit::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.

  • 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

Definition at line 5856 of file TSpectrum2Fit.cxx.

◆ GetVolumeErrors()

void TSpectrum2Fit::GetVolumeErrors ( Double_t volumeErrors)

This function gets errors of the volumes of fitted 2D peaks.

  • volumeErrors - gets vector of volumes errors of 2D peaks

Definition at line 5778 of file TSpectrum2Fit.cxx.

◆ GetVolumes()

void TSpectrum2Fit::GetVolumes ( Double_t volumes)

This function gets the volumes of fitted 2D peaks.

  • volumes - gets vector of volumes of 2D peaks

Definition at line 5767 of file TSpectrum2Fit.cxx.

◆ Ourpowl()

Double_t TSpectrum2Fit::Ourpowl ( Double_t  a,
Int_t  pw 
)
protected

power function

Definition at line 345 of file TSpectrum2Fit.cxx.

◆ SetBackgroundParameters()

void TSpectrum2Fit::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:

  • a0Init - initial value of a0 parameter (background 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)

Definition at line 5650 of file TSpectrum2Fit.cxx.

◆ SetFitParameters()

void TSpectrum2Fit::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:

  • 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.

Definition at line 5523 of file TSpectrum2Fit.cxx.

◆ SetPeakParameters()

void TSpectrum2Fit::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:

  • 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 - array 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 - array 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 - array of initial values of 2D peaks amplitudes
  • fixAmp - array of logical values which allow to fix appropriate amplitudes of 2D peaks (not fit). However they are present in the estimated functional
  • ampInitX1 - array of initial values of amplitudes of 1D ridges in x direction
  • fixAmpX1 - 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
  • ampInitY1 - array of initial values of amplitudes of 1D ridges in y direction
  • fixAmpY1 - 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

Definition at line 5581 of file TSpectrum2Fit.cxx.

◆ SetTailParameters()

void TSpectrum2Fit::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.

  • 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)

Definition at line 5679 of file TSpectrum2Fit.cxx.

◆ Shape2()

Double_t TSpectrum2Fit::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 
)
protected

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

Definition at line 435 of file TSpectrum2Fit.cxx.

◆ StiefelInversion()

void TSpectrum2Fit::StiefelInversion ( Double_t **  a,
Int_t  size 
)
protected

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

Definition at line 372 of file TSpectrum2Fit.cxx.

◆ Volume()

Double_t TSpectrum2Fit::Volume ( Double_t  a,
Double_t  sx,
Double_t  sy,
Double_t  ro 
)
protected

This function calculates volume of a peak.

Function parameters:

  • a-amplitude of the peak
  • sx,sy-sigmas of peak
  • ro-correlation coefficient

Definition at line 1527 of file TSpectrum2Fit.cxx.

Member Data Documentation

◆ fA0Calc

Double_t TSpectrum2Fit::fA0Calc
protected

calculated value of background a0 parameter

Definition at line 87 of file TSpectrum2Fit.h.

◆ fA0Err

Double_t TSpectrum2Fit::fA0Err
protected

error value of background a0 parameter

Definition at line 88 of file TSpectrum2Fit.h.

◆ fA0Init

Double_t TSpectrum2Fit::fA0Init
protected

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

Definition at line 86 of file TSpectrum2Fit.h.

◆ fAlpha

Double_t TSpectrum2Fit::fAlpha
protected

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

Definition at line 28 of file TSpectrum2Fit.h.

◆ fAlphaOptim

Int_t TSpectrum2Fit::fAlphaOptim
protected

optimization of convergence algorithm, possible values kFitAlphaHalving, kFitAlphaOptimal

Definition at line 25 of file TSpectrum2Fit.h.

◆ fAmpCalc

Double_t* TSpectrum2Fit::fAmpCalc
protected

[fNPeaks] array of calculated values of amplitudes of 2D peaks, output parameters

Definition at line 43 of file TSpectrum2Fit.h.

◆ fAmpCalcX1

Double_t* TSpectrum2Fit::fAmpCalcX1
protected

[fNPeaks] array of calculated values of amplitudes of 1D ridges in x direction, output parameters

Definition at line 46 of file TSpectrum2Fit.h.

◆ fAmpCalcY1

Double_t* TSpectrum2Fit::fAmpCalcY1
protected

[fNPeaks] array of calculated values of amplitudes of 1D ridges in y direction, output parameters

Definition at line 49 of file TSpectrum2Fit.h.

◆ fAmpErr

Double_t* TSpectrum2Fit::fAmpErr
protected

[fNPeaks] array of amplitudes errors of 2D peaks, output parameters

Definition at line 44 of file TSpectrum2Fit.h.

◆ fAmpErrX1

Double_t* TSpectrum2Fit::fAmpErrX1
protected

[fNPeaks] array of amplitudes errors of 1D ridges in x direction, output parameters

Definition at line 47 of file TSpectrum2Fit.h.

◆ fAmpErrY1

Double_t* TSpectrum2Fit::fAmpErrY1
protected

[fNPeaks] array of amplitudes errors of 1D ridges in y direction, output parameters

Definition at line 50 of file TSpectrum2Fit.h.

◆ fAmpInit

Double_t* TSpectrum2Fit::fAmpInit
protected

[fNPeaks] array of initial values of amplitudes of 2D peaks, input parameters

Definition at line 42 of file TSpectrum2Fit.h.

◆ fAmpInitX1

Double_t* TSpectrum2Fit::fAmpInitX1
protected

[fNPeaks] array of initial values of amplitudes of 1D ridges in x direction, input parameters

Definition at line 45 of file TSpectrum2Fit.h.

◆ fAmpInitY1

Double_t* TSpectrum2Fit::fAmpInitY1
protected

[fNPeaks] array of initial values of amplitudes of 1D ridges in y direction, input parameters

Definition at line 48 of file TSpectrum2Fit.h.

◆ fAxCalc

Double_t TSpectrum2Fit::fAxCalc
protected

calculated value of background ax parameter

Definition at line 90 of file TSpectrum2Fit.h.

◆ fAxErr

Double_t TSpectrum2Fit::fAxErr
protected

error value of background ax parameter

Definition at line 91 of file TSpectrum2Fit.h.

◆ fAxInit

Double_t TSpectrum2Fit::fAxInit
protected

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

Definition at line 89 of file TSpectrum2Fit.h.

◆ fAyCalc

Double_t TSpectrum2Fit::fAyCalc
protected

calculated value of background ay parameter

Definition at line 93 of file TSpectrum2Fit.h.

◆ fAyErr

Double_t TSpectrum2Fit::fAyErr
protected

error value of background ay parameter

Definition at line 94 of file TSpectrum2Fit.h.

◆ fAyInit

Double_t TSpectrum2Fit::fAyInit
protected

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

Definition at line 92 of file TSpectrum2Fit.h.

◆ fBxCalc

Double_t TSpectrum2Fit::fBxCalc
protected

calculated value of b parameter for 1D ridges in x direction

Definition at line 81 of file TSpectrum2Fit.h.

◆ fBxErr

Double_t TSpectrum2Fit::fBxErr
protected

error value of b parameter for 1D ridges in x direction

Definition at line 82 of file TSpectrum2Fit.h.

◆ fBxInit

Double_t TSpectrum2Fit::fBxInit
protected

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

Definition at line 80 of file TSpectrum2Fit.h.

◆ fByCalc

Double_t TSpectrum2Fit::fByCalc
protected

calculated value of b parameter for 1D ridges in y direction

Definition at line 84 of file TSpectrum2Fit.h.

◆ fByErr

Double_t TSpectrum2Fit::fByErr
protected

error value of b parameter for 1D ridges in y direction

Definition at line 85 of file TSpectrum2Fit.h.

◆ fByInit

Double_t TSpectrum2Fit::fByInit
protected

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

Definition at line 83 of file TSpectrum2Fit.h.

◆ fChi

Double_t TSpectrum2Fit::fChi
protected

here the fitting functions return resulting chi square

Definition at line 29 of file TSpectrum2Fit.h.

◆ fFitTaylor

Int_t TSpectrum2Fit::fFitTaylor
protected

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

Definition at line 27 of file TSpectrum2Fit.h.

◆ fFixA0

Bool_t TSpectrum2Fit::fFixA0
protected

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

Definition at line 113 of file TSpectrum2Fit.h.

◆ fFixAmp

Bool_t* TSpectrum2Fit::fFixAmp
protected

[fNPeaks] array of logical values which allow to fix appropriate amplitudes of 2D peaks (not fit). However they are present in the estimated functional

Definition at line 99 of file TSpectrum2Fit.h.

◆ fFixAmpX1

Bool_t* TSpectrum2Fit::fFixAmpX1
protected

[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

Definition at line 100 of file TSpectrum2Fit.h.

◆ fFixAmpY1

Bool_t* TSpectrum2Fit::fFixAmpY1
protected

[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

Definition at line 101 of file TSpectrum2Fit.h.

◆ fFixAx

Bool_t TSpectrum2Fit::fFixAx
protected

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

Definition at line 114 of file TSpectrum2Fit.h.

◆ fFixAy

Bool_t TSpectrum2Fit::fFixAy
protected

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

Definition at line 115 of file TSpectrum2Fit.h.

◆ fFixBx

Bool_t TSpectrum2Fit::fFixBx
protected

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

Definition at line 111 of file TSpectrum2Fit.h.

◆ fFixBy

Bool_t TSpectrum2Fit::fFixBy
protected

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

Definition at line 112 of file TSpectrum2Fit.h.

◆ fFixPositionX

Bool_t* TSpectrum2Fit::fFixPositionX
protected

[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

Definition at line 95 of file TSpectrum2Fit.h.

◆ fFixPositionX1

Bool_t* TSpectrum2Fit::fFixPositionX1
protected

[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

Definition at line 97 of file TSpectrum2Fit.h.

◆ fFixPositionY

Bool_t* TSpectrum2Fit::fFixPositionY
protected

[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

Definition at line 96 of file TSpectrum2Fit.h.

◆ fFixPositionY1

Bool_t* TSpectrum2Fit::fFixPositionY1
protected

[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

Definition at line 98 of file TSpectrum2Fit.h.

◆ fFixRo

Bool_t TSpectrum2Fit::fFixRo
protected

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

Definition at line 104 of file TSpectrum2Fit.h.

◆ fFixSigmaX

Bool_t TSpectrum2Fit::fFixSigmaX
protected

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

Definition at line 102 of file TSpectrum2Fit.h.

◆ fFixSigmaY

Bool_t TSpectrum2Fit::fFixSigmaY
protected

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

Definition at line 103 of file TSpectrum2Fit.h.

◆ fFixSx

Bool_t TSpectrum2Fit::fFixSx
protected

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

Definition at line 109 of file TSpectrum2Fit.h.

◆ fFixSxy

Bool_t TSpectrum2Fit::fFixSxy
protected

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

Definition at line 106 of file TSpectrum2Fit.h.

◆ fFixSy

Bool_t TSpectrum2Fit::fFixSy
protected

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

Definition at line 110 of file TSpectrum2Fit.h.

◆ fFixTx

Bool_t TSpectrum2Fit::fFixTx
protected

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

Definition at line 107 of file TSpectrum2Fit.h.

◆ fFixTxy

Bool_t TSpectrum2Fit::fFixTxy
protected

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

Definition at line 105 of file TSpectrum2Fit.h.

◆ fFixTy

Bool_t TSpectrum2Fit::fFixTy
protected

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

Definition at line 108 of file TSpectrum2Fit.h.

◆ fNPeaks

Int_t TSpectrum2Fit::fNPeaks
protected

number of peaks present in fit, input parameter, it should be > 0

Definition at line 18 of file TSpectrum2Fit.h.

◆ fNumberIterations

Int_t TSpectrum2Fit::fNumberIterations
protected

number of iterations in fitting procedure, input parameter, it should be > 0

Definition at line 19 of file TSpectrum2Fit.h.

◆ fPositionCalcX

Double_t* TSpectrum2Fit::fPositionCalcX
protected

[fNPeaks] array of calculated values of x positions of 2D peaks, output parameters

Definition at line 31 of file TSpectrum2Fit.h.

◆ fPositionCalcX1

Double_t* TSpectrum2Fit::fPositionCalcX1
protected

[fNPeaks] array of calculated x positions of 1D ridges, output parameters

Definition at line 37 of file TSpectrum2Fit.h.

◆ fPositionCalcY

Double_t* TSpectrum2Fit::fPositionCalcY
protected

[fNPeaks] array of calculated values of y positions of 2D peaks, output parameters

Definition at line 34 of file TSpectrum2Fit.h.

◆ fPositionCalcY1

Double_t* TSpectrum2Fit::fPositionCalcY1
protected

[fNPeaks] array of calculated y positions of 1D ridges, output parameters

Definition at line 40 of file TSpectrum2Fit.h.

◆ fPositionErrX

Double_t* TSpectrum2Fit::fPositionErrX
protected

[fNPeaks] array of error values of x positions of 2D peaks, output parameters

Definition at line 32 of file TSpectrum2Fit.h.

◆ fPositionErrX1

Double_t* TSpectrum2Fit::fPositionErrX1
protected

[fNPeaks] array of x positions errors of 1D ridges, output parameters

Definition at line 38 of file TSpectrum2Fit.h.

◆ fPositionErrY

Double_t* TSpectrum2Fit::fPositionErrY
protected

[fNPeaks] array of error values of y positions of 2D peaks, output parameters

Definition at line 35 of file TSpectrum2Fit.h.

◆ fPositionErrY1

Double_t* TSpectrum2Fit::fPositionErrY1
protected

[fNPeaks] array of y positions errors of 1D ridges, output parameters

Definition at line 41 of file TSpectrum2Fit.h.

◆ fPositionInitX

Double_t* TSpectrum2Fit::fPositionInitX
protected

[fNPeaks] array of initial values of x positions of 2D peaks, input parameters

Definition at line 30 of file TSpectrum2Fit.h.

◆ fPositionInitX1

Double_t* TSpectrum2Fit::fPositionInitX1
protected

[fNPeaks] array of initial x positions of 1D ridges, input parameters

Definition at line 36 of file TSpectrum2Fit.h.

◆ fPositionInitY

Double_t* TSpectrum2Fit::fPositionInitY
protected

[fNPeaks] array of initial values of y positions of 2D peaks, input parameters

Definition at line 33 of file TSpectrum2Fit.h.

◆ fPositionInitY1

Double_t* TSpectrum2Fit::fPositionInitY1
protected

[fNPeaks] array of initial y positions of 1D ridges, input parameters

Definition at line 39 of file TSpectrum2Fit.h.

◆ fPower

Int_t TSpectrum2Fit::fPower
protected

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

Definition at line 26 of file TSpectrum2Fit.h.

◆ fRoCalc

Double_t TSpectrum2Fit::fRoCalc
protected

calculated value of correlation coefficient

Definition at line 60 of file TSpectrum2Fit.h.

◆ fRoErr

Double_t TSpectrum2Fit::fRoErr
protected

error value of correlation coefficient

Definition at line 61 of file TSpectrum2Fit.h.

◆ fRoInit

Double_t TSpectrum2Fit::fRoInit
protected

initial value of correlation coefficient

Definition at line 59 of file TSpectrum2Fit.h.

◆ fSigmaCalcX

Double_t TSpectrum2Fit::fSigmaCalcX
protected

calculated value of sigma x parameter

Definition at line 54 of file TSpectrum2Fit.h.

◆ fSigmaCalcY

Double_t TSpectrum2Fit::fSigmaCalcY
protected

calculated value of sigma y parameter

Definition at line 57 of file TSpectrum2Fit.h.

◆ fSigmaErrX

Double_t TSpectrum2Fit::fSigmaErrX
protected

error value of sigma x parameter

Definition at line 55 of file TSpectrum2Fit.h.

◆ fSigmaErrY

Double_t TSpectrum2Fit::fSigmaErrY
protected

error value of sigma y parameter

Definition at line 58 of file TSpectrum2Fit.h.

◆ fSigmaInitX

Double_t TSpectrum2Fit::fSigmaInitX
protected

initial value of sigma x parameter

Definition at line 53 of file TSpectrum2Fit.h.

◆ fSigmaInitY

Double_t TSpectrum2Fit::fSigmaInitY
protected

initial value of sigma y parameter

Definition at line 56 of file TSpectrum2Fit.h.

◆ fStatisticType

Int_t TSpectrum2Fit::fStatisticType
protected

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

Definition at line 24 of file TSpectrum2Fit.h.

◆ fSxCalc

Double_t TSpectrum2Fit::fSxCalc
protected

calculated value of s parameter for 1D ridges in x direction

Definition at line 75 of file TSpectrum2Fit.h.

◆ fSxErr

Double_t TSpectrum2Fit::fSxErr
protected

error value of s parameter for 1D ridges in x direction

Definition at line 76 of file TSpectrum2Fit.h.

◆ fSxInit

Double_t TSpectrum2Fit::fSxInit
protected

initial value of s parameter for 1D ridges in x direction (relative amplitude of step), for details see html manual and references

Definition at line 74 of file TSpectrum2Fit.h.

◆ fSxyCalc

Double_t TSpectrum2Fit::fSxyCalc
protected

calculated value of s parameter for 2D peaks

Definition at line 66 of file TSpectrum2Fit.h.

◆ fSxyErr

Double_t TSpectrum2Fit::fSxyErr
protected

error value of s parameter for 2D peaks

Definition at line 67 of file TSpectrum2Fit.h.

◆ fSxyInit

Double_t TSpectrum2Fit::fSxyInit
protected

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

Definition at line 65 of file TSpectrum2Fit.h.

◆ fSyCalc

Double_t TSpectrum2Fit::fSyCalc
protected

calculated value of s parameter for 1D ridges in y direction

Definition at line 78 of file TSpectrum2Fit.h.

◆ fSyErr

Double_t TSpectrum2Fit::fSyErr
protected

error value of s parameter for 1D ridges in y direction

Definition at line 79 of file TSpectrum2Fit.h.

◆ fSyInit

Double_t TSpectrum2Fit::fSyInit
protected

initial value of s parameter for 1D ridges in y direction (relative amplitude of step), for details see html manual and references

Definition at line 77 of file TSpectrum2Fit.h.

◆ fTxCalc

Double_t TSpectrum2Fit::fTxCalc
protected

calculated value of t parameter for 1D ridges in x direction

Definition at line 69 of file TSpectrum2Fit.h.

◆ fTxErr

Double_t TSpectrum2Fit::fTxErr
protected

error value of t parameter for 1D ridges in x direction

Definition at line 70 of file TSpectrum2Fit.h.

◆ fTxInit

Double_t TSpectrum2Fit::fTxInit
protected

initial value of t parameter for 1D ridges in x direction (relative amplitude of tail), for details see html manual and references

Definition at line 68 of file TSpectrum2Fit.h.

◆ fTxyCalc

Double_t TSpectrum2Fit::fTxyCalc
protected

calculated value of t parameter for 2D peaks

Definition at line 63 of file TSpectrum2Fit.h.

◆ fTxyErr

Double_t TSpectrum2Fit::fTxyErr
protected

error value of t parameter for 2D peaks

Definition at line 64 of file TSpectrum2Fit.h.

◆ fTxyInit

Double_t TSpectrum2Fit::fTxyInit
protected

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

Definition at line 62 of file TSpectrum2Fit.h.

◆ fTyCalc

Double_t TSpectrum2Fit::fTyCalc
protected

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

Definition at line 72 of file TSpectrum2Fit.h.

◆ fTyErr

Double_t TSpectrum2Fit::fTyErr
protected

error value of t parameter for 1D ridges in y direction

Definition at line 73 of file TSpectrum2Fit.h.

◆ fTyInit

Double_t TSpectrum2Fit::fTyInit
protected

initial value of t parameter for 1D ridges in y direction (relative amplitude of tail), for details see html manual and references

Definition at line 71 of file TSpectrum2Fit.h.

◆ fVolume

Double_t* TSpectrum2Fit::fVolume
protected

[fNPeaks] array of calculated volumes of 2D peaks, output parameters

Definition at line 51 of file TSpectrum2Fit.h.

◆ fVolumeErr

Double_t* TSpectrum2Fit::fVolumeErr
protected

[fNPeaks] array of volumes errors of 2D peaks, output parameters

Definition at line 52 of file TSpectrum2Fit.h.

◆ fXmax

Int_t TSpectrum2Fit::fXmax
protected

last fitted channel in x direction

Definition at line 21 of file TSpectrum2Fit.h.

◆ fXmin

Int_t TSpectrum2Fit::fXmin
protected

first fitted channel in x direction

Definition at line 20 of file TSpectrum2Fit.h.

◆ fYmax

Int_t TSpectrum2Fit::fYmax
protected

last fitted channel in y direction

Definition at line 23 of file TSpectrum2Fit.h.

◆ fYmin

Int_t TSpectrum2Fit::fYmin
protected

first fitted channel in y direction

Definition at line 22 of file TSpectrum2Fit.h.

Libraries for TSpectrum2Fit:
[legend]

The documentation for this class was generated from the following files: