GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
TCal.cxx
Go to the documentation of this file.
1#include "TCal.h"
2
4{
5 /// Default constructor
6 InitTCal();
7}
8
9TCal::TCal(const char* name, const char* title)
10{
11 /// Constructor for naming the calibration
12 InitTCal();
13 SetNameTitle(name, title);
14 // fgraph->SetNameTitle(name,title);
15}
16
17TCal::TCal(const TCal& copy) : TGraphErrors(copy)
18{
19 /// Copy constructor
20 InitTCal();
21 copy.Copy(*this);
22}
23
24void TCal::SetNucleus(TNucleus* nuc, Option_t*)
25{
26 /// Sets the nucleus to be calibrated against
27 if(nuc == nullptr) {
28 Error("SetNucleus", "Nucleus does not exist");
29 return;
30 }
31 if(fNuc != nullptr) {
32 Warning("SetNucleus", "Overwriting nucleus: %s", fNuc->GetName());
33 }
34 fNuc = nuc;
35}
36
37void TCal::Copy(TObject& obj) const
38{
39 /// Copies the TCal.
40 static_cast<TCal&>(obj).fChan = fChan;
41 // Things to make deep copies of
42 if(fFitFunc != nullptr) {
43 *(static_cast<TCal&>(obj).fFitFunc) = *fFitFunc;
44 }
45
46 // Members to make shallow copies of
47 static_cast<TCal&>(obj).fNuc = fNuc;
48 TNamed::Copy(static_cast<TGraphErrors&>(obj));
49}
50
52{
53 /// Sets the channel being calibrated
54 if(chan == nullptr) {
55 Error("SetChannel", "TChannel does not exist");
56 return false;
57 }
58 // Set our TRef to point at the TChannel
59 fChan = chan;
60 return true;
61}
62
63void TCal::WriteToAllChannels(const std::string& mnemonic)
64{
65 /// Writes this calibration to all channels in the current TChannel Map
66 std::unordered_map<unsigned int, TChannel*>::iterator mapIt;
67 std::unordered_map<unsigned int, TChannel*>* chanMap = TChannel::GetChannelMap();
68 TChannel* origChan = GetChannel();
69 for(mapIt = chanMap->begin(); mapIt != chanMap->end(); mapIt++) {
70 if(mnemonic.empty() || (strncmp(mapIt->second->GetName(), mnemonic.c_str(), mnemonic.size()) == 0)) {
71 SetChannel(mapIt->second);
73 }
74 }
75 std::cout << std::endl;
76 if(origChan != nullptr) {
77 SetChannel(origChan);
78 }
79}
80
81std::vector<Double_t> TCal::GetParameters() const
82{
83 /// Returns all of the parameters in the current TCal.
84 std::vector<Double_t> paramList;
85 if(GetFitFunction() == nullptr) {
86 Error("GetParameters", "Function has not been fitted yet");
87 return paramList;
88 }
89
90 Int_t nParams = GetFitFunction()->GetNpar();
91
92 for(int i = 0; i < nParams; i++) {
93 paramList.push_back(GetParameter(i));
94 }
95
96 return paramList;
97}
98
99Double_t TCal::GetParameter(size_t parameter) const
100{
101 /// Returns the parameter at the index parameter
102 if(GetFitFunction() == nullptr) {
103 Error("GetParameter", "Function have not been fitted yet");
104 return 0;
105 }
106 return GetFitFunction()->GetParameter(parameter); // Root does all of the checking for us.
107}
108
109Bool_t TCal::SetChannel(UInt_t chanNum)
110{
111 /// Sets the channel for the calibration to the channel number channum. Returns
112 /// 0 if the channel does not exist
113 TChannel* chan = TChannel::GetChannelByNumber(static_cast<Int_t>(chanNum));
114 if(chan == nullptr) {
115 Error("SetChannel", "Channel Number %d does not exist in current memory.", chanNum);
116 return false;
117 }
118 return SetChannel(chan);
119}
120
122{
123 /// Gets the channel being pointed to by the TCal. Returns 0 if no channel
124 /// is set.
125 return static_cast<TChannel*>(fChan.GetObject()); // Gets the object pointed at by the TRef and casts it to a
126 // TChannel
127}
128
130{
131 /// Sets this histogram pointed to. TCal does NOT take ownership so you cannot delete this
132 /// histogram as long as you want to access the hist in the TCal/write it out. I will add this
133 /// functionality if I get annoyed enough with the way it is.
134 fHist = hist;
135}
136
137void TCal::Clear(Option_t*)
138{
139 /// Clears the calibration. Does not delete nuclei or channels.
140 fNuc = nullptr;
141 fChan = nullptr;
142 TGraphErrors::Clear();
143}
144
145void TCal::Print(Option_t*) const
146{
147 /// Prints calibration information
148 if(GetChannel() != nullptr) {
149 std::cout << "Channel Number: " << GetChannel()->GetNumber() << std::endl;
150 } else {
151 }
152
153 if(fFitFunc != nullptr) {
154 std::cout << std::endl
155 << "*******************************" << std::endl;
156 std::cout << " Fit:" << std::endl;
157 fFitFunc->Print();
158 std::cout << std::endl
159 << "*******************************" << std::endl;
160 } else {
161 std::cout << "Parameters: FIT NOT SET" << std::endl;
162 }
163
164 std::cout << "Nucleus: ";
165 if(GetNucleus() != nullptr) {
166 std::cout << GetNucleus()->GetName() << std::endl;
167 } else {
168 std::cout << "NOT SET" << std::endl;
169 }
170}
171
173{
174 /// Initiallizes the TCal.
175 /* fgraph = new TGraphErrors;*/
176 fFitFunc = nullptr;
177 fChan = nullptr;
178 fNuc = nullptr;
179 fHist = nullptr;
180 Clear();
181}
TH1D * hist
Definition UserFillObj.h:3
Definition TCal.h:44
void Print(Option_t *opt="") const override
Definition TCal.cxx:145
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
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
TChannel * GetChannel() const
Definition TCal.cxx:121
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
static std::unordered_map< unsigned int, TChannel * > * GetChannelMap()
Definition TChannel.h:72
int GetNumber() const
Definition TChannel.h:169
static TChannel * GetChannelByNumber(int temp_num)
Definition TChannel.cxx:482