GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
TCal.h
Go to the documentation of this file.
1#ifndef TCAL_H
2#define TCAL_H
3
4/** \addtogroup Calibration
5 * @{
6 */
7
8/////////////////////////////////////////////////////////////////
9///
10/// \class TCal
11///
12/// This is an abstract class that contains the basic info
13/// about a calibration. Calibrations here are TGraphErrors
14/// that are fit, with the resulting fit function being the
15/// calibrating function.
16///
17/////////////////////////////////////////////////////////////////
18
19#include <map>
20#include <vector>
21#include <utility>
22
23#include "Varargs.h"
24#include "TNamed.h"
25#include "TH1.h"
26#include "TF1.h"
27#include "TList.h"
28#include "TFitResult.h"
29#include "TFitResultPtr.h"
30#include "TRandom.h"
31#include "TSpectrum.h"
32#include "TVirtualFitter.h"
33#include "TMath.h"
34#include "TCanvas.h"
35#include "TROOT.h"
36#include "TMultiGraph.h"
37#include "TGraphErrors.h"
38#include "TRef.h"
39
40#include "TChannel.h"
41#include "TNucleus.h"
42#include "TGRSITransition.h"
43
44class TCal : public TGraphErrors {
45public:
46 TCal();
47 TCal(const char* name, const char* title);
48 TCal(const TCal&);
49 TCal(TCal&&) noexcept = default;
50 TCal& operator=(const TCal&) = default;
51 TCal& operator=(TCal&&) noexcept = default;
52 ~TCal() = default;
53
54 // pure virtual functions
55 virtual Bool_t IsGroupable() const = 0;
56
57 void Copy(TObject& obj) const override;
58 virtual void WriteToChannel() const { Error("WriteToChannel", "Not defined for %s", ClassName()); }
59 virtual TF1* GetFitFunction() const { return fFitFunc; }
60 virtual void SetFitFunction(TF1* func) { fFitFunc = func; }
61 virtual std::vector<Double_t> GetParameters() const;
62 virtual Double_t GetParameter(size_t parameter) const;
63
64 TChannel* GetChannel() const;
65 Bool_t SetChannel(TChannel* chan);
66 Bool_t SetChannel(UInt_t chanNum);
67 void Print(Option_t* opt = "") const override;
68 void Clear(Option_t* opt = "") override;
69 // virtual void Draw(Option_t* chopt = "");
70
71 virtual void WriteToAllChannels(const std::string& mnemonic = "");
72
73 virtual void SetHist(TH1* hist);
74 TH1* GetHist() const { return fHist; }
75 virtual void SetNucleus(TNucleus* nuc, Option_t* opt = "");
76 virtual TNucleus* GetNucleus() const { return fNuc; }
77
78protected:
79 void InitTCal();
80
81private:
82 TRef fChan{nullptr}; ///< This points at the TChannel
83 TF1* fFitFunc{nullptr}; ///< Fit function representing calibration
84 TH1* fHist{nullptr}; ///< Histogram that was fit by the TPeak.
85 TNucleus* fNuc{nullptr}; ///< Nucleus that we are calibrating against
86
87 /// \cond CLASSIMP
88 ClassDefOverride(TCal, 2) // NOLINT(readability-else-after-return)
89 /// \endcond
90};
91/*! @} */
92#endif
TH1D * hist
Definition UserFillObj.h:3
Definition TCal.h:44
void Print(Option_t *opt="") const override
Definition TCal.cxx:145
virtual void SetFitFunction(TF1 *func)
Definition TCal.h:60
virtual Double_t GetParameter(size_t parameter) const
Definition TCal.cxx:99
TH1 * fHist
Histogram that was fit by the TPeak.
Definition TCal.h:84
void InitTCal()
Definition TCal.cxx:172
TH1 * GetHist() const
Definition TCal.h:74
virtual void WriteToChannel() const
Definition TCal.h:58
TF1 * fFitFunc
Fit function representing calibration.
Definition TCal.h:83
void Clear(Option_t *opt="") override
Definition TCal.cxx:137
TRef fChan
This points at the TChannel.
Definition TCal.h:82
Bool_t SetChannel(TChannel *chan)
Definition TCal.cxx:51
void Copy(TObject &obj) const override
Definition TCal.cxx:37
TNucleus * fNuc
Nucleus that we are calibrating against.
Definition TCal.h:85
TCal()
Definition TCal.cxx:3
virtual TF1 * GetFitFunction() const
Definition TCal.h:59
virtual void SetNucleus(TNucleus *nuc, Option_t *opt="")
Definition TCal.cxx:24
virtual void WriteToAllChannels(const std::string &mnemonic="")
Definition TCal.cxx:63
virtual Bool_t IsGroupable() const =0
TChannel * GetChannel() const
Definition TCal.cxx:121
TCal(TCal &&) noexcept=default
virtual void SetHist(TH1 *hist)
Definition TCal.cxx:129
virtual std::vector< Double_t > GetParameters() const
Definition TCal.cxx:81
virtual TNucleus * GetNucleus() const
Definition TCal.h:76