ROOT 6.12/07 Reference Guide |
Advanced 1-dimensional spectra fitting functions.
Class for fitting 1D 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:
Definition at line 18 of file TSpectrumFit.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 | |
TSpectrumFit (void) | |
Default constructor. More... | |
TSpectrumFit (Int_t numberPeaks) | |
numberPeaks: number of fitted peaks (must be greater than zero) More... | |
virtual | ~TSpectrumFit () |
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... | |
Double_t * | GetAmplitudes () const |
Double_t * | GetAmplitudesErrors () const |
Double_t * | GetAreas () const |
Double_t * | GetAreasErrors () const |
void | GetBackgroundParameters (Double_t &a0, Double_t &a0Err, Double_t &a1, Double_t &a1Err, Double_t &a2, Double_t &a2Err) |
This function gets the background parameters and their errors. More... | |
Double_t | GetChi () const |
Double_t * | GetPositions () const |
Double_t * | GetPositionsErrors () const |
void | GetSigma (Double_t &sigma, Double_t &sigmaErr) |
This function gets the sigma parameter and its error. More... | |
void | GetTailParameters (Double_t &t, Double_t &tErr, Double_t &b, Double_t &bErr, Double_t &s, Double_t &sErr) |
This function gets the tail parameters and their errors. More... | |
void | SetBackgroundParameters (Double_t a0Init, Bool_t fixA0, Double_t a1Init, Bool_t fixA1, Double_t a2Init, Bool_t fixA2) |
This function sets the following fitting parameters of background: More... | |
void | SetFitParameters (Int_t xmin, Int_t xmax, Int_t numberIterations, Double_t alpha, Int_t statisticType, Int_t alphaOptim, Int_t power, Int_t fitTaylor) |
This function sets the following fitting parameters: More... | |
void | SetPeakParameters (Double_t sigma, Bool_t fixSigma, const Double_t *positionInit, const Bool_t *fixPosition, const Double_t *ampInit, const Bool_t *fixAmp) |
This function sets the following fitting parameters of peaks: More... | |
void | SetTailParameters (Double_t tInit, Bool_t fixT, Double_t bInit, Bool_t fixB, Double_t sInit, Bool_t fixS) |
This function sets the following fitting parameters of tails of peaks. 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 TObject * | Clone (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... | |
TNamed & | operator= (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 TObject * | DrawClone (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 TObject * | FindObject (const char *name) const |
Must be redefined in derived classes. More... | |
virtual TObject * | FindObject (const TObject *obj) const |
Must be redefined in derived classes. More... | |
virtual Option_t * | GetDrawOption () 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_t * | GetOption () 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... | |
void * | operator new (size_t sz) |
void * | operator new (size_t sz, void *vp) |
void * | operator new[] (size_t sz) |
void * | operator new[] (size_t sz, void *vp) |
TObject & | operator= (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 | Area (Double_t a, Double_t sigma, Double_t t, Double_t b) |
This function calculates area of a peak Function parameters: More... | |
Double_t | Dera1 (Double_t i) |
Derivative of background according to a1. More... | |
Double_t | Dera2 (Double_t i) |
Derivative of background according to a2. More... | |
Double_t | Deramp (Double_t i, Double_t i0, Double_t sigma, Double_t t, Double_t s, Double_t b) |
This function calculates derivative of peak shape function (see manual) according to amplitude of peak. More... | |
Double_t | Derb (Int_t num_of_fitted_peaks, Double_t i, const Double_t *parameter, Double_t sigma, Double_t t, Double_t b) |
This function calculates derivative of peaks shape function (see manual) according to slope b. More... | |
Double_t | Derderi0 (Double_t i, Double_t amp, Double_t i0, Double_t sigma) |
This function calculates second derivative of peak shape function (see manual) according to peak position. More... | |
Double_t | Derdersigma (Int_t num_of_fitted_peaks, Double_t i, const Double_t *parameter, Double_t sigma) |
This function calculates second derivative of peaks shape function (see manual) according to sigma of peaks. More... | |
Double_t | Derfc (Double_t x) |
This function calculates derivative of error function of x. More... | |
Double_t | Deri0 (Double_t i, Double_t amp, Double_t i0, Double_t sigma, Double_t t, Double_t s, Double_t b) |
This function calculates derivative of peak shape function (see manual) according to peak position. More... | |
Double_t | Derpa (Double_t sigma, Double_t t, Double_t b) |
This function calculates derivative of the area of peak according to its amplitude. More... | |
Double_t | Derpb (Double_t a, Double_t sigma, Double_t t, Double_t b) |
This function calculates derivative of the area of peak according to b parameter. More... | |
Double_t | Derpsigma (Double_t a, Double_t t, Double_t b) |
This function calculates derivative of the area of peak according to sigma of peaks. More... | |
Double_t | Derpt (Double_t a, Double_t sigma, Double_t b) |
This function calculates derivative of the area of peak according to t parameter. More... | |
Double_t | Ders (Int_t num_of_fitted_peaks, Double_t i, const Double_t *parameter, Double_t sigma) |
This function calculates derivative of peaks shape function (see manual) according to relative amplitude s. More... | |
Double_t | Dersigma (Int_t num_of_fitted_peaks, Double_t i, const Double_t *parameter, Double_t sigma, Double_t t, Double_t s, Double_t b) |
This function calculates derivative of peaks shape function (see manual) according to sigma of peaks. More... | |
Double_t | Dert (Int_t num_of_fitted_peaks, Double_t i, const Double_t *parameter, Double_t sigma, Double_t b) |
This function calculates derivative of peaks shape function (see manual) according to relative amplitude t. More... | |
Double_t | Erfc (Double_t x) |
Double_t | Ourpowl (Double_t a, Int_t pw) |
Power function. More... | |
Double_t | Shape (Int_t num_of_fitted_peaks, Double_t i, const Double_t *parameter, Double_t sigma, Double_t t, Double_t s, Double_t b, Double_t a0, Double_t a1, Double_t a2) |
This function calculates peaks shape function (see manual) Function parameters: More... | |
void | StiefelInversion (Double_t **a, Int_t rozmer) |
This function calculates solution of the system of linear equations. 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+a1*x+a2*x*x) More... | |
Double_t | fA1Calc |
calculated value of background a1 parameter More... | |
Double_t | fA1Err |
error value of background a1 parameter More... | |
Double_t | fA1Init |
initial value of background a1 parameter(backgroud is estimated as a0+a1*x+a2*x*x) More... | |
Double_t | fA2Calc |
calculated value of background a2 parameter More... | |
Double_t | fA2Err |
error value of background a2 parameter More... | |
Double_t | fA2Init |
initial value of background a2 parameter(backgroud is estimated as a0+a1*x+a2*x*x) 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_t * | fAmpCalc |
[fNPeaks] array of calculated values of fitted amplitudes, output parameters More... | |
Double_t * | fAmpErr |
[fNPeaks] array of amplitude errors More... | |
Double_t * | fAmpInit |
[fNPeaks] array of initial values of peaks amplitudes, input parameters More... | |
Double_t * | fArea |
[fNPeaks] array of calculated areas of peaks More... | |
Double_t * | fAreaErr |
[fNPeaks] array of errors of peak areas More... | |
Double_t | fBCalc |
calculated value of b parameter More... | |
Double_t | fBErr |
error value of b parameter More... | |
Double_t | fBInit |
initial value of b parameter (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_t | fFixA1 |
logical value of a1 parameter, which allows to fix the parameter (not to fit). More... | |
Bool_t | fFixA2 |
logical value of a2 parameter, which allows to fix the parameter (not to fit). More... | |
Bool_t * | fFixAmp |
[fNPeaks] array of logical values which allow to fix appropriate amplitudes (not fit). However they are present in the estimated functional More... | |
Bool_t | fFixB |
logical value of b parameter, which allows to fix the parameter (not to fit). More... | |
Bool_t * | fFixPosition |
[fNPeaks] array of logical values which allow to fix appropriate positions (not fit). However they are present in the estimated functional More... | |
Bool_t | fFixS |
logical value of s parameter, which allows to fix the parameter (not to fit). More... | |
Bool_t | fFixSigma |
logical value of sigma parameter, which allows to fix the parameter (not to fit). More... | |
Bool_t | fFixT |
logical value of t parameter, 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_t * | fPositionCalc |
[fNPeaks] array of calculated values of fitted positions, output parameters More... | |
Double_t * | fPositionErr |
[fNPeaks] array of position errors More... | |
Double_t * | fPositionInit |
[fNPeaks] array of initial values of peaks positions, 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 | fSCalc |
calculated value of s parameter More... | |
Double_t | fSErr |
error value of s parameter More... | |
Double_t | fSigmaCalc |
calculated value of sigma parameter More... | |
Double_t | fSigmaErr |
error value of sigma parameter More... | |
Double_t | fSigmaInit |
initial value of sigma parameter More... | |
Double_t | fSInit |
initial value of s parameter (relative amplitude of step), for details see html manual and references 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 | fTCalc |
calculated value of t parameter More... | |
Double_t | fTErr |
error value of t parameter More... | |
Double_t | fTInit |
initial value of t parameter (relative amplitude of tail), for details see html manual and references More... | |
Int_t | fXmax |
last fitted channel More... | |
Int_t | fXmin |
first fitted channel 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 <TSpectrumFit.h>
anonymous enum |
Definition at line 70 of file TSpectrumFit.h.
TSpectrumFit::TSpectrumFit | ( | void | ) |
Default constructor.
Definition at line 35 of file TSpectrumFit.cxx.
TSpectrumFit::TSpectrumFit | ( | 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 is
where a represents vector of fitted parameters (positions p(j), amplitudes A(j), sigma, relative amplitudes T, S and slope B).
Definition at line 101 of file TSpectrumFit.cxx.
|
virtual |
Destructor.
Definition at line 160 of file TSpectrumFit.cxx.
This function calculates area of a peak Function parameters:
Definition at line 569 of file TSpectrumFit.cxx.
Derivative of background according to a1.
Definition at line 492 of file TSpectrumFit.cxx.
Derivative of background according to a2.
Definition at line 500 of file TSpectrumFit.cxx.
|
protected |
This function calculates derivative of peak shape function (see manual) according to amplitude of peak.
Function parameters:
Definition at line 229 of file TSpectrumFit.cxx.
|
protected |
This function calculates derivative of peaks shape function (see manual) according to slope b.
Function parameters:
Definition at line 463 of file TSpectrumFit.cxx.
This function calculates second derivative of peak shape function (see manual) according to peak position.
Function parameters:
Definition at line 303 of file TSpectrumFit.cxx.
|
protected |
This function calculates second derivative of peaks shape function (see manual) according to sigma of peaks.
Function parameters:
Definition at line 376 of file TSpectrumFit.cxx.
This function calculates derivative of error function of x.
Definition at line 200 of file TSpectrumFit.cxx.
|
protected |
This function calculates derivative of peak shape function (see manual) according to peak position.
Function parameters:
Definition at line 266 of file TSpectrumFit.cxx.
This function calculates derivative of the area of peak according to its amplitude.
Function parameters:
Definition at line 592 of file TSpectrumFit.cxx.
This function calculates derivative of the area of peak according to b parameter.
Function parameters:
Definition at line 658 of file TSpectrumFit.cxx.
This function calculates derivative of the area of peak according to sigma of peaks.
Function parameters:
Definition at line 614 of file TSpectrumFit.cxx.
This function calculates derivative of the area of peak according to t parameter.
Function parameters:
Definition at line 636 of file TSpectrumFit.cxx.
|
protected |
This function calculates derivative of peaks shape function (see manual) according to relative amplitude s.
Function parameters:
Definition at line 437 of file TSpectrumFit.cxx.
|
protected |
This function calculates derivative of peaks shape function (see manual) according to sigma of peaks.
Function parameters:
Definition at line 331 of file TSpectrumFit.cxx.
|
protected |
This function calculates derivative of peaks shape function (see manual) according to relative amplitude t.
Function parameters:
Definition at line 409 of file TSpectrumFit.cxx.
Definition at line 177 of file TSpectrumFit.cxx.
This function fits the source spectrum.
The calling program should fill in input parameters of the TSpectrumFit class The fitted parameters are written into TSpectrumFit class output parameters and fitted data are written into source spectrum.
Function parameters:
Goal: to estimate simultaneously peak shape parameters in spectra with large number of peaks
where T and S are relative amplitudes and B is slope.
[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. Morhac, J. Kliman, M. Jandel, L. 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 to illustrate fitting function using AWMI algorithm. To execute this example, do:
root > .x FitAwmi.C
Definition at line 820 of file TSpectrumFit.cxx.
This function fits the source spectrum.
The calling program should fill in input parameters The fitted parameters are written into output parameters and fitted data are written into source spectrum.
Function parameters:
Example to illustrate fitting function using Stiefel-Hestens method. To execute this example, do:
root > .x FitStiefel.C
Definition at line 1857 of file TSpectrumFit.cxx.
|
inline |
Definition at line 116 of file TSpectrumFit.h.
|
inline |
Definition at line 117 of file TSpectrumFit.h.
|
inline |
Definition at line 118 of file TSpectrumFit.h.
|
inline |
Definition at line 119 of file TSpectrumFit.h.
void TSpectrumFit::GetBackgroundParameters | ( | Double_t & | a0, |
Double_t & | a0Err, | ||
Double_t & | a1, | ||
Double_t & | a1Err, | ||
Double_t & | a2, | ||
Double_t & | a2Err | ||
) |
This function gets the background parameters and their errors.
Definition at line 2740 of file TSpectrumFit.cxx.
|
inline |
Definition at line 121 of file TSpectrumFit.h.
|
inline |
Definition at line 122 of file TSpectrumFit.h.
|
inline |
Definition at line 123 of file TSpectrumFit.h.
This function gets the sigma parameter and its error.
Definition at line 2725 of file TSpectrumFit.cxx.
void TSpectrumFit::GetTailParameters | ( | Double_t & | t, |
Double_t & | tErr, | ||
Double_t & | b, | ||
Double_t & | bErr, | ||
Double_t & | s, | ||
Double_t & | sErr | ||
) |
This function gets the tail parameters and their errors.
Definition at line 2760 of file TSpectrumFit.cxx.
Power function.
Definition at line 674 of file TSpectrumFit.cxx.
void TSpectrumFit::SetBackgroundParameters | ( | Double_t | a0Init, |
Bool_t | fixA0, | ||
Double_t | a1Init, | ||
Bool_t | fixA1, | ||
Double_t | a2Init, | ||
Bool_t | fixA2 | ||
) |
This function sets the following fitting parameters of background:
Definition at line 2691 of file TSpectrumFit.cxx.
void TSpectrumFit::SetFitParameters | ( | Int_t | xmin, |
Int_t | xmax, | ||
Int_t | numberIterations, | ||
Double_t | alpha, | ||
Int_t | statisticType, | ||
Int_t | alphaOptim, | ||
Int_t | power, | ||
Int_t | fitTaylor | ||
) |
This function sets the following fitting parameters:
Definition at line 2608 of file TSpectrumFit.cxx.
void TSpectrumFit::SetPeakParameters | ( | Double_t | sigma, |
Bool_t | fixSigma, | ||
const Double_t * | positionInit, | ||
const Bool_t * | fixPosition, | ||
const Double_t * | ampInit, | ||
const Bool_t * | fixAmp | ||
) |
This function sets the following fitting parameters of peaks:
Definition at line 2656 of file TSpectrumFit.cxx.
void TSpectrumFit::SetTailParameters | ( | Double_t | tInit, |
Bool_t | fixT, | ||
Double_t | bInit, | ||
Bool_t | fixB, | ||
Double_t | sInit, | ||
Bool_t | fixS | ||
) |
This function sets the following fitting parameters of tails of peaks.
Definition at line 2710 of file TSpectrumFit.cxx.
|
protected |
This function calculates peaks shape function (see manual) Function parameters:
Definition at line 516 of file TSpectrumFit.cxx.
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:
Definition at line 1721 of file TSpectrumFit.cxx.
|
protected |
calculated value of background a0 parameter
Definition at line 51 of file TSpectrumFit.h.
|
protected |
error value of background a0 parameter
Definition at line 52 of file TSpectrumFit.h.
|
protected |
initial value of background a0 parameter(backgroud is estimated as a0+a1*x+a2*x*x)
Definition at line 50 of file TSpectrumFit.h.
|
protected |
calculated value of background a1 parameter
Definition at line 54 of file TSpectrumFit.h.
|
protected |
error value of background a1 parameter
Definition at line 55 of file TSpectrumFit.h.
|
protected |
initial value of background a1 parameter(backgroud is estimated as a0+a1*x+a2*x*x)
Definition at line 53 of file TSpectrumFit.h.
|
protected |
calculated value of background a2 parameter
Definition at line 57 of file TSpectrumFit.h.
|
protected |
error value of background a2 parameter
Definition at line 58 of file TSpectrumFit.h.
|
protected |
initial value of background a2 parameter(backgroud is estimated as a0+a1*x+a2*x*x)
Definition at line 56 of file TSpectrumFit.h.
|
protected |
convergence coefficient, input parameter, it should be positive number and <=1, for details see references
Definition at line 28 of file TSpectrumFit.h.
|
protected |
optimization of convergence algorithm, possible values kFitAlphaHalving, kFitAlphaOptimal
Definition at line 25 of file TSpectrumFit.h.
|
protected |
[fNPeaks] array of calculated values of fitted amplitudes, output parameters
Definition at line 34 of file TSpectrumFit.h.
|
protected |
[fNPeaks] array of amplitude errors
Definition at line 35 of file TSpectrumFit.h.
|
protected |
[fNPeaks] array of initial values of peaks amplitudes, input parameters
Definition at line 33 of file TSpectrumFit.h.
|
protected |
[fNPeaks] array of calculated areas of peaks
Definition at line 36 of file TSpectrumFit.h.
|
protected |
[fNPeaks] array of errors of peak areas
Definition at line 37 of file TSpectrumFit.h.
|
protected |
calculated value of b parameter
Definition at line 45 of file TSpectrumFit.h.
|
protected |
error value of b parameter
Definition at line 46 of file TSpectrumFit.h.
|
protected |
initial value of b parameter (slope), for details see html manual and references
Definition at line 44 of file TSpectrumFit.h.
|
protected |
here the fitting functions return resulting chi square
Definition at line 29 of file TSpectrumFit.h.
|
protected |
order of Taylor expansion, possible values kFitTaylorOrderFirst, kFitTaylorOrderSecond. It applies only for Awmi fitting function.
Definition at line 27 of file TSpectrumFit.h.
|
protected |
logical value of a0 parameter, which allows to fix the parameter (not to fit).
Definition at line 65 of file TSpectrumFit.h.
|
protected |
logical value of a1 parameter, which allows to fix the parameter (not to fit).
Definition at line 66 of file TSpectrumFit.h.
|
protected |
logical value of a2 parameter, which allows to fix the parameter (not to fit).
Definition at line 67 of file TSpectrumFit.h.
|
protected |
[fNPeaks] array of logical values which allow to fix appropriate amplitudes (not fit). However they are present in the estimated functional
Definition at line 60 of file TSpectrumFit.h.
|
protected |
logical value of b parameter, which allows to fix the parameter (not to fit).
Definition at line 63 of file TSpectrumFit.h.
|
protected |
[fNPeaks] array of logical values which allow to fix appropriate positions (not fit). However they are present in the estimated functional
Definition at line 59 of file TSpectrumFit.h.
|
protected |
logical value of s parameter, which allows to fix the parameter (not to fit).
Definition at line 64 of file TSpectrumFit.h.
|
protected |
logical value of sigma parameter, which allows to fix the parameter (not to fit).
Definition at line 61 of file TSpectrumFit.h.
|
protected |
logical value of t parameter, which allows to fix the parameter (not to fit).
Definition at line 62 of file TSpectrumFit.h.
|
protected |
number of peaks present in fit, input parameter, it should be > 0
Definition at line 20 of file TSpectrumFit.h.
|
protected |
number of iterations in fitting procedure, input parameter, it should be > 0
Definition at line 21 of file TSpectrumFit.h.
|
protected |
[fNPeaks] array of calculated values of fitted positions, output parameters
Definition at line 31 of file TSpectrumFit.h.
|
protected |
[fNPeaks] array of position errors
Definition at line 32 of file TSpectrumFit.h.
|
protected |
[fNPeaks] array of initial values of peaks positions, input parameters
Definition at line 30 of file TSpectrumFit.h.
|
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 TSpectrumFit.h.
|
protected |
calculated value of s parameter
Definition at line 48 of file TSpectrumFit.h.
|
protected |
error value of s parameter
Definition at line 49 of file TSpectrumFit.h.
|
protected |
calculated value of sigma parameter
Definition at line 39 of file TSpectrumFit.h.
|
protected |
error value of sigma parameter
Definition at line 40 of file TSpectrumFit.h.
|
protected |
initial value of sigma parameter
Definition at line 38 of file TSpectrumFit.h.
|
protected |
initial value of s parameter (relative amplitude of step), for details see html manual and references
Definition at line 47 of file TSpectrumFit.h.
|
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 TSpectrumFit.h.
|
protected |
calculated value of t parameter
Definition at line 42 of file TSpectrumFit.h.
|
protected |
error value of t parameter
Definition at line 43 of file TSpectrumFit.h.
|
protected |
initial value of t parameter (relative amplitude of tail), for details see html manual and references
Definition at line 41 of file TSpectrumFit.h.
|
protected |
last fitted channel
Definition at line 23 of file TSpectrumFit.h.
|
protected |
first fitted channel
Definition at line 22 of file TSpectrumFit.h.