GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
TChannel.h
Go to the documentation of this file.
1#ifndef TCHANNEL_H
2#define TCHANNEL_H
3
4/** \addtogroup Sorting
5 * @{
6 */
7
8/*
9 * Author: P.C. Bender, <pcbend@gmail.com>
10 *
11 * Please indicate changes with your initials.
12 *
13 *
14 */
15
16//////////////////////////////////////////////////////////////////////////
17///
18/// \class TChannel
19///
20/// The TChannel is designed to hold all non-essential
21/// information of a TFragment (now the same as a hit; name, energy coeff, etc..)
22/// that would otherwise clog up the FragmentTree. The TChannel class
23/// contains a static map to every channel make retrieval fairly
24/// easy. The TChannel class also contains the ability to
25/// read and write a custom calibration file to set or
26/// save the TChannel information. Most of the information is read out
27/// of the ODB information in the MIDAS file.
28///
29//////////////////////////////////////////////////////////////////////////
30
31#include <string>
32#include <cmath>
33#include <utility>
34#include <unordered_map>
35
36#include "TNamed.h"
37#include "TTree.h"
38#include "TGraph.h"
39#include "TClassRef.h"
40#include "TPriorityValue.h"
41
42enum class EDigitizer : char;
43class TMnemonic;
44
45class TChannel : public TNamed {
46public:
47 static TChannel* GetChannel(unsigned int temp_address, bool warn = false);
48 static TChannel* GetChannelByNumber(int temp_num);
49 static TChannel* FindChannelByName(const char* ccName);
50 static std::vector<TChannel*> FindChannelByRegEx(const char* ccName);
51
52 TChannel();
53 explicit TChannel(const char*);
54 explicit TChannel(TChannel*);
55 TChannel(const TChannel&);
56 TChannel(TChannel&&) noexcept;
57
58 TChannel& operator=(const TChannel& rhs);
59 TChannel& operator=(TChannel&&) noexcept;
60
61 ~TChannel();
62
63 static size_t GetNumberOfChannels() { return fChannelMap->size(); }
64 static void AddChannel(TChannel*, Option_t* opt = "");
65 static int UpdateChannel(TChannel*, Option_t* opt = "");
66
67 static std::unordered_map<unsigned int, TChannel*>* GetChannelMap() { return fChannelMap; }
68 static std::unordered_map<unsigned int, int>* GetMissingChannelMap() { return fMissingChannelMap; }
69 static void DeleteAllChannels();
70
71 static bool CompareChannels(const TChannel*, const TChannel*);
72
74
75 static void SetMnemonicClass(const TClassRef& cls) { fMnemonicClass = cls; }
76 static TClassRef GetMnemonicClass() { return fMnemonicClass; }
77
78private:
79 unsigned int fAddress{0}; // The address of the digitizer
80 TPriorityValue<int> fIntegration{1}; // The charge integration setting
87
88 mutable int fDetectorNumber{-1};
89 mutable int fSegmentNumber{-1};
90 mutable int fCrystalNumber{-1};
91
93 TPriorityValue<double> fTimeDrift; ///< Time drift factor
95 static TClassRef fMnemonicClass;
96
97 TPriorityValue<std::vector<std::vector<Float_t>>> fENGCoefficients; ///< Energy calibration coeffs (low to high order)
99 TPriorityValue<std::vector<double>> fENGChi2; ///< Chi2 of the energy calibration
100 TPriorityValue<std::vector<Float_t>> fENGDriftCoefficents; ///< Energy drift coefficents (applied after energy calibration has been applied)
101 TPriorityValue<std::vector<double>> fCFDCoefficients; ///< CFD calibration coeffs (low to high order)
102 TPriorityValue<double> fCFDChi2; ///< Chi2 of the CFD calibration
103 TPriorityValue<std::vector<double>> fLEDCoefficients; ///< LED calibration coeffs (low to high order)
104 TPriorityValue<double> fLEDChi2; ///< Chi2 of LED calibration
105 TPriorityValue<std::vector<double>> fTIMECoefficients; ///< Time calibration coeffs (low to high order)
106 TPriorityValue<double> fTIMEChi2; ///< Chi2 of the Time calibration
107 TPriorityValue<std::vector<double>> fEFFCoefficients; ///< Efficiency calibration coeffs (low to high order)
108 TPriorityValue<double> fEFFChi2; ///< Chi2 of Efficiency calibration
110 TPriorityValue<TGraph> fEnergyNonlinearity; ///< Energy nonlinearity as TGraph, is used as E=E+GetEnergyNonlinearity(E), so y should be E(source)-calibration(peak)
112
114 bool InUse;
115 double BaseLine;
116 double TauDecay;
117 double TauRise;
118 };
119
121
122 static std::unordered_map<unsigned int, TChannel*>* fChannelMap; // A map to all of the channels based on address
123 static std::unordered_map<unsigned int, int>* fMissingChannelMap; // A map to all of the missing channels based on address
124 static std::unordered_map<int, TChannel*>* fChannelNumberMap; // A map of TChannels based on channel number
126 void AppendChannel(TChannel*);
127
128 void SetupEnergyNonlinearity(); // sort energy nonlinearity graph and set name/title
130
131 static std::vector<TChannel*> SortedChannels();
132
133public:
134 void SetName(const char* tmpName) override;
135 void SetAddress(unsigned int tmpadd);
137 {
138 if(fNumber == tmp) { return; }
139 // channel number has changed so we need to delete the old one and insert the new one
141 fNumber = tmp;
142 if((fNumber.Value() != 0) && (fChannelNumberMap->count(fNumber.Value()) == 0)) {
143 fChannelNumberMap->insert(std::make_pair(fNumber.Value(), this));
144 }
145 }
147 static void SetIntegration(const std::string& mnemonic, int tmpint, EPriority pr);
148 void SetStream(const TPriorityValue<int>& tmp) { fStream = tmp; }
150 static void SetDigitizerType(const std::string& mnemonic, const char* tmpstr, EPriority prio);
153
154 void SetDetectorNumber(int tempint) { fDetectorNumber = tempint; }
155 void SetSegmentNumber(int tempint) { fSegmentNumber = tempint; }
156 void SetCrystalNumber(int tempint) { fCrystalNumber = tempint; }
157
158 double GetTime(Long64_t timestamp, Float_t cfd, double energy) const;
159 int GetDetectorNumber() const;
160 int GetSegmentNumber() const;
161 int GetCrystalNumber() const;
162 const TMnemonic* GetMnemonic() const;
163 TClass* GetClassType() const;
164 void SetClassType(TClass* cl_type);
165
166 int GetNumber() const { return fNumber.Value(); }
167 unsigned int GetAddress() const { return fAddress; }
168 int GetIntegration() const { return fIntegration.Value(); }
169 int GetStream() const { return fStream.Value(); }
170 const char* GetDigitizerTypeString() const { return fDigitizerTypeString.c_str(); }
172 int GetTimeStampUnit() const { return fTimeStampUnit.Value(); }
173 Long64_t GetTimeOffset() const { return fTimeOffset.Value(); }
174 double GetTimeDrift() const { return fTimeDrift.Value(); }
175 // write the rest of the gettters/setters...
176
177 std::vector<double> GetAllENGChi2() const { return fENGChi2.Value(); }
178 double GetENGChi2(size_t range = 0) const { return fENGChi2.Value()[range]; }
179 double GetCFDChi2() const { return fCFDChi2.Value(); }
180 double GetLEDChi2() const { return fLEDChi2.Value(); }
181 double GetTIMEChi2() const { return fTIMEChi2.Value(); }
182 double GetEFFChi2() const { return fEFFChi2.Value(); }
183
185 static void SetUseCalFileIntegration(const std::string& mnemonic, bool flag, EPriority pr);
187
188 std::vector<std::vector<Float_t>> GetAllENGCoeff() const { return fENGCoefficients.Value(); }
189 std::vector<Float_t> GetENGCoeff(size_t range = 0) const { return fENGCoefficients.Value()[range]; }
190 std::vector<double> GetCFDCoeff() const { return fCFDCoefficients.Value(); }
191 std::vector<double> GetLEDCoeff() const { return fLEDCoefficients.Value(); }
192 std::vector<double> GetTIMECoeff() const { return fTIMECoefficients.Value(); }
193 std::vector<double> GetEFFCoeff() const { return fEFFCoefficients.Value(); }
194 std::vector<double> GetCTCoeff() const { return fCTCoefficients.Value(); }
196 double GetEnergyNonlinearity(double eng) const;
197 TGraph GetTimeNonlinearity() const { return fTimeNonlinearity.Value(); }
198 double GetTimeNonlinearity(Long64_t mytimestamp) const;
199 std::vector<std::pair<double, double>> GetENGRanges() const { return fENGRanges.Value(); }
200 std::pair<double, double> GetENGRange(size_t range) const { return fENGRanges.Value()[range]; }
201 std::vector<Float_t> GetENGDriftCoefficents() const { return fENGDriftCoefficents.Value(); }
202
203 void AddENGCoefficient(Float_t temp, size_t range = 0)
204 {
205 if(range >= fENGCoefficients.size()) { fENGCoefficients.resize(range + 1); }
206 fENGCoefficients.Address()->at(range).push_back(temp);
207 }
208 void AddENGDriftCoefficent(Float_t temp) { fENGDriftCoefficents.Address()->push_back(temp); }
209 void AddCFDCoefficient(double temp) { fCFDCoefficients.Address()->push_back(temp); }
210 void AddLEDCoefficient(double temp) { fLEDCoefficients.Address()->push_back(temp); }
211 void AddTIMECoefficient(double temp) { fTIMECoefficients.Address()->push_back(temp); }
212 void AddEFFCoefficient(double temp) { fEFFCoefficients.Address()->push_back(temp); }
213 void AddCTCoefficient(double temp) { fCTCoefficients.Address()->push_back(temp); }
214 void AddEnergyNonlinearityPoint(double x, double y) { fEnergyNonlinearity.Address()->SetPoint(fEnergyNonlinearity.Address()->GetN(), x, y); }
215 void AddTimeNonlinearityPoint(double x, double y) { fTimeNonlinearity.Address()->SetPoint(fTimeNonlinearity.Address()->GetN(), x, y); }
216
217 void ResizeENG(size_t size)
218 {
219 fENGCoefficients.resize(size);
220 fENGChi2.resize(size);
221 fENGRanges.resize(size);
222 }
223
224 void SetAllENGCoefficients(const TPriorityValue<std::vector<std::vector<Float_t>>>& tmp) { fENGCoefficients = tmp; }
225 void SetENGCoefficients(const std::vector<Float_t>& tmp, size_t range = 0)
226 {
227 if(range >= fENGCoefficients.size()) { fENGCoefficients.resize(range + 1); }
228 fENGCoefficients.Address()->at(range) = tmp;
229 }
230 void SetENGRanges(const TPriorityValue<std::vector<std::pair<double, double>>>& tmp) { fENGRanges = tmp; }
231 void SetENGRange(const std::pair<double, double>& tmp, const size_t& range)
232 {
233 if(range >= fENGRanges.size()) { fENGRanges.resize(range + 1); }
234 fENGRanges.Address()->at(range) = tmp;
235 }
236 void SetENGDriftCoefficents(const TPriorityValue<std::vector<Float_t>>& tmp) { fENGDriftCoefficents = tmp; }
237 void SetCFDCoefficients(const TPriorityValue<std::vector<double>>& tmp) { fCFDCoefficients = tmp; }
238 void SetLEDCoefficients(const TPriorityValue<std::vector<double>>& tmp) { fLEDCoefficients = tmp; }
239 void SetTIMECoefficients(const TPriorityValue<std::vector<double>>& tmp) { fTIMECoefficients = tmp; }
240 void SetEFFCoefficients(const TPriorityValue<std::vector<double>>& tmp) { fEFFCoefficients = tmp; }
241 void SetCTCoefficients(const TPriorityValue<std::vector<double>>& tmp) { fCTCoefficients = tmp; }
244
245 void SetAllENGChi2(const TPriorityValue<std::vector<double>>& tmp) { fENGChi2 = tmp; }
246 void SetENGChi2(const TPriorityValue<double>& tmp, const size_t& range = 0)
247 {
248 if(tmp.Priority() >= fENGChi2.Priority()) {
249 if(range >= fENGChi2.size()) { fENGChi2.resize(range + 1); }
250 fENGChi2.Address()->at(range) = tmp.Value();
251 }
252 }
253 void SetCFDChi2(const TPriorityValue<double>& tmp) { fCFDChi2 = tmp; }
254 void SetLEDChi2(const TPriorityValue<double>& tmp) { fLEDChi2 = tmp; }
255 void SetTIMEChi2(const TPriorityValue<double>& tmp) { fTIMEChi2 = tmp; }
256 void SetEFFChi2(const TPriorityValue<double>& tmp) { fEFFChi2 = tmp; }
257
258 void SetWaveRise(const double& temp)
259 {
260 WaveFormShape.TauRise = temp;
262 }
263 void SetWaveDecay(const double& temp)
264 {
265 WaveFormShape.TauDecay = temp;
267 }
268 void SetWaveBaseLine(double temp)
269 {
270 WaveFormShape.BaseLine = temp;
272 }
273
274 void SetUseWaveParam(bool temp = true) { WaveFormShape.InUse = temp; }
276
277 double GetWaveRise() const { return WaveFormShape.TauRise; }
278 double GetWaveDecay() const { return WaveFormShape.TauDecay; }
279 double GetWaveBaseLine() const { return WaveFormShape.BaseLine; }
280 bool UseWaveParam() const { return WaveFormShape.InUse; }
282
283 double CalibrateENG(double) const;
284 double CalibrateENG(double, int temp_int) const;
285 double CalibrateENG(int, int temp_int = 0) const;
286
287 double CalibrateCFD(double) const;
288 double CalibrateCFD(int) const;
289
290 double CalibrateLED(double) const;
291 double CalibrateLED(int) const;
292
293 double CalibrateTIME(double) const;
294 double CalibrateTIME(int) const;
295 double GetTZero(double tempd) const { return CalibrateTIME(tempd); }
296 double GetTZero(int tempi) const { return CalibrateTIME(tempi); }
297
298 double CalibrateEFF(double) const;
299
300 void DestroyCalibrations();
301
302 void DestroyENGCal();
303 void DestroyCFDCal();
304 void DestroyLEDCal();
305 void DestroyTIMECal();
306 void DestroyEFFCal();
307 void DestroyCTCal();
310
311 static Int_t ReadCalFromCurrentFile(Option_t* opt = "overwrite");
312 static Int_t ReadCalFromTree(TTree*, Option_t* opt = "overwrite");
313 static Int_t ReadCalFromFile(TFile* tempf, Option_t* opt = "overwrite");
314 static Int_t ReadCalFile(std::ifstream& infile);
315 static Int_t ReadCalFile(const char* filename = "");
316 static Int_t ParseInputData(const char* inputdata = "", Option_t* opt = "", EPriority prio = EPriority::kUser);
317 static void WriteCalFile(const std::string& outfilename = "");
318 static void WriteCTCorrections(const std::string& outfilename = "");
319 static void WriteCalBuffer(Option_t* opt = "");
320 static void ReadEnergyNonlinearities(TFile*, const char* graphName = "EnergyNonlinearity0x", bool all = false);
321 static void ReadTimeNonlinearities(TFile*, const char* graphName = "TimeNonlinearity0x", bool all = false);
322
323 void Print(Option_t* opt = "") const override;
324 void Clear(Option_t* opt = "") override;
325 // static void PrintAll(Option_t* opt = "");
326 std::string PrintToString(Option_t* opt = "") const;
327 std::string PrintCTToString(Option_t* opt = "") const;
328 void PrintCTCoeffs(Option_t* opt = "") const;
329
330 static int WriteToRoot(TFile* fileptr = nullptr);
331
332private:
333 // the follow is to make the custom streamer
334 // stuff play nice. pcb.
335 static std::string fFileData;
336 static void InitChannelInput();
337 static void SaveToSelf();
338 static void SaveToSelf(const char*);
339
340 static Int_t ReadFile(TFile* tempf);
341
342 /// \cond CLASSIMP
343 ClassDefOverride(TChannel, 6) // NOLINT(readability-else-after-return)
344 /// \endcond
345};
346/*! @} */
347#endif
EDigitizer
TPriorityValue< double > fCFDChi2
Chi2 of the CFD calibration.
Definition TChannel.h:102
TPriorityValue< std::string > fDigitizerTypeString
Definition TChannel.h:81
TPriorityValue< std::vector< double > > fCFDCoefficients
CFD calibration coeffs (low to high order)
Definition TChannel.h:101
void SetWaveDecay(const double &temp)
Definition TChannel.h:263
TPriorityValue< bool > fUseCalFileInt
Definition TChannel.h:86
void SetWaveRise(const double &temp)
Definition TChannel.h:258
void SetENGRanges(const TPriorityValue< std::vector< std::pair< double, double > > > &tmp)
Definition TChannel.h:230
static TChannel * GetDefaultChannel()
Definition TChannel.cxx:424
static void SaveToSelf()
double CalibrateENG(double) const
Definition TChannel.cxx:668
void AddLEDCoefficient(double temp)
Definition TChannel.h:210
TPriorityValue< std::vector< double > > fENGChi2
Chi2 of the energy calibration.
Definition TChannel.h:99
std::vector< double > GetEFFCoeff() const
Definition TChannel.h:193
static void SetMnemonicClass(const TClassRef &cls)
Definition TChannel.h:75
int GetStream() const
Definition TChannel.h:169
void AddCFDCoefficient(double temp)
Definition TChannel.h:209
void DestroyTimeNonlinearity()
Definition TChannel.cxx:615
static int WriteToRoot(TFile *fileptr=nullptr)
static std::unordered_map< unsigned int, int > * fMissingChannelMap
Definition TChannel.h:123
int GetDetectorNumber() const
void SetTimeNonlinearity(const TPriorityValue< TGraph > &tmp)
Definition TChannel.h:243
std::string PrintToString(Option_t *opt="") const
Definition TChannel.cxx:886
void SetTIMEChi2(const TPriorityValue< double > &tmp)
Definition TChannel.h:255
static std::unordered_map< unsigned int, TChannel * > * GetChannelMap()
Definition TChannel.h:67
std::vector< double > GetTIMECoeff() const
Definition TChannel.h:192
int GetCrystalNumber() const
void AddCTCoefficient(double temp)
Definition TChannel.h:213
static int UpdateChannel(TChannel *, Option_t *opt="")
Definition TChannel.cxx:409
double GetWaveDecay() const
Definition TChannel.h:278
void SetENGDriftCoefficents(const TPriorityValue< std::vector< Float_t > > &tmp)
Definition TChannel.h:236
static void WriteCalFile(const std::string &outfilename="")
static bool CompareChannels(const TChannel *, const TChannel *)
Definition TChannel.cxx:273
std::vector< double > GetCTCoeff() const
Definition TChannel.h:194
static TClassRef GetMnemonicClass()
Definition TChannel.h:76
TPriorityValue< Long64_t > fTimeOffset
Definition TChannel.h:92
TPriorityValue< TGraph > fTimeNonlinearity
Definition TChannel.h:111
static void WriteCalBuffer(Option_t *opt="")
TPriorityValue< std::vector< double > > fCTCoefficients
Cross talk coefficients.
Definition TChannel.h:109
TPriorityValue< TMnemonic * > fMnemonic
Definition TChannel.h:94
static TClassRef fMnemonicClass
Definition TChannel.h:95
void AddEFFCoefficient(double temp)
Definition TChannel.h:212
void SetCTCoefficients(const TPriorityValue< std::vector< double > > &tmp)
Definition TChannel.h:241
std::vector< Float_t > GetENGCoeff(size_t range=0) const
Definition TChannel.h:189
void SetTimeDrift(const TPriorityValue< double > &tmp)
Definition TChannel.h:152
static std::unordered_map< unsigned int, TChannel * > * fChannelMap
Definition TChannel.h:122
static void DeleteAllChannels()
Definition TChannel.cxx:282
void SetCFDChi2(const TPriorityValue< double > &tmp)
Definition TChannel.h:253
static Int_t ParseInputData(const char *inputdata="", Option_t *opt="", EPriority prio=EPriority::kUser)
void SetAllENGCoefficients(const TPriorityValue< std::vector< std::vector< Float_t > > > &tmp)
Definition TChannel.h:224
void SetLEDCoefficients(const TPriorityValue< std::vector< double > > &tmp)
Definition TChannel.h:238
void OverWriteChannel(TChannel *)
Definition TChannel.cxx:336
void DestroyLEDCal()
Definition TChannel.cxx:586
std::vector< double > GetLEDCoeff() const
Definition TChannel.h:191
void SetStream(const TPriorityValue< int > &tmp)
Definition TChannel.h:148
void SetLEDChi2(const TPriorityValue< double > &tmp)
Definition TChannel.h:254
std::pair< double, double > GetENGRange(size_t range) const
Definition TChannel.h:200
int fSegmentNumber
Definition TChannel.h:89
static void SaveToSelf(const char *)
void SetSegmentNumber(int tempint)
Definition TChannel.h:155
double CalibrateLED(double) const
Definition TChannel.cxx:752
void SetAddress(unsigned int tmpadd)
Definition TChannel.cxx:557
void SetTimeOffset(const TPriorityValue< Long64_t > &tmp)
Definition TChannel.h:151
double GetTime(Long64_t timestamp, Float_t cfd, double energy) const
void SetClassType(TClass *cl_type)
void SetTIMECoefficients(const TPriorityValue< std::vector< double > > &tmp)
Definition TChannel.h:239
void DestroyENGCal()
Definition TChannel.cxx:571
TPriorityValue< std::vector< double > > fTIMECoefficients
Time calibration coeffs (low to high order)
Definition TChannel.h:105
int GetSegmentNumber() const
unsigned int GetAddress() const
Definition TChannel.h:167
void SetEFFCoefficients(const TPriorityValue< std::vector< double > > &tmp)
Definition TChannel.h:240
void AddTimeNonlinearityPoint(double x, double y)
Definition TChannel.h:215
const TMnemonic * GetMnemonic() const
TPriorityValue< std::vector< Float_t > > fENGDriftCoefficents
Energy drift coefficents (applied after energy calibration has been applied)
Definition TChannel.h:100
static Int_t ReadCalFromCurrentFile(Option_t *opt="overwrite")
static TChannel * FindChannelByName(const char *ccName)
Definition TChannel.cxx:511
static Int_t ReadCalFile(std::ifstream &infile)
void SetUseCalFileIntegration(const TPriorityValue< bool > &tmp=TPriorityValue< bool >(true, EPriority::kUser))
Definition TChannel.h:184
TGraph GetEnergyNonlinearity() const
Definition TChannel.h:195
static std::unordered_map< int, TChannel * > * fChannelNumberMap
Definition TChannel.h:124
std::vector< Float_t > GetENGDriftCoefficents() const
Definition TChannel.h:201
WaveFormShapePar WaveFormShape
Definition TChannel.h:120
static Int_t ReadCalFromFile(TFile *tempf, Option_t *opt="overwrite")
static void ReadTimeNonlinearities(TFile *, const char *graphName="TimeNonlinearity0x", bool all=false)
void SetENGChi2(const TPriorityValue< double > &tmp, const size_t &range=0)
Definition TChannel.h:246
std::string PrintCTToString(Option_t *opt="") const
Definition TChannel.cxx:861
unsigned int fAddress
Definition TChannel.h:79
bool UseCalFileIntegration()
Definition TChannel.h:186
double CalibrateTIME(double) const
Definition TChannel.cxx:776
TPriorityValue< EDigitizer > fDigitizerType
Definition TChannel.h:82
std::vector< double > GetAllENGChi2() const
Definition TChannel.h:177
void AddENGDriftCoefficent(Float_t temp)
Definition TChannel.h:208
static std::vector< TChannel * > SortedChannels()
static TChannel * GetChannel(unsigned int temp_address, bool warn=false)
Definition TChannel.cxx:476
int fDetectorNumber
Definition TChannel.h:88
TChannel & operator=(const TChannel &rhs)
Definition TChannel.cxx:172
std::vector< double > GetCFDCoeff() const
Definition TChannel.h:190
double CalibrateEFF(double) const
Definition TChannel.cxx:792
static std::unordered_map< unsigned int, int > * GetMissingChannelMap()
Definition TChannel.h:68
double GetTZero(int tempi) const
Definition TChannel.h:296
double GetWaveRise() const
Definition TChannel.h:277
void SetupEnergyNonlinearity()
void AppendChannel(TChannel *)
Definition TChannel.cxx:375
void DestroyTIMECal()
Definition TChannel.cxx:592
Long64_t GetTimeOffset() const
Definition TChannel.h:173
double GetENGChi2(size_t range=0) const
Definition TChannel.h:178
TGraph GetTimeNonlinearity() const
Definition TChannel.h:197
TClass * GetClassType() const
double GetTimeDrift() const
Definition TChannel.h:174
TPriorityValue< int > fStream
Definition TChannel.h:85
TPriorityValue< double > fEFFChi2
Chi2 of Efficiency calibration.
Definition TChannel.h:108
static Int_t ReadCalFromTree(TTree *, Option_t *opt="overwrite")
TPriorityValue< double > fTimeDrift
Time drift factor.
Definition TChannel.h:93
TPriorityValue< std::vector< std::pair< double, double > > > fENGRanges
Range of energy calibrations.
Definition TChannel.h:98
static void WriteCTCorrections(const std::string &outfilename="")
void DestroyCFDCal()
Definition TChannel.cxx:580
void SetupTimeNonlinearity()
void DestroyEFFCal()
Definition TChannel.cxx:598
TPriorityValue< int > fIntegration
Definition TChannel.h:80
void Print(Option_t *opt="") const override
Definition TChannel.cxx:855
void DestroyCalibrations()
Definition TChannel.cxx:620
void PrintCTCoeffs(Option_t *opt="") const
Definition TChannel.cxx:840
void SetIntegration(const TPriorityValue< int > &tmp)
Definition TChannel.h:146
void ResizeENG(size_t size)
Definition TChannel.h:217
TPriorityValue< int > fTimeStampUnit
Definition TChannel.h:83
EDigitizer GetDigitizerType() const
Definition TChannel.h:171
void Clear(Option_t *opt="") override
Definition TChannel.cxx:432
TPriorityValue< double > fTIMEChi2
Chi2 of the Time calibration.
Definition TChannel.h:106
double CalibrateCFD(double) const
Definition TChannel.cxx:727
double GetTIMEChi2() const
Definition TChannel.h:181
void SetNumber(const TPriorityValue< int > &tmp)
Definition TChannel.h:136
void AddTIMECoefficient(double temp)
Definition TChannel.h:211
int GetNumber() const
Definition TChannel.h:166
int fCrystalNumber
Definition TChannel.h:90
TPriorityValue< double > fLEDChi2
Chi2 of LED calibration.
Definition TChannel.h:104
std::vector< std::vector< Float_t > > GetAllENGCoeff() const
Definition TChannel.h:188
TPriorityValue< std::vector< std::vector< Float_t > > > fENGCoefficients
Energy calibration coeffs (low to high order)
Definition TChannel.h:97
double GetLEDChi2() const
Definition TChannel.h:180
static TChannel * GetChannelByNumber(int temp_num)
Definition TChannel.cxx:499
void AddENGCoefficient(Float_t temp, size_t range=0)
Definition TChannel.h:203
static Int_t ReadFile(TFile *tempf)
static void InitChannelInput()
Definition TChannel.cxx:263
const char * GetDigitizerTypeString() const
Definition TChannel.h:170
int GetTimeStampUnit() const
Definition TChannel.h:172
void SetENGCoefficients(const std::vector< Float_t > &tmp, size_t range=0)
Definition TChannel.h:225
void SetName(const char *tmpName) override
Definition TChannel.cxx:253
void SetCrystalNumber(int tempint)
Definition TChannel.h:156
void DestroyCTCal()
Definition TChannel.cxx:604
TPriorityValue< TGraph > fEnergyNonlinearity
Energy nonlinearity as TGraph, is used as E=E+GetEnergyNonlinearity(E), so y should be E(source)-cali...
Definition TChannel.h:110
double GetWaveBaseLine() const
Definition TChannel.h:279
void SetDetectorNumber(int tempint)
Definition TChannel.h:154
void DestroyEnergyNonlinearity()
Definition TChannel.cxx:610
void AddEnergyNonlinearityPoint(double x, double y)
Definition TChannel.h:214
TPriorityValue< std::vector< double > > fLEDCoefficients
LED calibration coeffs (low to high order)
Definition TChannel.h:103
TPriorityValue< int > fNumber
Definition TChannel.h:84
TPriorityValue< std::vector< double > > fEFFCoefficients
Efficiency calibration coeffs (low to high order)
Definition TChannel.h:107
bool UseWaveParam() const
Definition TChannel.h:280
double GetCFDChi2() const
Definition TChannel.h:179
double GetEFFChi2() const
Definition TChannel.h:182
void SetUseWaveParam(bool temp=true)
Definition TChannel.h:274
static void AddChannel(TChannel *, Option_t *opt="")
Definition TChannel.cxx:294
void SetAllENGChi2(const TPriorityValue< std::vector< double > > &tmp)
Definition TChannel.h:245
static size_t GetNumberOfChannels()
Definition TChannel.h:63
void SetWaveParam(WaveFormShapePar temp)
Definition TChannel.h:275
static void ReadEnergyNonlinearities(TFile *, const char *graphName="EnergyNonlinearity0x", bool all=false)
void SetEnergyNonlinearity(const TPriorityValue< TGraph > &tmp)
Definition TChannel.h:242
static std::vector< TChannel * > FindChannelByRegEx(const char *ccName)
Definition TChannel.cxx:535
void SetENGRange(const std::pair< double, double > &tmp, const size_t &range)
Definition TChannel.h:231
void SetWaveBaseLine(double temp)
Definition TChannel.h:268
int GetIntegration() const
Definition TChannel.h:168
void SetDigitizerType(const TPriorityValue< std::string > &tmp)
static std::string fFileData
Definition TChannel.h:335
void SetCFDCoefficients(const TPriorityValue< std::vector< double > > &tmp)
Definition TChannel.h:237
double GetTZero(double tempd) const
Definition TChannel.h:295
std::vector< std::pair< double, double > > GetENGRanges() const
Definition TChannel.h:199
void SetEFFChi2(const TPriorityValue< double > &tmp)
Definition TChannel.h:256
WaveFormShapePar GetWaveParam() const
Definition TChannel.h:281
EPriority Priority() const
const T & Value() const
EPriority