GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
TPeakFitter.h
Go to the documentation of this file.
1#ifndef TPEAKFITTER_H
2#define TPEAKFITTER_H
3
4/** \addtogroup Fitting Fitting & Analysis
5 * @{
6 */
7
8#include <cstdarg>
9
10#include "TF1.h"
11#include "TFitResultPtr.h"
12#include "TGraph.h"
13
14#include "Globals.h"
15#include "TSinglePeak.h"
16
17/////////////////////////////////////////////////////////////////
18///
19/// \class TPeakFitter
20///
21/// This class is used to fit things that resemble "peaks" in data
22///
23/////////////////////////////////////////////////////////////////
24class TSinglePeak;
25
26class TPeakFitter : public TObject {
27public:
28 // ctors and dtors
29 TPeakFitter() : TPeakFitter(0., 0.) {}
30 TPeakFitter(const Double_t& rangeLow, const Double_t& rangeHigh);
31 TPeakFitter(const TPeakFitter&) = default;
32 TPeakFitter(TPeakFitter&&) noexcept = default;
33 TPeakFitter& operator=(const TPeakFitter&) = default;
34 TPeakFitter& operator=(TPeakFitter&&) noexcept = default;
35 ~TPeakFitter();
36
37 void AddPeak(TSinglePeak* peak)
38 {
39 fPeaksToFit.push_back(peak);
41 }
43 {
44 fPeaksToFit.remove(peak);
46 }
48 {
49 fPeaksToFit.clear();
51 }
52 size_t NumberOfPeaks() { return fPeaksToFit.size(); }
53 std::list<TSinglePeak*>& Peaks() { return fPeaksToFit; }
54 TSinglePeak* Peak(const size_t& index)
55 {
56 if(index >= fPeaksToFit.size()) {
57 std::cerr << "Only " << fPeaksToFit.size() << " peaks in this peak fitter, can't access peak #" << index << std::endl;
58 return nullptr;
59 }
60 auto it = fPeaksToFit.begin();
61 std::advance(it, index);
62 return *it;
63 }
64
65 void SetBackground(TF1* bg_to_fit) { fBGToFit = bg_to_fit; }
66 void InitializeParameters(TH1* fit_hist);
67 void InitializeBackgroundParameters(TH1* fit_hist);
68
69 void Print(Option_t* opt = "") const override;
70 void PrintParameters() const;
71
72 TF1* GetBackground() const { return fBGToFit; }
73 TF1* GetFitFunction() const { return fTotalFitFunction; }
74 void SetRange(const Double_t& low, const Double_t& high);
75 Int_t GetNParameters() const;
76 TFitResultPtr Fit(TH1* fit_hist, Option_t* opt = "");
77 void DrawPeaks(Option_t* = "") const;
78
79 void ResetInitFlag() { fInitFlag = false; }
80
81 void SetColorIndex(const int& index) { fColorIndex = index; }
82
83 static void VerboseLevel(EVerbosity val) { fVerboseLevel = val; }
85
86private:
88 void UpdatePeakParameters(const TFitResultPtr& fit_res, TH1* fit_hist);
89 Double_t DefaultBackgroundFunction(Double_t* dim, Double_t* par);
91 {
92 delete fTotalFitFunction;
93 fTotalFitFunction = nullptr;
94 }
95
96 std::list<TSinglePeak*> fPeaksToFit;
97 TF1* fBGToFit{nullptr};
98
99 TF1* fTotalFitFunction{nullptr};
100
101 Double_t fRangeLow{0.};
102 Double_t fRangeHigh{0.};
103
104 Double_t FitFunction(Double_t* dim, Double_t* par);
105 Double_t BackgroundFunction(Double_t* dim, Double_t* par);
106
107 bool fInitFlag{false};
108
109 TH1* fLastHistFit{nullptr};
110
111 int fColorIndex{0}; ///< this index is added to the colors kRed for the total function and kMagenta for the individual peaks
112
113 static EVerbosity fVerboseLevel; ///< Changes verbosity of code.
114
115 /// \cond CLASSIMP
116 ClassDefOverride(TPeakFitter, 2) // NOLINT(readability-else-after-return)
117 /// \endcond
118};
119/*! @} */
120#endif
EVerbosity
Definition Globals.h:143
static void VerboseLevel(EVerbosity val)
Definition TPeakFitter.h:83
void UpdatePeakParameters(const TFitResultPtr &fit_res, TH1 *fit_hist)
void RemovePeak(TSinglePeak *peak)
Definition TPeakFitter.h:42
TF1 * fTotalFitFunction
Definition TPeakFitter.h:99
void RemoveAllPeaks()
Definition TPeakFitter.h:47
void PrintParameters() const
TSinglePeak * Peak(const size_t &index)
Definition TPeakFitter.h:54
void SetRange(const Double_t &low, const Double_t &high)
std::list< TSinglePeak * > & Peaks()
Definition TPeakFitter.h:53
void Print(Option_t *opt="") const override
Int_t GetNParameters() const
void InitializeBackgroundParameters(TH1 *fit_hist)
void SetColorIndex(const int &index)
Definition TPeakFitter.h:81
Double_t FitFunction(Double_t *dim, Double_t *par)
Double_t fRangeLow
void InitializeParameters(TH1 *fit_hist)
int fColorIndex
this index is added to the colors kRed for the total function and kMagenta for the individual peaks
TF1 * GetBackground() const
Definition TPeakFitter.h:72
TPeakFitter(TPeakFitter &&) noexcept=default
void DrawPeaks(Option_t *="") const
static EVerbosity VerboseLevel()
Definition TPeakFitter.h:84
std::list< TSinglePeak * > fPeaksToFit
Definition TPeakFitter.h:96
void ResetInitFlag()
Definition TPeakFitter.h:79
Double_t BackgroundFunction(Double_t *dim, Double_t *par)
void SetBackground(TF1 *bg_to_fit)
Definition TPeakFitter.h:65
void ResetTotalFitFunction()
Definition TPeakFitter.h:90
TF1 * GetFitFunction() const
Definition TPeakFitter.h:73
Double_t fRangeHigh
void AddPeak(TSinglePeak *peak)
Definition TPeakFitter.h:37
static EVerbosity fVerboseLevel
Changes verbosity of code.
TFitResultPtr Fit(TH1 *fit_hist, Option_t *opt="")
size_t NumberOfPeaks()
Definition TPeakFitter.h:52
TPeakFitter(const TPeakFitter &)=default
Double_t DefaultBackgroundFunction(Double_t *dim, Double_t *par)
TH1 * fLastHistFit
TF1 * fBGToFit
Definition TPeakFitter.h:97
void UpdateFitterParameters()