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