GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
TDecay.h
Go to the documentation of this file.
1// Author: Ryan Dunlop 09/15
2#ifndef TDECAY_H
3#define TDECAY_H
4
5/** \addtogroup Fitting Fitting & Analysis
6 * @{
7 */
8
9#include "TNamed.h"
10#include "TFitResult.h"
11#include "TFitResultPtr.h"
12#include "TGraph.h"
13#include "TF1.h"
14#include "TH1.h"
15
16class TVirtualDecay;
17class TSingleDecay;
18class TDecay;
19class TDecayChain;
20
21class TDecayFit : public TF1 {
22public:
23 TDecayFit() = default;
24 TDecayFit(const char* name, const char* formula, Double_t xmin = 0, Double_t xmax = 1)
25 : TF1(name, formula, xmin, xmax)
26 {
28 }
29 TDecayFit(const char* name, Double_t xmin, Double_t xmax, Int_t npar) : TF1(name, xmin, xmax, npar)
30 {
32 }
33 TDecayFit(const char* name, const ROOT::Math::ParamFunctor& f, Double_t xmin = 0, Double_t xmax = 1, Int_t npar = 0)
34 : TF1(name, f, xmin, xmax, npar)
35 {
37 }
38#if !defined(__CINT__) && !defined(__CLING__)
39 TDecayFit(const char* name, Double_t (*fcn)(Double_t*, Double_t*), Double_t xmin = 0, Double_t xmax = 1, Int_t npar = 0)
40 : TF1(name, fcn, xmin, xmax, npar)
41 {
43 }
44 TDecayFit(const char* name, Double_t (*fcn)(const Double_t*, const Double_t*), Double_t xmin = 0, Double_t xmax = 1, Int_t npar = 0)
45 : TF1(name, fcn, xmin, xmax, npar)
46 {
48 }
49#endif
50
51 template <class PtrObj, typename MemFn>
52 TDecayFit(const char* name, const PtrObj& p, MemFn memFn, Double_t xmin, Double_t xmax, Int_t npar,
53 const char* className = nullptr, const char* methodName = nullptr)
54 : TF1(name, p, memFn, xmin, xmax, npar, className, methodName)
55 {
57 }
58
59 template <typename Func>
60 TDecayFit(const char* name, Func f, Double_t xmin, Double_t xmax, Int_t npar, const char* className = nullptr)
61 : TF1(name, f, xmin, xmax, npar, className)
62 {
64 }
65 TDecayFit(const TDecayFit&) = default;
66 TDecayFit(TDecayFit&&) noexcept = default;
67 TDecayFit& operator=(const TDecayFit&) = default;
68 TDecayFit& operator=(TDecayFit&&) noexcept = default;
69 ~TDecayFit() = default;
70
71 void SetDecay(TVirtualDecay* decay);
72 TVirtualDecay* GetDecay() const;
73 void DrawComponents() const; // *MENU*
74
75 void Print(Option_t* opt = "") const override;
76 void UpdateResiduals(TH1* hist);
77 void DrawResiduals(); // *MENU*
78 TGraph* GetResiduals() { return &fResiduals; }
79 TFitResultPtr Fit(TH1* hist, Option_t* opt = "");
80
81private:
82 void DefaultGraphs();
83
84 TVirtualDecay* fDecay{nullptr}; // VirtualDecay that made this fit
85 TGraph fResiduals; // Last histogram fit by this function
86
87 /// \cond CLASSIMP
88 ClassDefOverride(TDecayFit, 1); // NOLINT(readability-else-after-return)
89 /// \endcond
90};
91
92class TVirtualDecay : public TNamed {
93public:
94 TVirtualDecay() = default;
95 TVirtualDecay(const TVirtualDecay&) = default;
96 TVirtualDecay(TVirtualDecay&&) noexcept = default;
97 TVirtualDecay& operator=(const TVirtualDecay&) = default;
98 TVirtualDecay& operator=(TVirtualDecay&&) noexcept = default;
99 ~TVirtualDecay() = default;
100
101 virtual void DrawComponents(Option_t* opt = "", Bool_t color_flag = true);
102 void Print(Option_t* opt = "") const override = 0;
103
104private:
105 virtual const TDecayFit* GetFitFunction() = 0;
106
107 /// \cond CLASSIMP
108 ClassDefOverride(TVirtualDecay, 1) // NOLINT(readability-else-after-return)
109 /// \endcond
110};
111
113 friend class TDecayChain;
114 friend class TDecayFit;
115 friend class TDecay;
116 // friend class TDecay;
117
118public:
120 : fDetectionEfficiency(1.0)
121 {
122 }
123 TSingleDecay(UInt_t generation, TSingleDecay* parent, Double_t tlow = 0, Double_t thigh = 10);
124 explicit TSingleDecay(TSingleDecay* parent, Double_t tlow = 0, Double_t thigh = 10);
125 TSingleDecay(const TSingleDecay&) = default;
126 TSingleDecay(TSingleDecay&&) noexcept = default;
127 TSingleDecay& operator=(const TSingleDecay&) = default;
128 TSingleDecay& operator=(TSingleDecay&&) noexcept = default;
129 ~TSingleDecay();
130
131 ///// TF1 Helpers ////
132 UInt_t GetGeneration() const { return fGeneration; }
133 Double_t GetHalfLife() const
134 {
135 return (fDecayFunc->GetParameter(1) > 0.) ? std::log(2.0) / fDecayFunc->GetParameter(1) : 0.0;
136 }
137 Double_t GetDecayRate() const { return fDecayFunc->GetParameter(1); }
138 Double_t GetIntensity() const { return fDecayFunc->GetParameter(0); }
139 Double_t GetEfficiency() const { return fDetectionEfficiency; }
140 Double_t GetHalfLifeError() const
141 {
142 return (GetDecayRate() > 0.) ? GetHalfLife() * GetDecayRateError() / GetDecayRate() : 0.0;
143 }
144 Double_t GetDecayRateError() const { return fDecayFunc->GetParError(1); }
145 Double_t GetIntensityError() const { return fDecayFunc->GetParError(0); }
146 void SetHalfLife(const Double_t& halflife)
147 {
148 fDecayFunc->SetParameter(1, std::log(2.0) / halflife);
149 UpdateDecays();
150 }
151 void SetDecayRate(const Double_t& decayrate)
152 {
153 fDecayFunc->SetParameter(1, decayrate);
154 UpdateDecays();
155 }
156 void SetIntensity(const Double_t& intens)
157 {
158 fDecayFunc->SetParameter(0, intens);
159 UpdateDecays();
160 }
161 void SetEfficiency(const Double_t& eff) { fDetectionEfficiency = eff; }
162 void FixHalfLife(const Double_t& halflife)
163 {
164 fDecayFunc->FixParameter(1, std::log(2) / halflife);
165 UpdateDecays();
166 }
168 {
169 fDecayFunc->FixParameter(1, GetHalfLife());
170 UpdateDecays();
171 }
172 void FixDecayRate(const Double_t& decayrate)
173 {
174 fDecayFunc->FixParameter(1, decayrate);
175 UpdateDecays();
176 }
178 {
179 fDecayFunc->FixParameter(0, GetDecayRate());
180 UpdateDecays();
181 }
182 void FixIntensity(const Double_t& intensity) { fDecayFunc->FixParameter(0, intensity); }
183 void FixIntensity() { fDecayFunc->FixParameter(0, GetIntensity()); }
184 void SetHalfLifeLimits(const Double_t& low, const Double_t& high);
185 void SetIntensityLimits(const Double_t& low, const Double_t& high);
186 void SetDecayRateLimits(const Double_t& low, const Double_t& high);
187 void GetHalfLifeLimits(Double_t& low, Double_t& high) const;
188 void GetIntensityLimits(Double_t& low, Double_t& high) const;
189 void GetDecayRateLimits(Double_t& low, Double_t& high) const;
190 void ReleaseHalfLife() { fDecayFunc->ReleaseParameter(1); }
191 void ReleaseDecayRate() { fDecayFunc->ReleaseParameter(1); }
192 void ReleaseIntensity() { fDecayFunc->ReleaseParameter(0); }
193 void Draw(Option_t* option = "") override;
194 Double_t Eval(Double_t t);
195 Double_t EvalPar(const Double_t* x, const Double_t* par = nullptr);
196 TFitResultPtr Fit(TH1* fithist, Option_t* opt = "");
197 void Fix();
198 void Release();
199 void SetRange(Double_t tlow, Double_t thigh);
200 void SetName(const char* name) override;
201 void SetLineColor(Color_t color) { fTotalDecayFunc->SetLineColor(color); }
202 Color_t GetLineColor() const { return fTotalDecayFunc->GetLineColor(); }
203 void SetMinimum(Double_t min)
204 {
205 fTotalDecayFunc->SetMinimum(min);
206 fDecayFunc->SetMinimum(min);
207 }
208 void SetMaximum(Double_t max)
209 {
210 fTotalDecayFunc->SetMaximum(max);
211 fDecayFunc->SetMaximum(max);
212 }
213
214 void SetDaughterDecay(TSingleDecay* daughter) { fDaughter = daughter; }
215 void SetParentDecay(TSingleDecay* parent) { fParent = parent; }
216 void SetTotalDecayParameters();
217 void SetDecayId(Int_t Id) { fUnId = Id; }
218 Int_t GetDecayId() const { return fUnId; }
219 Int_t GetChainId() const { return fChainId; }
220
221 const TDecayFit* GetDecayFunc() const { return fDecayFunc; }
223 {
224 SetTotalDecayParameters();
225 return fTotalDecayFunc;
226 }
227
228 TSingleDecay* GetParentDecay();
229 TSingleDecay* GetDaughterDecay();
230
231 Double_t ActivityFunc(Double_t* dim, Double_t* par);
232
233 void Print(Option_t* option = "") const override;
234
235private:
236 void SetDecayRateError(Double_t err) { fDecayFunc->SetParError(1, err); }
237 void SetIntensityError(Double_t err) { fDecayFunc->SetParError(0, err); }
238
239 void UpdateDecays();
240
241 void SetChainId(Int_t id) { fChainId = id; }
242
243 const TDecayFit* GetFitFunction() override
244 {
245 SetTotalDecayParameters();
246 return fTotalDecayFunc;
247 }
248
249 UInt_t fGeneration{0}; // Generation from the primary
250 Double_t fDetectionEfficiency{0.}; // The probability that this decay can be detected
251 TDecayFit* fDecayFunc{nullptr}; // Function describing decay
252 TDecayFit* fTotalDecayFunc{nullptr}; // Function used to access other fits
253 TSingleDecay* fParent{nullptr}; // Parent Decay
254 TSingleDecay* fDaughter{nullptr}; // Daughter Decay
255 TSingleDecay* fFirstParent{nullptr}; // FirstParent in the decay
256 Int_t fUnId{0}; // The Unique ID of the Decay
257 static UInt_t fCounter; // Helps set unique Id's
258 Int_t fChainId{-1}; // The chain that the single decay belongs to
259
260 /// \cond CLASSIMP
261 ClassDefOverride(TSingleDecay, 1) // NOLINT(readability-else-after-return)
262 /// \endcond
263};
264
266public:
267 TDecayChain() = default;
268 explicit TDecayChain(UInt_t generations);
269 TDecayChain(const TDecayChain&) = default;
270 TDecayChain(TDecayChain&&) noexcept = default;
271 TDecayChain& operator=(const TDecayChain&) = default;
272 TDecayChain& operator=(TDecayChain&&) noexcept = default;
273 ~TDecayChain();
274
275 TSingleDecay* GetDecay(UInt_t generation);
276 Double_t Eval(Double_t t) const;
277 void Draw(Option_t* opt = "") override;
278 Int_t Size() const { return static_cast<Int_t>(fDecayChain.size()); }
279
280 void Print(Option_t* option = "") const override;
281
282 void SetChainParameters();
283 void SetRange(Double_t xlow, Double_t xhigh);
285 {
286 SetChainParameters();
287 return fChainFunc;
288 }
289 void DrawComponents(Option_t* opt = "", Bool_t color_flag = true) override;
290 TFitResultPtr Fit(TH1* fithist, Option_t* opt = "");
291 Double_t EvalPar(const Double_t* x, const Double_t* par = nullptr);
292
293 Int_t GetChainId() const { return fChainId; }
294
295private:
296 void AddToChain(TSingleDecay* decay);
297 Double_t ChainActivityFunc(Double_t* dim, Double_t* par);
298 static UInt_t fChainCounter;
299 const TDecayFit* GetFitFunction() override
300 {
301 SetChainParameters();
302 return fChainFunc;
303 }
304
305 std::vector<TSingleDecay*> fDecayChain; // The Decays in the Decay Chain
306 TDecayFit* fChainFunc{nullptr}; // Function describing the total chain activity
307 Int_t fChainId{-1};
308
309 /// \cond CLASSIMP
310 ClassDefOverride(TDecayChain, 1) // NOLINT(readability-else-after-return)
311 /// \endcond
312};
313
314////////////////////////////////////////////////////////////////////////
315///
316/// \class TDecay
317///
318/// TDecay is a class for fitting halflives during nuclear decay
319/// A TDecay consists of multiple TDecayChains, where a TDecayChain
320/// is starts at a specific nucleus which has a population before the
321/// decay fit takes place. This could be a nucleus with a daughter.
322/// One TDecayChain would consist of just the daughter while the
323/// the other decay chain would be the parent and daughter.
324/// TDecayChains are made up of multiple TSingleDecays which holds
325/// the nucleus specific information such as name, id, halflife and
326/// intensity. When any of the above classes are fit to a histogram,
327/// they use a TDecayFit. The TDecayFit is a a TF1 with extra information
328/// such as the class that was used to create the TDecayFit. Furthermore,
329/// the function DrawComponents() can be used to draw the activites of the
330/// individual nuclei involved in the TDecayFit.
331///
332////////////////////////////////////////////////////////////////////////
333
334class TDecay : public TVirtualDecay {
335public:
336 TDecay() = default;
337 explicit TDecay(std::vector<TDecayChain*> chainList);
338 TDecay(const TDecay&) = default;
339 TDecay(TDecay&&) noexcept = default;
340 TDecay& operator=(const TDecay&) = default;
341 TDecay& operator=(TDecay&&) noexcept = default;
342 ~TDecay() = default;
343
344 void AddChain(TDecayChain* chain) { fChainList.push_back(chain); }
345 Double_t DecayFit(Double_t* dim, Double_t* par);
346 TDecayChain* GetChain(UInt_t idx);
347
348 void SetHalfLife(Int_t Id, Double_t halflife);
349 void SetHalfLifeLimits(Int_t Id, Double_t low, Double_t high);
350 void SetDecayRateLimits(Int_t Id, Double_t low, Double_t high);
351 void FixHalfLife(Int_t Id, Double_t halflife)
352 {
353 SetHalfLife(Id, halflife);
354 SetHalfLifeLimits(Id, halflife, halflife);
355 }
356 TFitResultPtr Fit(TH1* fithist, Option_t* opt = "");
357
358 void Print(Option_t* opt = "") const override;
359 void PrintMap() const;
360 TDecayFit* GetFitFunc() { return fFitFunc; }
361 void SetBackground(Double_t background) { fFitFunc->SetParameter(0, background); }
362 Double_t GetBackground() const { return fFitFunc->GetParameter(0); }
363 Double_t GetBackgroundError() const { return fFitFunc->GetParError(0); }
364 void SetRange(Double_t xlow, Double_t xhigh);
365 void DrawComponents(Option_t* opt = "", Bool_t color_flag = true) override;
366 void Draw(Option_t* opt = "") override;
367 void DrawBackground(Option_t* opt = "");
368 void FixBackground(const Double_t& background) { fFitFunc->FixParameter(0, background); }
369 void FixBackground() { fFitFunc->FixParameter(0, GetBackground()); }
370 void SetBackgroundLimits(const Double_t& low, const Double_t& high) { fFitFunc->SetParLimits(0, low, high); }
371 void ReleaseBackground() { fFitFunc->ReleaseParameter(0); }
372
373 TGraph* GetResiduals() { return fFitFunc->GetResiduals(); }
374
375private:
376 void RemakeMap();
377 void SetParameters();
378 Double_t ComponentFunc(Double_t* dim, Double_t* par);
379 const TDecayFit* GetFitFunction() override { return fFitFunc; }
380
381 std::vector<TDecayChain*> fChainList;
382 TDecayFit* fFitFunc{nullptr};
383 std::map<Int_t, std::vector<TSingleDecay*>> fDecayMap; //
384
385 /// \cond CLASSIMP
386 ClassDefOverride(TDecay, 1) // NOLINT(readability-else-after-return)
387 /// \endcond
388};
389/*! @} */
390#endif
TH1D * hist
Definition UserFillObj.h:3
std::vector< TSingleDecay * > fDecayChain
Definition TDecay.h:305
TDecayChain(const TDecayChain &)=default
static UInt_t fChainCounter
Definition TDecay.h:298
TDecayChain()=default
const TDecayFit * GetFitFunction() override
Definition TDecay.h:299
TDecayChain(TDecayChain &&) noexcept=default
Int_t GetChainId() const
Definition TDecay.h:293
const TDecayFit * GetChainFunc()
Definition TDecay.h:284
TDecayFit(const char *name, Double_t(*fcn)(Double_t *, Double_t *), Double_t xmin=0, Double_t xmax=1, Int_t npar=0)
Definition TDecay.h:39
TDecayFit()=default
TVirtualDecay * fDecay
Definition TDecay.h:84
void DefaultGraphs()
Definition TDecay.cxx:54
void Print(Option_t *opt="") const override
Definition TDecay.cxx:27
TDecayFit(const char *name, Double_t xmin, Double_t xmax, Int_t npar)
Definition TDecay.h:29
TFitResultPtr Fit(TH1 *hist, Option_t *opt="")
Definition TDecay.cxx:39
TDecayFit(TDecayFit &&) noexcept=default
TGraph * GetResiduals()
Definition TDecay.h:78
TDecayFit(const char *name, const PtrObj &p, MemFn memFn, Double_t xmin, Double_t xmax, Int_t npar, const char *className=nullptr, const char *methodName=nullptr)
Definition TDecay.h:52
TGraph fResiduals
Definition TDecay.h:85
void SetDecay(TVirtualDecay *decay)
Definition TDecay.cxx:21
TDecayFit(const char *name, const char *formula, Double_t xmin=0, Double_t xmax=1)
Definition TDecay.h:24
void DrawComponents() const
Definition TDecay.cxx:15
TDecayFit(const char *name, Double_t(*fcn)(const Double_t *, const Double_t *), Double_t xmin=0, Double_t xmax=1, Int_t npar=0)
Definition TDecay.h:44
void UpdateResiduals(TH1 *hist)
Definition TDecay.cxx:61
TDecayFit(const char *name, const ROOT::Math::ParamFunctor &f, Double_t xmin=0, Double_t xmax=1, Int_t npar=0)
Definition TDecay.h:33
void DrawResiduals()
Definition TDecay.cxx:93
TDecayFit(const char *name, Func f, Double_t xmin, Double_t xmax, Int_t npar, const char *className=nullptr)
Definition TDecay.h:60
TDecayFit(const TDecayFit &)=default
TVirtualDecay * GetDecay() const
Definition TDecay.cxx:34
TDecay(TDecay &&) noexcept=default
Double_t GetBackgroundError() const
Definition TDecay.h:363
std::map< Int_t, std::vector< TSingleDecay * > > fDecayMap
Definition TDecay.h:383
std::vector< TDecayChain * > fChainList
Definition TDecay.h:381
void FixBackground()
Definition TDecay.h:369
void SetBackgroundLimits(const Double_t &low, const Double_t &high)
Definition TDecay.h:370
TGraph * GetResiduals()
Definition TDecay.h:373
Double_t GetBackground() const
Definition TDecay.h:362
const TDecayFit * GetFitFunction() override
Definition TDecay.h:379
TDecayFit * GetFitFunc()
Definition TDecay.h:360
TDecay()=default
TDecay(const TDecay &)=default
void FixHalfLife(Int_t Id, Double_t halflife)
Definition TDecay.h:351
void FixBackground(const Double_t &background)
Definition TDecay.h:368
void ReleaseBackground()
Definition TDecay.h:371
void SetBackground(Double_t background)
Definition TDecay.h:361
void SetDecayRate(const Double_t &decayrate)
Definition TDecay.h:151
Int_t GetChainId() const
Definition TDecay.h:219
void FixHalfLife(const Double_t &halflife)
Definition TDecay.h:162
void SetDaughterDecay(TSingleDecay *daughter)
Definition TDecay.h:214
Color_t GetLineColor() const
Definition TDecay.h:202
Int_t GetDecayId() const
Definition TDecay.h:218
static UInt_t fCounter
Definition TDecay.h:257
Double_t GetHalfLife() const
Definition TDecay.h:133
void ReleaseDecayRate()
Definition TDecay.h:191
Double_t GetDecayRateError() const
Definition TDecay.h:144
Double_t GetIntensityError() const
Definition TDecay.h:145
void SetMaximum(Double_t max)
Definition TDecay.h:208
void FixDecayRate(const Double_t &decayrate)
Definition TDecay.h:172
void SetChainId(Int_t id)
Definition TDecay.h:241
void FixIntensity()
Definition TDecay.h:183
Double_t GetEfficiency() const
Definition TDecay.h:139
void FixHalfLife()
Definition TDecay.h:167
const TDecayFit * GetDecayFunc() const
Definition TDecay.h:221
void SetDecayRateError(Double_t err)
Definition TDecay.h:236
void SetMinimum(Double_t min)
Definition TDecay.h:203
void SetIntensity(const Double_t &intens)
Definition TDecay.h:156
void SetLineColor(Color_t color)
Definition TDecay.h:201
const TDecayFit * GetFitFunction() override
Definition TDecay.h:243
void FixIntensity(const Double_t &intensity)
Definition TDecay.h:182
Double_t GetDecayRate() const
Definition TDecay.h:137
void SetDecayId(Int_t Id)
Definition TDecay.h:217
void ReleaseIntensity()
Definition TDecay.h:192
void SetHalfLife(const Double_t &halflife)
Definition TDecay.h:146
Double_t GetHalfLifeError() const
Definition TDecay.h:140
void SetParentDecay(TSingleDecay *parent)
Definition TDecay.h:215
Double_t GetIntensity() const
Definition TDecay.h:138
void FixDecayRate()
Definition TDecay.h:177
void SetIntensityError(Double_t err)
Definition TDecay.h:237
void ReleaseHalfLife()
Definition TDecay.h:190
const TDecayFit * GetTotalDecayFunc()
Definition TDecay.h:222
TSingleDecay(TSingleDecay &&) noexcept=default
void SetEfficiency(const Double_t &eff)
Definition TDecay.h:161
TSingleDecay(const TSingleDecay &)=default
virtual const TDecayFit * GetFitFunction()=0
void Print(Option_t *opt="") const override=0
virtual void DrawComponents(Option_t *opt="", Bool_t color_flag=true)
Definition TDecay.cxx:103
TVirtualDecay(TVirtualDecay &&) noexcept=default
TVirtualDecay()=default
TVirtualDecay(const TVirtualDecay &)=default