FairRoot/PandaRoot
Public Member Functions | List of all members
PndHypGePeakFitFunction Class Reference

#include <PndHypGePeakFitFunction.h>

Public Member Functions

 PndHypGePeakFitFunction ()
 
Double_t PeakFunc (Double_t *x, Double_t *par)
 
Double_t PeakFuncGaussian (Double_t *x, Double_t *par)
 
Double_t PeakFuncFreeSkewedPosition (Double_t *x, Double_t *par)
 
Double_t PeakFuncFreeSkewedPositionWithoutFirstGausian (Double_t *x, Double_t *par)
 
Double_t PeakFuncDoubleGausian (Double_t *x, Double_t *par)
 
Double_t NewFunction (Double_t *x, Double_t *par)
 
Double_t GausOnly (Double_t *x, Double_t *par)
 
Double_t SmoothedStepOnly (Double_t *x, Double_t *par)
 
Double_t LinearOnly (Double_t *x, Double_t *par)
 
Double_t SkewedOnly (Double_t *x, Double_t *par)
 
Double_t FreeSkewedOnly (Double_t *x, Double_t *par)
 
Double_t SecondGausOnly (Double_t *x, Double_t *par)
 
Double_t NewSkewedOnly (Double_t *x, Double_t *par)
 

Detailed Description

Definition at line 30 of file PndHypGePeakFitFunction.h.

Constructor & Destructor Documentation

PndHypGePeakFitFunction::PndHypGePeakFitFunction ( )

Definition at line 28 of file PndHypGePeakFitFunction.cxx.

29 {
30 
31 }

Member Function Documentation

Double_t PndHypGePeakFitFunction::FreeSkewedOnly ( Double_t x,
Double_t par 
)

Definition at line 160 of file PndHypGePeakFitFunction.cxx.

References Double_t, and CAMath::Sqrt().

Referenced by PndHypGeSpectrumAnalyser::FitPeaks(), and PeakFuncFreeSkewedPositionWithoutFirstGausian().

161 {
162  Double_t SkewedGausian = par[6]
163  * TMath::Exp((x[0]-par[8])/par[7])
164  * TMath::Erfc( ((x[0] - par[8]) /(TMath::Sqrt(2)*par[9]) ) + par[9]/TMath::Sqrt(2)/par[7]); //par[6]: amplitude skewed gaus, par[7]: beta, par[8]: x0_s, par[9]: sigma_s
165  return SkewedGausian;
166 }
static T Sqrt(const T &x)
Definition: PndCAMath.h:37
Double_t par[3]
Double_t
Double_t x
Double_t PndHypGePeakFitFunction::GausOnly ( Double_t x,
Double_t par 
)

Definition at line 133 of file PndHypGePeakFitFunction.cxx.

References Double_t.

Referenced by PndHypGeSpectrumAnalyser::FitPeaks(), and PeakFuncGaussian().

134 {
135 
136  Double_t Gausian = par[0] * TMath::Exp(- 0.5*pow((x[0] - par[1])/par[2],2)); // par[0]: amplitude, par[1]: x[0]0; par[2]: sigma
137 
138  return Gausian ;
139 }
Double_t par[3]
Double_t
Double_t x
Double_t PndHypGePeakFitFunction::LinearOnly ( Double_t x,
Double_t par 
)

Definition at line 147 of file PndHypGePeakFitFunction.cxx.

References Double_t.

Referenced by PndHypGeSpectrumAnalyser::FitPeaks(), PeakFuncFreeSkewedPositionWithoutFirstGausian(), and PeakFuncGaussian().

148 {
149  //Double_t Constant = par[4]*x[0]+par[5]; // par[4]: gradient of lin bg, par[5]: const bg
150  Double_t Constant = par[5]; // par[4]: gradient of lin bg, par[5]: const bg
151  return Constant;
152 }
Double_t par[3]
Double_t
Double_t PndHypGePeakFitFunction::NewFunction ( Double_t x,
Double_t par 
)

Definition at line 182 of file PndHypGePeakFitFunction.cxx.

References Double_t, and CAMath::Sqrt().

Referenced by PndHypGeSpectrumAnalyser::FitPeaks().

183 {
184  Double_t value = par[3]*par[0]/2*TMath::Exp(-par[0]/2*(2*par[1]+par[0]*par[2]*par[2]-2*x[0])) * TMath::Erfc(-(par[1]+par[0]*par[2]*par[2]-x[0])/TMath::Sqrt(2)/par[2])+ par[4]*TMath::Gaus(x[0],par[5],par[6]) +par[7];
185  return value;
186 }
static T Sqrt(const T &x)
Definition: PndCAMath.h:37
Double_t par[3]
Double_t
Double_t x
Double_t PndHypGePeakFitFunction::NewSkewedOnly ( Double_t x,
Double_t par 
)

Definition at line 174 of file PndHypGePeakFitFunction.cxx.

References Double_t, and sqrt().

Referenced by PeakFuncFreeSkewedPosition().

175 {
176  Double_t NewSkewedGausian = par[6]
177  * TMath::Exp(-1/2*pow((x[0]-par[8])/par[9],2))
178  * (1 + TMath::Erf(1/sqrt(2)*par[7]*(x[0]-par[8])/par[9]));
179  return NewSkewedGausian;
180 }
friend F32vec4 sqrt(const F32vec4 &a)
Definition: P4_F32vec4.h:29
Double_t par[3]
Double_t
Double_t x
Double_t PndHypGePeakFitFunction::PeakFunc ( Double_t x,
Double_t par 
)

Definition at line 59 of file PndHypGePeakFitFunction.cxx.

References PeakFuncGaussian(), and SkewedOnly().

Referenced by PndHypGeSpectrumAnalyser::FitPeaks().

60 {
61  // Peak shape is taken from gf2 http://www.phy.anl.gov/gammasphere/doc/gf2.hlp
62  // composed by the sum of:
63  // (1) a Gaussian,
64  // (2) a skewed Gaussian,
65  // (3) a smoothed step function to increase the background on the low-energy
66  // side of the peak. Components (2) and/or (3) can easily be set to zero if not
67  // required, and
68  // (4) a constant // not in gf2
69  //
70  // (2) y = constant
71  // * EXP( (x-c)/beta )
72  // * ERFC( (x-c)/(SQRT(2)*sigma) + sigma/(SQRT(2)*beta) )
73  // (3) y = constant * ERFC( (x-c)/(SQRT(2)*sigma) )
74 
75  // par[0]: amplitude gaus, par[1]: x0; par[2]: sigma, par[3]: amplitude step, par[4]: gradient of lin bg
76  // par[5]: amplitude skewed gaus, par[6]: beta
77 
78  return PeakFuncGaussian (x,par) + SkewedOnly(x,par);
79 }
Double_t par[3]
Double_t SkewedOnly(Double_t *x, Double_t *par)
Double_t x
Double_t PeakFuncGaussian(Double_t *x, Double_t *par)
Double_t PndHypGePeakFitFunction::PeakFuncDoubleGausian ( Double_t x,
Double_t par 
)

Definition at line 113 of file PndHypGePeakFitFunction.cxx.

References PeakFuncGaussian(), and SecondGausOnly().

Referenced by PndHypGeSpectrumAnalyser::FitPeaks().

114 {
115  // Peak shape is taken from gf2 http://www.phy.anl.gov/gammasphere/doc/gf2.hlp
116  // composed by the sum of:
117  // (1) a Gaussian,
118  // (2) a second Gaussian,
119  // (3) a smoothed step function to increase the background on the low-energy
120  // side of the peak. Components (2) and/or (3) can easily be set to zero if not
121  // required, and
122  // (4) a constant // not in gf2
123  //
124  // (3) y = constant * ERFC( (x-c)/(SQRT(2)*sigma) )
125 
126  // par[0]: amplitude gaus, par[1]: x0; par[2]: sigma, par[3]: amplitude step, par[4]: gradient of lin bg
127  // par[5]: amplitude second gaus, par[6]: x0_2, par[7]: sigma_2
128 
130 }
Double_t par[3]
Double_t SecondGausOnly(Double_t *x, Double_t *par)
Double_t x
Double_t PeakFuncGaussian(Double_t *x, Double_t *par)
Double_t PndHypGePeakFitFunction::PeakFuncFreeSkewedPosition ( Double_t x,
Double_t par 
)

Definition at line 81 of file PndHypGePeakFitFunction.cxx.

References NewSkewedOnly(), and PeakFuncGaussian().

82 {
83  // Peak shape is taken from gf2 http://www.phy.anl.gov/gammasphere/doc/gf2.hlp
84  // composed by the sum of:
85  // (1) a Gaussian,
86  // (2) a skewed Gaussian,
87  // (3) a smoothed step function to increase the background on the low-energy
88  // side of the peak. Components (2) and/or (3) can easily be set to zero if not
89  // required, and
90  // (4) a constant // not in gf2
91  //
92  // (2) y = constant
93  // * EXP( (x-c)/beta )
94  // * ERFC( (x-c)/(SQRT(2)*sigma) + sigma/(SQRT(2)*beta) )
95  // (3) y = constant * ERFC( (x-c)/(SQRT(2)*sigma) )
96 
97  // par[0]: amplitude gaus, par[1]: x0; par[2]: sigma, par[3]: amplitude step, par[4]: gradient of lin bg
98  // par[5]: amplitude skewed gaus, par[6]: beta, par[7]: x0_s, par[8]: sigma_s
99 
100  //Double_t SkewedGausian = par[5]
101  // * TMath::Exp((x[0]-par[7])/par[6])
102  // * TMath::Erfc( ((x[0] - par[7]) /(TMath::Sqrt(2)*par[8]) ) + par[8]/TMath::Sqrt(2)/par[6]); //par[5]: amplitude skewed gaus, par[6]: beta, par[7]: x0_s, par[8]: sigma_s
103  //return PeakFuncGaussian (x,par) + FreeSkewedOnly(x,par);//SkewedGausian;
104  return PeakFuncGaussian (x,par) + NewSkewedOnly(x,par);
105 
106 }
Double_t par[3]
Double_t x
Double_t PeakFuncGaussian(Double_t *x, Double_t *par)
Double_t NewSkewedOnly(Double_t *x, Double_t *par)
Double_t PndHypGePeakFitFunction::PeakFuncFreeSkewedPositionWithoutFirstGausian ( Double_t x,
Double_t par 
)

Definition at line 108 of file PndHypGePeakFitFunction.cxx.

References FreeSkewedOnly(), LinearOnly(), and SmoothedStepOnly().

Referenced by PndHypGeSpectrumAnalyser::FitPeaks().

109 {
111 }
Double_t LinearOnly(Double_t *x, Double_t *par)
Double_t par[3]
Double_t FreeSkewedOnly(Double_t *x, Double_t *par)
Double_t SmoothedStepOnly(Double_t *x, Double_t *par)
Double_t x
Double_t PndHypGePeakFitFunction::PeakFuncGaussian ( Double_t x,
Double_t par 
)

Definition at line 33 of file PndHypGePeakFitFunction.cxx.

References GausOnly(), LinearOnly(), and SmoothedStepOnly().

Referenced by PndHypGeSpectrumAnalyser::FitPeaks(), PeakFunc(), PeakFuncDoubleGausian(), and PeakFuncFreeSkewedPosition().

34 {
35  // Peak shape is taken from gf2 http://www.phy.anl.gov/gammasphere/doc/gf2.hlp
36  // composed by the sum of:
37  // (1) a Gaussian,
38  // (2) a smoothed step function to increase the background on the low-energy
39  // side of the peak. Components (2) and/or (3) can easily be set to zero if not
40  // required, and
41  // (3) a constant // not in gf2
42  //
43 
44  // (2) y = constant * ERFC( (x-c)/(SQRT(2)*sigma) )
45 
46  // par[0]: amplitude gaus, par[1]: x0; par[2]: sigma, par[3]: amplitude step, par[4]: gradient of lin bg
47 
48  // with sigma in amplitude
49  //Double_t Gausian = par[0]*TMath::Sqrt(2 * TMath::Pi())/par[2] * TMath::Exp(- 0.5*pow((x[0] - par[1])/par[2],2)); // par[0]: amplitude, par[1]: x[0]0; par[2]: sigma
50  // or without sigma in amplitude
51 
52  //Double_t Gausian = par[0]/*TMath::Sqrt(2 * TMath::Pi())/par[2]*/ * TMath::Exp(- 0.5*pow((x[0] - par[1])/par[2],2)); // par[0]: amplitude, par[1]: x[0]0; par[2]: sigma
53  //Double_t SmoothedStep = par[3] * TMath::Erfc( (x[0] - par[1]) /(TMath::Sqrt(2)*par[2])); // par[3]: amplitude
54  //Double_t Constant = par[4]*x[0]; // par[4]: gradient of lin b
55  //return Gausian + SmoothedStep + Constant;
57 }
Double_t LinearOnly(Double_t *x, Double_t *par)
Double_t par[3]
Double_t SmoothedStepOnly(Double_t *x, Double_t *par)
Double_t x
Double_t GausOnly(Double_t *x, Double_t *par)
Double_t PndHypGePeakFitFunction::SecondGausOnly ( Double_t x,
Double_t par 
)

Definition at line 168 of file PndHypGePeakFitFunction.cxx.

References Double_t.

Referenced by PndHypGeSpectrumAnalyser::FitPeaks(), and PeakFuncDoubleGausian().

169 {
170  Double_t SecondGausian = par[6] * TMath::Exp(- 0.5*pow((x[0] - par[7])/par[8],2)); // par[6]: amplitude, par[7]: x0_2; par[8]: sigma_2
171  return SecondGausian;
172 }
Double_t par[3]
Double_t
Double_t x
Double_t PndHypGePeakFitFunction::SkewedOnly ( Double_t x,
Double_t par 
)

Definition at line 153 of file PndHypGePeakFitFunction.cxx.

References Double_t, and CAMath::Sqrt().

Referenced by PndHypGeSpectrumAnalyser::FitPeaks(), and PeakFunc().

154 {
155  Double_t SkewedGausian = par[6]
156  * TMath::Exp((x[0]-par[1])/par[7])
157  * TMath::Erfc( ((x[0] - par[1]) /(TMath::Sqrt(2)*par[2]) ) + par[2]/TMath::Sqrt(2)/par[7]); // par[6]: amplitude skewed gaus, par[7]: beta
158  return SkewedGausian;
159 }
static T Sqrt(const T &x)
Definition: PndCAMath.h:37
Double_t par[3]
Double_t
Double_t x
Double_t PndHypGePeakFitFunction::SmoothedStepOnly ( Double_t x,
Double_t par 
)

Definition at line 141 of file PndHypGePeakFitFunction.cxx.

References Double_t, and CAMath::Sqrt().

Referenced by PndHypGeSpectrumAnalyser::FitPeaks(), PeakFuncFreeSkewedPositionWithoutFirstGausian(), and PeakFuncGaussian().

142 {
143  Double_t SmoothedStep = par[3] * TMath::Erfc( (x[0] - par[1]) /(TMath::Sqrt(2)*par[2])); // par[3]: amplitude
144  return SmoothedStep;
145 }
static T Sqrt(const T &x)
Definition: PndCAMath.h:37
Double_t par[3]
Double_t
Double_t x

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