ATLAS Offline Software
Loading...
Searching...
No Matches
T0MTHistos.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3*/
4
6
7#include <cmath>
8
11#include "GaudiKernel/MsgStream.h"
14#include "TDirectory.h"
15#include "TLine.h"
16#include "TRandom3.h"
17#include "list"
18#include "boost/thread/tss.hpp"
19
20namespace MuonCalib {
21
22 TRandom3* getTLSRandomGen()
23 {
24 static boost::thread_specific_ptr<TRandom3> rnd ATLAS_THREAD_SAFE;
25 TRandom3* random = rnd.get();
26 if (!random) {
27 random = new TRandom3();
28 rnd.reset(random);
29 }
30 return random;
31 }
32
34 inline Double_t mt_t0_fermi(Double_t *x, Double_t *par) {
35 // more convenient parameters
36 const Double_t &t(x[0]);
37 const Double_t &t_0(par[T0MTHistos::T0_PAR_NR_T0]), &T(par[T0MTHistos::T0_PAR_NR_T]), &back(par[T0MTHistos::T0_PAR_NR_BACK]),
39 // the formula
40 return (back + A / (1 + std::exp(-(t - t_0) / T)));
41 }
42
44 inline Double_t mt_tmax_fermi(Double_t *x, Double_t *par) {
45 // more convenient parameters
46 Double_t &t(x[0]);
49 // the formula
50 return (back + (std::exp(a + b * (t - t_0))) / (1 + std::exp((t - t_max) / T)));
51 }
52
54 // Fill.. //
56 void T0MTHistos::FillT(double t) { m_time->Fill(static_cast<Axis_t>(t)); }
57
59 // Initialize(int id, const T0MTSettings & settings) //
61 void T0MTHistos::Initialize(int id, const T0MTSettings *settings, const char *hname) {
62 m_settings = settings;
63#ifndef NDEBUG
64 MsgStream log(Athena::getMessageSvc(), "T0MTHistos");
65 if (log.level() <= MSG::VERBOSE) log << MSG::VERBOSE << "T0MTHistos::Initialize: called" << endmsg;
66#endif
67 char buf[100];
68 if (hname == nullptr)
69 snprintf(buf, 100, "t_spec_%d", id);
70 else
71 snprintf(buf, 100, "t_spec_%s", hname);
72#ifndef NDEBUG
73 if (log.level() <= MSG::VERBOSE) log << MSG::VERBOSE << "directory=" << gDirectory->GetName() << endmsg;
74#endif
75 m_time = std::make_unique<TH1F>(buf, "", settings->NBinsTime(), settings->TimeMin(), settings->TimeMax());
76 m_id = id;
77 if (settings->DrawDebugGraphs()) {
78#ifndef NDEBUG
79 if (log.level() <= MSG::VERBOSE) log << MSG::VERBOSE << "T0MTHistos::Initialize: debug directory created" << endmsg;
80#endif
81 TDirectory *cwd = gDirectory;
82 snprintf(buf, 100, "t0_tmax_dir_%d", id);
83 m_dir = gDirectory->mkdir(buf, buf);
84 cwd->cd();
85 } else {
86 m_dir = nullptr;
87 }
88 m_t0_ok = false;
89 m_tmax_ok = false;
90 } // end T0MTHistos::Initialize
91
93 // SetTSpec(int id, TH1F *spec) //
95 void T0MTHistos::SetTSpec(int id, TH1F *spec, const T0MTSettings *settings, bool copy_spec) {
96 m_settings = settings;
97 TDirectory *cwd = gDirectory;
98 if (copy_spec)
99 m_time = std::make_unique<TH1F>(*spec);
100 else
101 m_time.reset(spec);
102 m_id = id;
103 if (settings->DrawDebugGraphs()) {
104 char buffer[100];
105 snprintf(buffer, 100, "t0_tmax_dir_%d", id);
106 m_dir = gDirectory->mkdir(buffer, buffer);
107 cwd->cd();
108 }
109 } // end T0MTHistos::SetTSpec
110
112 // FitT0() //
115 if (m_time->GetEntries() < 1000) {
116 m_status_code = 2;
117 return false;
118 }
119 if (!NormalFit()) return false;
120 if (m_time->GetEntries() < 10000) {
121 m_status_code = 0;
122 return true;
123 }
124 if (m_settings->T0Settings()->ScrambleThreshold() > 0 && m_chi2 > m_settings->T0Settings()->ScrambleThreshold()) {
125 if (!T0Scramble()) {
126 m_status_code = 3;
127 return false;
128 }
129 }
130 if (m_settings->T0Settings()->SlicingThreshold() > 0 && m_chi2 > m_settings->T0Settings()->SlicingThreshold()) { TopSlicing(); }
131#ifndef NDEBUG
132 MsgStream log(Athena::getMessageSvc(), "T0MTHistos");
133 if (log.level() <= MSG::VERBOSE) log << MSG::VERBOSE << m_time->GetName() << " " << m_chi2 << endmsg;
134#endif
135 m_status_code = 0;
136 return true;
137 } // end T0MTHistos::FitT0
138
140 // FitTmax() //
143 TDirectory *cwd = gDirectory;
144 if (m_dir != nullptr) m_dir->cd();
145 if (!m_time) {
146 MsgStream log(Athena::getMessageSvc(), "T0MTHistos");
147 log << MSG::WARNING << "T0MTHistos::FitTmax: Class is not initialized!" << endmsg;
148 m_tmax_ok = false;
149 cwd->cd();
150 return false;
151 }
152 // check if t0-fit was successfull t0 is needed for tmax-pattern recognition
153 if (!m_t0_ok || m_t0_fermi == nullptr) {
154 MsgStream log(Athena::getMessageSvc(), "T0MTHistos");
155 log << MSG::WARNING << "T0MTHistos::FitTmax for tube " << m_id << ": No valid t0-value!" << endmsg;
156 m_tmax_ok = false;
157 cwd->cd();
158 return false;
159 }
160 // Create pattern Recognition Class
162 // perform pattern recognition
163 if (!rec.Initialize(m_time.get(), m_t0_fermi->GetParameter(T0_PAR_NR_T0), m_settings)) {
164 MsgStream log(Athena::getMessageSvc(), "T0MTHistos");
165 log << MSG::WARNING << "T0MTHistos::FitTmax for tube " << m_id << ": Pattern recognition failed!" << endmsg;
166 m_tmax_ok = false;
167 cwd->cd();
168 return false;
169 }
170 // create function object
171 char buffer[100];
172 snprintf(buffer, 100, "mt_tmax_fermi");
173 if (!m_tmax_fermi) {
174 m_tmax_fermi = std::make_unique<TF1>(buffer, mt_tmax_fermi, rec.GetFitRangeMin(), rec.GetFitRangeMax(), N_TMAX_FIT_PAR);
175 // set parameter names
176 m_tmax_fermi->SetParName(TMAX_PAR_NR_TMAX, "t_{max}");
177 m_tmax_fermi->SetParName(TMAX_PAR_NR_T, "T");
178 m_tmax_fermi->SetParName(TMAX_PAR_NR_BACK, "r_{b}");
179 m_tmax_fermi->SetParName(TMAX_PAR_NR_A, "a");
180 m_tmax_fermi->SetParName(TMAX_PAR_NR_B, "b");
181 // set fixed values
182 m_tmax_fermi->FixParameter(TMAX_PAR_NR_BACK, rec.GetBackground());
183 m_tmax_fermi->FixParameter(TMAX_PAR_NR_A, rec.GetA());
184 m_tmax_fermi->FixParameter(TMAX_PAR_NR_B, rec.GetB());
185 m_tmax_fermi->FixParameter(TMAX_PAR_NR_T0, m_t0_fermi->GetParameter(T0_PAR_NR_T0));
186 // set start values
187 m_tmax_fermi->SetParameter(TMAX_PAR_NR_TMAX, rec.GetEstimatedTMax());
188 m_tmax_fermi->SetParameter(TMAX_PAR_NR_T, 3.0);
189 }
190 // perform fit
191 if (m_dir != nullptr) {
192 m_tmax_fermi->SetLineColor(3);
193 m_tmax_fermi->Write();
194 }
195 m_tmax_fermi->SetLineColor(4);
196 std::string fitopt("LR");
197 if (m_settings->VerboseLevel() == 0) fitopt += "Q";
198 if (m_settings->AddFitfun()) {
199 fitopt += "+";
200 } else {
201 fitopt += "0";
202 }
203 m_time->Fit(m_tmax_fermi.get(), fitopt.c_str());
204 m_tmax_ok = true;
205 cwd->cd();
206 return true;
207 } // end T0MTHistos::FitTmax
208
210 // FitT0() //
213#ifndef NDEBUG
214 MsgStream log(Athena::getMessageSvc(), "T0MTHistos");
215 if (log.level() <= MSG::VERBOSE) log << MSG::VERBOSE << "T0MTHistos::FitT0(): called" << endmsg;
216#endif
217 TDirectory *cwd = gDirectory;
218 if (m_dir) m_dir->cd();
219 // check if class is initialized
220 if (!m_time) {
221#ifdef NDEBUG
222 MsgStream log(Athena::getMessageSvc(), "T0MTHistos");
223#endif
224 log << MSG::WARNING << "T0MTHistos::FitT0: Class is not initialized!" << endmsg;
225 m_t0_ok = false;
226 cwd->cd();
227 m_status_code = 3;
228 return false;
229 }
230 // create pattern recognition class
232 // perform pattern recognition
233 if (!rec.Initialize(m_time.get(), m_settings)) {
234 m_t0_ok = false;
235#ifdef NDEBUG
236 MsgStream log(Athena::getMessageSvc(), "T0MTHistos");
237#endif
238 log << MSG::WARNING << "T0MTHistos::FitT0 for tube " << m_id << ": Pattern recognition failed!" << endmsg;
239 cwd->cd();
240 m_status_code = 3;
241 return false;
242 }
243 m_t0_ok = true;
244 // create function class
245 char buffer[100];
246 // sprintf(buffer,"mt_t0_fermi_%d", m_id);
247 snprintf(buffer, 100, "mt_t0_fermi");
248 m_t0_fermi = std::make_unique<TF1>(buffer, mt_t0_fermi, rec.GetFitRangeMin(), rec.GetFitRangeMax(), N_T0_FIT_PAR);
249 // set parameter names
250 m_t0_fermi->SetParName(T0_PAR_NR_T0, "t_{0}");
251 m_t0_fermi->SetParName(T0_PAR_NR_T, "T");
252 m_t0_fermi->SetParName(T0_PAR_NR_BACK, "r_{b}");
253 m_t0_fermi->SetParName(T0_PAR_NR_A, "A");
254 // set fixed values
255 m_t0_fermi->FixParameter(T0_PAR_NR_BACK, rec.GetBackground());
256 m_t0_fermi->SetParameter(T0_PAR_NR_A, rec.GetHeight() - rec.GetBackground());
257 // set estimates as start values
258 m_t0_fermi->SetParameter(T0_PAR_NR_T0, rec.GetEstimatedT0());
259 // set resonable start value for T
260 m_t0_fermi->SetParameter(T0_PAR_NR_T, 3.0);
261 // perform fit - NOTE: The return value of the Fit function is not documented!
262 if (m_dir != nullptr) {
263 std::unique_ptr<TLine> ln = std::make_unique<TLine>(rec.GetFitRangeMin(), 0, rec.GetFitRangeMin(), m_time->GetMaximum());
264 ln->Write("t0_range_min");
265 ln = std::make_unique<TLine>(rec.GetFitRangeMax(), 0, rec.GetFitRangeMax(), m_time->GetMaximum());
266 ln->Write("t0_range_max");
267 m_t0_fermi->SetLineColor(3);
268 m_t0_fermi->Write();
269 }
270 m_t0_fermi->SetLineColor(2);
271 std::string fitopt("BLR");
272 if (m_settings->VerboseLevel() == 0) fitopt += "Q";
273 if (m_settings->AddFitfun()) {
274 fitopt += "+";
275 } else {
276 fitopt += "0";
277 }
278 m_time->Fit(m_t0_fermi.get(), fitopt.c_str(), "", rec.GetFitRangeMin(), rec.GetFitRangeMax());
279 if (m_settings->T0Settings()->UseTopChi2())
280 TopChi2();
281 else
282 m_chi2 = m_time->GetFunction("mt_t0_fermi")->GetChisquare() / m_time->GetFunction("mt_t0_fermi")->GetNDF();
283 m_t0_ok = true;
284 cwd->cd();
285 m_status_code = 0;
286 return true;
287 } // end T0MTHistos::NormalFit
288
290 // T0Scramble() //
293#ifndef NDEBUG
294 MsgStream log(Athena::getMessageSvc(), "T0MTHistos");
295 if (log.level() <= MSG::VERBOSE) log << MSG::VERBOSE << "Scrambling for " << m_time->GetName() << endmsg;
296#endif
297 std::string fitopt("BLR");
298 if (m_settings->VerboseLevel() == 0) fitopt += "Q";
299 if (m_settings->AddFitfun()) {
300 fitopt += "+";
301 } else {
302 fitopt += "0";
303 }
304 // create scrambled histogram
305 char scramhistname[100];
306 snprintf(scramhistname, 100, "%s_scram", m_time->GetName());
307 std::unique_ptr<TH1F> scramhist = std::make_unique<TH1F>(scramhistname, "scrambled histogram", m_time->GetSize() - 2,
308 m_time->GetXaxis()->GetXmin(), m_time->GetXaxis()->GetXmax());
309
310 for (int binnr = 0; binnr < m_time->GetSize(); binnr++) {
311 scramhist->SetBinContent(binnr, m_time->GetBinContent(binnr) + getTLSRandomGen()->Gaus(0, m_time->GetBinError(binnr)));
312 scramhist->SetBinError(binnr, m_time->GetBinError(binnr) * 1.41421356);
313 if (scramhist->GetBinContent(binnr) < 0) scramhist->SetBinContent(binnr, 0);
314 }
315 TDirectory *cwd = gDirectory;
316 if (m_dir) m_dir->cd();
317 MTT0PatternRecognition scramrec;
318 // perform pattern recognition
319 if (!scramrec.Initialize(scramhist.get(), m_settings)) {
320#ifdef NDEBUG
321 MsgStream log(Athena::getMessageSvc(), "T0MTHistos");
322#endif
323 log << MSG::WARNING << "T0MTHistos::FitT0 for tube " << m_id << ": Scrambed pattern recognition failed!" << endmsg;
324 cwd->cd();
325 return false;
326 }
327 char scrambuffer[100];
328 snprintf(scrambuffer, 100, "scrammt_t0_fermi");
329 std::unique_ptr<TF1> scramm_t0_fermi = std::make_unique<TF1>();
330 m_t0_fermi->Copy(*scramm_t0_fermi);
331 scramm_t0_fermi->SetName(scrambuffer);
332 scramm_t0_fermi->SetRange(scramrec.GetFitRangeMin(), scramrec.GetFitRangeMax());
333 // set parameter names
334 scramm_t0_fermi->SetParName(T0_PAR_NR_T0, "t_{0}");
335 scramm_t0_fermi->SetParName(T0_PAR_NR_T, "T");
336 scramm_t0_fermi->SetParName(T0_PAR_NR_BACK, "r_{b}");
337 scramm_t0_fermi->SetParName(T0_PAR_NR_A, "A");
338 // set fixed values
339 scramm_t0_fermi->FixParameter(T0_PAR_NR_BACK, scramrec.GetBackground());
340 scramm_t0_fermi->SetParameter(T0_PAR_NR_A, scramrec.GetHeight() - scramrec.GetBackground());
341 // set estimates as start values
342 scramm_t0_fermi->SetParameter(T0_PAR_NR_T0, scramrec.GetEstimatedT0());
343 // set resonable start value for T
344 scramm_t0_fermi->SetParameter(T0_PAR_NR_T, 3.0);
345 // perform fit - NOTE: The return value of the Fit function is not documented!
346 scramhist->Fit(scrambuffer, fitopt.c_str(), "", scramrec.GetFitRangeMin(), scramrec.GetFitRangeMax());
347 // set parameter for the new fit of the original histogram
348 m_time->GetListOfFunctions()->Clear();
349 // set fixed values
350 m_t0_fermi->FixParameter(T0_PAR_NR_BACK, scramm_t0_fermi->GetParameter(T0_PAR_NR_BACK));
351 m_t0_fermi->SetParameter(T0_PAR_NR_A, scramm_t0_fermi->GetParameter(T0_PAR_NR_A));
352 // set estimates as start values
353 m_t0_fermi->SetParameter(T0_PAR_NR_T0, scramm_t0_fermi->GetParameter(T0_PAR_NR_T0));
354 // set resonable start value for T
355 m_t0_fermi->SetParameter(T0_PAR_NR_T, scramm_t0_fermi->GetParameter(T0_PAR_NR_T));
356 m_time->GetListOfFunctions()->Clear();
357 m_time->Fit("mt_t0_fermi", fitopt.c_str(), "", scramrec.GetFitRangeMin(), scramrec.GetFitRangeMax());
358 if (m_settings->T0Settings()->UseTopChi2())
359 TopChi2();
360 else
361 m_chi2 = m_time->GetFunction("mt_t0_fermi")->GetChisquare() / m_time->GetFunction("mt_t0_fermi")->GetNDF();
362 cwd->cd();
363 return true;
364 } // end T0MTHistos::T0Scramble
365
367 // calculate topchi2
368 m_chi2 = 0;
369 int topndf = 0;
370 TF1 *t0_fermi = m_time->GetFunction("mt_t0_fermi");
371 Double_t min, max;
372 t0_fermi->GetRange(min, max);
373 int startbin = m_time->FindBin(min) + 1;
374 int endbin = m_time->FindBin(max) - 1;
375 for (int bin = startbin; bin < endbin; bin++) {
376 float measval = m_time->GetBinContent(bin);
377 float funcval = t0_fermi->Eval(m_time->GetBinCenter(bin));
378 float errval = m_time->GetBinError(bin);
379 // take only chi2 from top part or if the bin content is min 10 or if the function >10
380 if (measval < 10 && funcval < 10 &&
381 m_time->GetBinCenter(bin) < t0_fermi->GetParameter(T0_PAR_NR_T0) + 2 * t0_fermi->GetParameter(T0_PAR_NR_T))
382 continue;
383 if (errval == 0)
384 m_chi2 += (measval - funcval) * (measval - funcval);
385 else
386 m_chi2 += (measval - funcval) * (measval - funcval) / (errval * errval);
387 topndf++;
388 }
389 if (topndf != 0) m_chi2 = m_chi2 / topndf;
390 } // end T0MTHistos::TopChi2
391
393#ifndef NDEBUG
394 MsgStream log(Athena::getMessageSvc(), "T0MTHistos");
395 if (log.level() <= MSG::VERBOSE) log << MSG::VERBOSE << "Slicing for " << m_time->GetName() << endmsg;
396#endif
397 // vector with slice chi2
398 std::list<Slice> slice_chi2;
399 Slice current;
400 current.chi_2 = 0.0;
401 current.n_bins = 0;
402 TF1 *t0_fermi = m_time->GetFunction("mt_t0_fermi");
403 current.min_bin = m_time->FindBin(t0_fermi->GetParameter(T0_PAR_NR_T0) + 2 * t0_fermi->GetParameter(T0_PAR_NR_T));
404 Double_t min, max;
405 t0_fermi->GetRange(min, max);
406#ifndef NDEBUG
407 if (log.level() <= MSG::VERBOSE) log << MSG::VERBOSE << current.min_bin << " " << max << endmsg;
408#endif
409 for (int bin = m_time->FindBin(t0_fermi->GetParameter(T0_PAR_NR_T0) + 2 * t0_fermi->GetParameter(T0_PAR_NR_T));
410 bin < m_time->FindBin(max) - 1; bin++) {
411 if (current.n_bins == 10) {
412#ifndef NDEBUG
413 if (log.level() <= MSG::VERBOSE) log << MSG::VERBOSE << current.chi_2 / current.n_bins << endmsg;
414#endif
415 current.max_bin = bin;
416 slice_chi2.push_back(current);
417 current.chi_2 = 0.0;
418 current.min_bin = bin;
419 current.n_bins = 0;
420 }
421 double measval = m_time->GetBinContent(bin);
422 double funcval = t0_fermi->Eval(m_time->GetBinCenter(bin));
423 double errval = m_time->GetBinError(bin);
424 if (errval == 0)
425 current.chi_2 += std::pow(measval - funcval, 2.0);
426 else
427 current.chi_2 += std::pow(measval - funcval, 2.0) / std::pow(errval, 2);
428 current.n_bins++;
429 }
430#ifndef NDEBUG
431 if (log.level() <= MSG::VERBOSE) log << MSG::VERBOSE << "number of slices: " << slice_chi2.size() << endmsg;
432#endif
433 std::list<Slice>::iterator it = slice_chi2.end();
434 do {
435 --it;
436 if (it == slice_chi2.begin()) {
437#ifndef NDEBUG
438 if (log.level() <= MSG::VERBOSE) log << MSG::VERBOSE << "No gain in slicing!" << endmsg;
439#endif
440 return;
441 }
442 } while (it->chi_2 / static_cast<double>(it->n_bins) > 3);
443 max = m_time->GetBinCenter(it->min_bin);
444 m_time->GetListOfFunctions()->Clear();
445 std::string fitopt("BLR");
446 if (m_settings->VerboseLevel() == 0) fitopt += "Q";
447 if (m_settings->AddFitfun()) {
448 fitopt += "+";
449 } else {
450 fitopt += "0";
451 }
452 m_time->Fit("mt_t0_fermi", fitopt.c_str(), "", min, max);
453 if (m_settings->T0Settings()->UseTopChi2())
454 TopChi2();
455 else
456 m_chi2 = m_time->GetFunction("mt_t0_fermi")->GetChisquare() / m_time->GetFunction("mt_t0_fermi")->GetNDF();
457 } // end T0MTHistos::TopSlicing
458
459} // namespace MuonCalib
#define endmsg
static Double_t a
static TRandom * rnd
#define x
#define min(a, b)
Definition cfImp.cxx:40
#define max(a, b)
Definition cfImp.cxx:41
Define macros for attributes used to control the static checker.
#define ATLAS_THREAD_SAFE
double GetEstimatedT0() const
get estimated t0
bool Initialize(TH1F *hist, const T0MTSettings *settings)
Initialize class - returns true if pattern recognition was successfull.
double GetBackground() const
get the background level
double GetFitRangeMax() const
get fit range
double GetFitRangeMin() const
get fit range
double GetB() const
get parameter a in exp-function representing the end of the spectrum
double GetEstimatedTMax() const
get estimated t0
bool Initialize(TH1F *hist, double t0, const T0MTSettings *settings)
Initialize class.
double GetBackground() const
get the background level
double GetA() const
get parameter a in exp-function representing the end of the spectrum
static constexpr int T0_PAR_NR_T0
parameter numbers in t0 fit
Definition T0MTHistos.h:57
static constexpr int N_T0_FIT_PAR
number of parameters in t0 fit
Definition T0MTHistos.h:55
static constexpr int TMAX_PAR_NR_TMAX
parameters numbers for tmax fit
Definition T0MTHistos.h:61
static constexpr int TMAX_PAR_NR_T0
Definition T0MTHistos.h:62
bool m_t0_ok
is true if t0 fit was successful
Definition T0MTHistos.h:125
static constexpr int T0_PAR_NR_T
Definition T0MTHistos.h:57
double m_tmax_ok
is true if tmax fit was successful
Definition T0MTHistos.h:131
static constexpr int TMAX_PAR_NR_T
Definition T0MTHistos.h:61
void Initialize(int id, const T0MTSettings *settings, const char *hname=nullptr)
Initialize class.
std::unique_ptr< TF1 > m_t0_fermi
function fitted to the riding edghe of the spectrum
Definition T0MTHistos.h:123
static constexpr int TMAX_PAR_NR_B
Definition T0MTHistos.h:61
int m_status_code
status code for t0 fit (0 ok, 1 not fitted, 2 low statistics, 3 failed)
Definition T0MTHistos.h:127
bool T0Scramble()
try to get better start values from a scrambled histogram
static constexpr int TMAX_PAR_NR_BACK
Definition T0MTHistos.h:61
bool FitT0()
Perform t0-fit Returns true if fit is successfull.
bool FitTmax()
Performs tmax-fit Returns true if fit is successfull.
void SetTSpec(int id, TH1F *spec, const T0MTSettings *settings, bool copy_spec=true)
set the pointer of the drift-time spectrum to an existing spectrum.
static constexpr int T0_PAR_NR_A
Definition T0MTHistos.h:57
static constexpr int N_TMAX_FIT_PAR
number of parameters for tmax fit
Definition T0MTHistos.h:59
double m_chi2
chi2/NDF value
Definition T0MTHistos.h:137
void TopChi2()
top chi2 calculation
bool NormalFit()
normal t0 fit
TDirectory * m_dir
TDirectory where debug and result histograms are stored.
Definition T0MTHistos.h:133
std::unique_ptr< TH1F > m_time
time spectrum
Definition T0MTHistos.h:119
std::unique_ptr< TF1 > m_tmax_fermi
function fitted to the falling edge of the spectrum
Definition T0MTHistos.h:129
const T0MTSettings * m_settings
Pointer to settings class.
Definition T0MTHistos.h:135
static constexpr int TMAX_PAR_NR_A
Definition T0MTHistos.h:61
void TopSlicing()
top slicing metyhod
void FillT(double t)
fill drift time spectrum
static constexpr int T0_PAR_NR_BACK
Definition T0MTHistos.h:57
Settings for the T0 calibration (histogram booking and fitting parameters) Parameters for pattern rec...
const bool & DrawDebugGraphs() const
If set to true for every tube a TDirectory will be created.
const int & NBinsTime() const
Number of bins for time histogram and range.
double TimeMax() const
double TimeMin() const
singleton-like access to IMessageSvc via open function and helper
std::string cwd
Definition listroot.cxx:38
IMessageSvc * getMessageSvc(bool quiet=false)
CscCalcPed - algorithm that finds the Cathode Strip Chamber pedestals from an RDO.
TRandom3 * getTLSRandomGen()
Double_t mt_t0_fermi(Double_t *x, Double_t *par)
The fermi function to be fitted at the rising edge of the spectrum.
Double_t mt_tmax_fermi(Double_t *x, Double_t *par)
The fermi function to be fitted at the trailing slope of the spectrum.
hold the test vectors and ease the comparison