GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
TPeak.h
Go to the documentation of this file.
1#ifndef TPEAK_H
2#define TPEAK_H
3
4/** \addtogroup Fitting Fitting & Analysis
5 * @{
6 */
7
8#include "TF1.h"
9#include "TFitResult.h"
10#include "TGraph.h"
11
12#include "TGRSIFit.h"
13
14/////////////////////////////////////////////////////////////////
15///
16/// \class TPeak
17///
18/// This Class is used to represent fitted data that is
19/// Gaussian like in nature (ie centroid and area).
20///
21/////////////////////////////////////////////////////////////////
22
23class TPeak : public TGRSIFit {
24 friend class TMultiPeak;
25
26public:
27 // ctors and dtors
28 TPeak(Double_t cent, Double_t xlow, Double_t xhigh, TF1* background = nullptr);
29 TPeak(); // I might make it so if you call this ctor, the TPeak yells at you since it's a fairly useless call anyway
30 TPeak(const TPeak& copy);
31 TPeak(TPeak&&) noexcept = default;
32 TPeak& operator=(const TPeak&) = default;
33 TPeak& operator=(TPeak&&) noexcept = default;
34 ~TPeak();
35
36protected:
37 void InitNames();
38
39public:
40 void Copy(TObject& obj) const override;
41 void SetCentroid(Double_t cent) { SetParameter("centroid", cent); }
42
43 Bool_t Fit(TH1* fitHist, Option_t* opt = ""); // Might switch this to TFitResultPtr
44 // Bool_t Fit(TH1* fithist = 0);
45
46 Double_t GetCentroid() const { return GetParameter("centroid"); }
47 Double_t GetCentroidErr() const { return GetParError(GetParNumber("centroid")); }
48 Double_t GetArea() const { return fArea; }
49 Double_t GetAreaErr() const { return fDArea; }
50 Double_t GetFWHM() const { return GetParameter("sigma") * 2.3548; }
51 Double_t GetFWHMErr() const { return GetParError(GetParNumber("sigma")) * 2.3548; }
52 Double_t GetIntegralArea();
53 Double_t GetIntegralArea(Double_t int_low, Double_t int_high);
54 Double_t GetIntegralAreaErr();
55 Double_t GetIntegralAreaErr(Double_t int_low, Double_t int_high);
56 /*
57 Double_t Fit(Option_t* opt = "");
58 Double_t Fit(TH1* hist, Option_t* opt = "");
59 Double_t Fit(const char* histname, Option_t* opt);
60 */
61 const TF1* GetFitFunction() const
62 {
63 return static_cast<const TF1*>(this);
64 } // I might move the fit functions to TGRSIFit, it's just a little tricky to initilize the function
65
66 Double_t Centroid() const { return GetCentroid(); }
67 Double_t CentroidErr() const { return GetCentroidErr(); }
68 Double_t Area() const { return GetArea(); }
69 Double_t AreaErr() const { return GetAreaErr(); }
70 Double_t FWHM() const { return GetFWHM(); }
71 Double_t FWHMErr() const { return GetFWHMErr(); }
72 Double_t IntegralArea() { return GetIntegralArea(); }
73 Double_t IntegralArea(Double_t int_low, Double_t int_high) { return GetIntegralArea(int_low, int_high); }
74 Double_t IntegralAreaErr() { return GetIntegralAreaErr(); }
75 Double_t IntegralAreaErr(Double_t int_low, Double_t int_high) { return GetIntegralAreaErr(int_low, int_high); }
76 const TF1* FitFunction() const { return GetFitFunction(); }
77
78protected:
79 void SetArea(Double_t area) { fArea = area; }
80 void SetAreaErr(Double_t areaErr) { fDArea = areaErr; }
81 void SetArea(Double_t area, Double_t areaErr)
82 {
83 SetArea(area);
84 SetAreaErr(areaErr);
85 }
86 void SetChi2(Double_t chi2) { fChi2 = chi2; }
87 void SetNdf(Double_t Ndf) { fNdf = Ndf; }
88
89public:
90 Bool_t InitParams(TH1* fitHist = nullptr) override;
91 TF1* Background() const { return fBackground; }
92 void DrawBackground(Option_t* opt = "SAME") const; // *MENU*
93 void DrawResiduals(); // *MENU*
94 void CheckArea();
95 void CheckArea(Double_t int_low, Double_t int_high);
96
97 static void SetLogLikelihoodFlag(Bool_t flag = true) { fLogLikelihoodFlag = flag; }
98 static Bool_t GetLogLikelihoodFlag() { return fLogLikelihoodFlag; }
99
100 static Bool_t CompareEnergy(const TPeak& lhs, const TPeak& rhs) { return lhs.GetCentroid() < rhs.GetCentroid(); }
101 static Bool_t CompareArea(const TPeak& lhs, const TPeak& rhs) { return lhs.GetArea() < rhs.GetArea(); }
102 static Bool_t CompareEnergy(const TPeak* lhs, const TPeak* rhs) { return lhs->GetCentroid() < rhs->GetCentroid(); }
103 static Bool_t CompareArea(const TPeak* lhs, const TPeak* rhs) { return lhs->GetArea() < rhs->GetArea(); }
104
105 static TPeak* GetLastFit() { return fLastFit; }
106
107 void Print(Option_t* opt = "") const override;
108 void Clear(Option_t* opt = "") override;
109
110private:
111 bool GoodStatus();
112 // Centroid will eventually be read from parameters
113 Double_t fArea{0.};
114 Double_t fDArea{0.};
115 Double_t fChi2{0.};
116 Double_t fNdf{0.};
117 Bool_t fOwnBgFlag{false};
118
119 static bool fLogLikelihoodFlag; //!<!
120 static TPeak* fLastFit; //!<!
121
122 TF1* fBackground{nullptr};
123 TGraph* fResiduals{nullptr};
124
125 /// \cond CLASSIMP
126 ClassDefOverride(TPeak, 2) // NOLINT(readability-else-after-return)
127 /// \endcond
128};
129/*! @} */
130#endif
Definition TPeak.h:23
Bool_t Fit(TH1 *fitHist, Option_t *opt="")
Definition TPeak.cxx:196
Double_t GetIntegralArea()
Definition TPeak.cxx:444
void SetNdf(Double_t Ndf)
Definition TPeak.h:87
void SetArea(Double_t area)
Definition TPeak.h:79
const TF1 * GetFitFunction() const
Definition TPeak.h:61
TF1 * Background() const
Definition TPeak.h:91
void InitNames()
Definition TPeak.cxx:97
Double_t GetIntegralAreaErr()
Definition TPeak.cxx:505
Double_t fDArea
Definition TPeak.h:114
Double_t GetFWHMErr() const
Definition TPeak.h:51
TPeak(TPeak &&) noexcept=default
void DrawResiduals()
Definition TPeak.cxx:393
Double_t IntegralAreaErr(Double_t int_low, Double_t int_high)
Definition TPeak.h:75
Double_t GetAreaErr() const
Definition TPeak.h:49
static void SetLogLikelihoodFlag(Bool_t flag=true)
Definition TPeak.h:97
void CheckArea()
Definition TPeak.cxx:538
static Bool_t CompareArea(const TPeak *lhs, const TPeak *rhs)
Definition TPeak.h:103
static Bool_t CompareEnergy(const TPeak &lhs, const TPeak &rhs)
Definition TPeak.h:100
Double_t Centroid() const
Definition TPeak.h:66
Double_t fChi2
Definition TPeak.h:115
Double_t CentroidErr() const
Definition TPeak.h:67
Bool_t fOwnBgFlag
Definition TPeak.h:117
TGraph * fResiduals
Definition TPeak.h:123
static Bool_t CompareEnergy(const TPeak *lhs, const TPeak *rhs)
Definition TPeak.h:102
Double_t GetCentroidErr() const
Definition TPeak.h:47
void Print(Option_t *opt="") const override
Definition TPeak.cxx:374
Double_t FWHMErr() const
Definition TPeak.h:71
static TPeak * GetLastFit()
Definition TPeak.h:105
Double_t IntegralArea()
Definition TPeak.h:72
static Bool_t CompareArea(const TPeak &lhs, const TPeak &rhs)
Definition TPeak.h:101
Double_t IntegralArea(Double_t int_low, Double_t int_high)
Definition TPeak.h:73
TF1 * fBackground
Definition TPeak.h:122
void DrawBackground(Option_t *opt="SAME") const
Definition TPeak.cxx:388
bool GoodStatus()
Definition TPeak.cxx:431
void Copy(TObject &obj) const override
Definition TPeak.cxx:111
Double_t fNdf
Definition TPeak.h:116
void SetAreaErr(Double_t areaErr)
Definition TPeak.h:80
Double_t Area() const
Definition TPeak.h:68
Double_t fArea
Definition TPeak.h:113
Double_t IntegralAreaErr()
Definition TPeak.h:74
const TF1 * FitFunction() const
Definition TPeak.h:76
Bool_t InitParams(TH1 *fitHist=nullptr) override
Definition TPeak.cxx:136
static TPeak * fLastFit
!
Definition TPeak.h:120
void SetCentroid(Double_t cent)
Definition TPeak.h:41
void Clear(Option_t *opt="") override
Definition TPeak.cxx:361
void SetChi2(Double_t chi2)
Definition TPeak.h:86
Double_t FWHM() const
Definition TPeak.h:70
static bool fLogLikelihoodFlag
!
Definition TPeak.h:119
static Bool_t GetLogLikelihoodFlag()
Definition TPeak.h:98
Double_t GetArea() const
Definition TPeak.h:48
TPeak()
Definition TPeak.cxx:74
Double_t GetFWHM() const
Definition TPeak.h:50
Double_t AreaErr() const
Definition TPeak.h:69
void SetArea(Double_t area, Double_t areaErr)
Definition TPeak.h:81
Double_t GetCentroid() const
Definition TPeak.h:46