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