GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
TSinglePeak.h
Go to the documentation of this file.
1#ifndef TSINGLEPEAK_H
2#define TSINGLEPEAK_H
3
4/** \addtogroup Fitting Fitting & Analysis
5 * @{
6 */
7
8#include <vector>
9#include <cstdarg>
10
11#include "TObject.h"
12#include "TF1.h"
13#include "TGraph.h"
14#include "TPeakFitter.h"
15
16/////////////////////////////////////////////////////////////////
17///
18/// \class TSinglePeak
19///
20/// This is the base class for anything that is used to fit
21/// things that resemble "peaks" in data.
22/// Any derived classes must implement functions to
23/// - initialize parameters and parameter names,
24/// - set the total function and peak function,
25/// The centroid of the peak should always be parameter 1.
26///
27/////////////////////////////////////////////////////////////////
28class TPeakFitter;
29
30class TSinglePeak : public TObject {
31public:
32 friend class TPeakFitter;
33 // ctors and dtors
34 TSinglePeak() = default;
35 TSinglePeak(const TSinglePeak&) = default;
36 TSinglePeak(TSinglePeak&&) noexcept = default;
37 TSinglePeak& operator=(const TSinglePeak&) = default;
38 TSinglePeak& operator=(TSinglePeak&&) noexcept = default;
39 ~TSinglePeak() = default;
40
41 virtual void InitParNames() {}
42 virtual void InitializeParameters(TH1*, const double&, const double&) {}
43 bool IsBackgroundParameter(const Int_t& par) const;
44 bool IsPeakParameter(const Int_t& par) const;
45 void SetListOfBGPar(const std::vector<bool>& list_of_bg_par) { fListOfBGPars = list_of_bg_par; }
46 Int_t GetNParameters() const;
47
48 void SetArea(const Double_t& area) { fArea = area; }
49 void SetAreaErr(const Double_t& area_err) { fAreaErr = area_err; }
50
51 Double_t Area() const { return fArea; }
52 Double_t AreaErr() const { return fAreaErr; }
53
54 virtual void Centroid(const Double_t& centroid) = 0;
55 virtual Double_t Centroid() const = 0;
56 virtual Double_t CentroidErr() const = 0;
57 virtual Double_t Width() const = 0;
58 virtual Double_t Sigma() const = 0;
59 virtual Double_t SigmaErr() const = 0;
60 virtual Double_t FWHM(); // not constant because we have to update the parameters of the peak function
61 virtual Double_t FWHMErr(); // not constant because we have to update the parameters of the peak function
62
63 void Print(Option_t* = "") const override;
64 void Draw(Option_t* opt = "") override;
65 virtual void DrawBackground(Option_t* opt = "")
66 {
67 if(fGlobalBackground != nullptr) { fGlobalBackground->Draw(opt); }
68 }
69 virtual void DrawComponents(Option_t* opt = "");
70 virtual void PrintParameters() const;
71
72 TF1* GetFitFunction() const { return fTotalFunction; }
73 TF1* GetPeakFunction() const { return fPeakFunction; }
74 TF1* GetGlobalBackground() const { return fGlobalBackground; }
76 void SetGlobalBackground(TF1* background)
77 {
78 fGlobalBackground = background;
79 fGlobalBackground->SetLineStyle(kDashed);
80 }
81
84
85 Double_t GetChi2() const { return fChi2; }
86 Double_t GetNDF() const { return fNDF; }
87 Double_t GetReducedChi2() const { return fChi2 / fNDF; }
88
89 bool ParameterSetByUser(int par);
90
91 void SetLineColor(Color_t color) { fTotalFunction->SetLineColor(color); }
92 void SetLineStyle(Style_t style) { fTotalFunction->SetLineStyle(style); }
93
94 Color_t GetLineColor() { return fTotalFunction->GetLineColor(); }
95 Style_t GetLineStyle() { return fTotalFunction->GetLineStyle(); }
96
97protected:
98 Double_t TotalFunction(Double_t* dim, Double_t* par);
99 virtual Double_t BackgroundFunction(Double_t*, Double_t*) { return 0.0; }
100 virtual Double_t PeakFunction(Double_t*, Double_t*) { return 0.0; }
101 virtual Double_t PeakOnGlobalFunction(Double_t* dim, Double_t* par);
102
103 void SetChi2(const Double_t& chi2) { fChi2 = chi2; }
104 void SetNDF(const Int_t& ndf) { fNDF = ndf; }
105
106 void SetFitFunction(TF1* function) { fTotalFunction = function; }
107 void SetPeakFunction(TF1* function) { fPeakFunction = function; }
108
109private:
110 TF1* fTotalFunction{nullptr};
111 TF1* fBackgroundFunction{nullptr};
112 TF1* fGlobalBackground{nullptr};
113 TF1* fPeakOnGlobal{nullptr};
114 TF1* fPeakFunction{nullptr};
115
116 std::vector<bool> fListOfBGPars;
117 Double_t fArea{-0.1};
118 Double_t fAreaErr{0.0};
119 Double_t fChi2{std::numeric_limits<Double_t>::quiet_NaN()};
120 Int_t fNDF{0};
121
122public:
123 /// \cond CLASSIMP
124 ClassDefOverride(TSinglePeak, 2) // NOLINT(readability-else-after-return)
125 /// \endcond
126};
127/*! @} */
128#endif
Int_t GetNParameters() const
TF1 * GetFitFunction() const
Definition TSinglePeak.h:72
virtual Double_t CentroidErr() const =0
void Draw(Option_t *opt="") override
TSinglePeak()=default
Double_t fChi2
void SetPeakFunction(TF1 *function)
void UpdateBackgroundParameters()
TSinglePeak(const TSinglePeak &)=default
virtual Double_t FWHMErr()
Double_t fAreaErr
bool ParameterSetByUser(int par)
virtual Double_t Sigma() const =0
Color_t GetLineColor()
Definition TSinglePeak.h:94
Double_t GetReducedChi2() const
Definition TSinglePeak.h:87
TF1 * fGlobalBackground
void SetLineStyle(Style_t style)
Definition TSinglePeak.h:92
virtual Double_t BackgroundFunction(Double_t *, Double_t *)
Definition TSinglePeak.h:99
TF1 * fPeakFunction
TF1 * fTotalFunction
void SetGlobalBackground(TF1 *background)
Definition TSinglePeak.h:76
void SetLineColor(Color_t color)
Definition TSinglePeak.h:91
bool IsPeakParameter(const Int_t &par) const
Double_t fArea
TSinglePeak(TSinglePeak &&) noexcept=default
std::vector< bool > fListOfBGPars
Double_t AreaErr() const
Definition TSinglePeak.h:52
virtual void DrawComponents(Option_t *opt="")
virtual Double_t Centroid() const =0
Double_t GetChi2() const
Definition TSinglePeak.h:85
TF1 * GetPeakFunction() const
Definition TSinglePeak.h:73
virtual Double_t PeakOnGlobalFunction(Double_t *dim, Double_t *par)
bool IsBackgroundParameter(const Int_t &par) const
virtual Double_t PeakFunction(Double_t *, Double_t *)
Style_t GetLineStyle()
Definition TSinglePeak.h:95
Double_t TotalFunction(Double_t *dim, Double_t *par)
virtual Double_t SigmaErr() const =0
virtual void Centroid(const Double_t &centroid)=0
virtual Double_t FWHM()
void SetAreaErr(const Double_t &area_err)
Definition TSinglePeak.h:49
void SetArea(const Double_t &area)
Definition TSinglePeak.h:48
void SetFitFunction(TF1 *function)
TF1 * GetGlobalBackground() const
Definition TSinglePeak.h:74
TF1 * fBackgroundFunction
virtual void PrintParameters() const
virtual Double_t Width() const =0
void SetNDF(const Int_t &ndf)
TF1 * fPeakOnGlobal
TF1 * GetBackgroundFunction()
virtual void InitializeParameters(TH1 *, const double &, const double &)
Definition TSinglePeak.h:42
Double_t GetNDF() const
Definition TSinglePeak.h:86
virtual void DrawBackground(Option_t *opt="")
Definition TSinglePeak.h:65
void Print(Option_t *="") const override
void SetListOfBGPar(const std::vector< bool > &list_of_bg_par)
Definition TSinglePeak.h:45
void SetChi2(const Double_t &chi2)
virtual void InitParNames()
Definition TSinglePeak.h:41
void UpdatePeakParameters()
Double_t Area() const
Definition TSinglePeak.h:51