ATLAS Offline Software
ZDCPulseAnalyzer.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3 */
4 
6 
7 #include "TFitResult.h"
8 #include "TFitResultPtr.h"
9 #include "TVirtualFitter.h"
10 #include "TList.h"
11 #include "TMinuit.h"
12 
13 #include <algorithm>
14 #include <sstream>
15 #include <cmath>
16 #include <numeric>
17 #include <iomanip>
18 #include <stdexcept>
19 
21 //
22 // List of allowed JSON configuration parameters
23 //
24 // For each parameter we have name, JSON value type, whether it can be set per channel, and whether it is required
25 //
26 // if the type is -1, then there's no value, the presence of the parameter itself is a boolean -- i.e. enabling
27 //
29  {"tag", {JSON::value_t::string, 1, true, true}},
30  {"enabled", {JSON::value_t::boolean, 1, true, true}},
31  {"LGMode", {JSON::value_t::number_unsigned, 1, false, true}},
32  {"Nsample", {JSON::value_t::number_integer, 1, false, true}},
33  {"FADCFreqMHz", {JSON::value_t::number_integer, 1, false, true}},
34  {"preSampleIdx", {JSON::value_t::number_integer, 1, true, true}},
35  {"nominalPedestal", {JSON::value_t::number_integer, 1, false, true}},
36  {"fitFunction", {JSON::value_t::string, 1, true, true}},
37  {"peakSample", {JSON::value_t::number_integer, 1, true, true}},
38  {"peakTolerance", {JSON::value_t::number_integer, 1, true, false}},
39  {"2ndDerivThreshHG", {JSON::value_t::number_integer, 1, true, true}},
40  {"2ndDerivThreshLG", {JSON::value_t::number_integer, 1, true, true}},
41  {"2ndDerivStep", {JSON::value_t::number_integer, 1, false, false}},
42  {"HGOverflowADC", {JSON::value_t::number_integer, 1, true, true}},
43  {"HGUnderflowADC", {JSON::value_t::number_integer, 1, true, true}},
44  {"LGOverflowADC", {JSON::value_t::number_integer, 1, true, true}},
45  {"nominalT0HG", {JSON::value_t::number_float, 1, true, true}},
46  {"nominalT0LG", {JSON::value_t::number_float, 1, true, true}},
47  {"nominalTau1", {JSON::value_t::number_float, 1, true, false}},
48  {"nominalTau2", {JSON::value_t::number_float, 1, true, false}},
49  {"fixTau1", {JSON::value_t::boolean, 1, true, false}},
50  {"fixTau2", {JSON::value_t::boolean, 1, true, false}},
51  {"T0CutsHG", {JSON::value_t::array, 2, true, true}},
52  {"T0CutsLG", {JSON::value_t::array, 2, true, true}},
53  {"chisqDivAmpCutHG", {JSON::value_t::number_float, 1, true, true}},
54  {"chisqDivAmpCutLG", {JSON::value_t::number_float, 1, true, true}},
55  {"chisqDivAmpOffsetHG", {JSON::value_t::number_float, 1, true, true}},
56  {"chisqDivAmpOffsetLG", {JSON::value_t::number_float, 1, true, true}},
57  {"chisqDivAmpPowerHG", {JSON::value_t::number_float, 1, true, true}},
58  {"chisqDivAmpPowerLG", {JSON::value_t::number_float, 1, true, true}},
59  {"gainFactorHG", {JSON::value_t::number_float, 1, true, true}},
60  {"gainFactorLG", {JSON::value_t::number_float, 1, true, true}},
61  {"noiseSigmaHG", {JSON::value_t::number_float, 1, true, true}},
62  {"noiseSigmaLG", {JSON::value_t::number_float, 1, true, true}},
63  {"enableRepass", {JSON::value_t::boolean, 1, false, false}},
64  {"Repass2ndDerivThreshHG", {JSON::value_t::number_integer, 1, true, true}},
65  {"Repass2ndDerivThreshLG", {JSON::value_t::number_integer, 1, true, true}},
66  {"fitAmpMinMaxHG", {JSON::value_t::array, 2, true, false}},
67  {"fitAmpMinMaxLG", {JSON::value_t::array, 2, true, false}},
68  {"ampMinSignifHGLG", {JSON::value_t::array, 2, true, false}},
69  {"enablePreExclusion", {JSON::value_t::array, 3, false, false}},
70  {"enablePostExclusion", {JSON::value_t::array, 3, false, false}},
71  {"enableTimingCorrection", {JSON::value_t::array, 2, false, false}},
72  {"timeCorrCoeffHG", {JSON::value_t::array, 0, true, false}},
73  {"timeCorrCoeffLG", {JSON::value_t::array, 0, true, false}},
74  {"enableADCNLCorrection", {JSON::value_t::array, 3, false, false}},
75  {"ADCNLCorrCoeffs", {JSON::value_t::array, 0, true, false}},
76  {"useDelayed", {JSON::value_t::boolean, 1, false, false}},
77  {"delayDeltaT", {JSON::value_t::number_float, 1, true, false}}
78 };
79 
80 
85 float ZDCPulseAnalyzer::s_combinedFitTMin = -0.5; // add to allow switch to high gain by skipping early samples
86 std::vector<float> ZDCPulseAnalyzer::s_pullValues;
87 
88 void ZDCPulseAnalyzer::CombinedPulsesFCN(int& /*numParam*/, double*, double& f, double* par, int flag)
89 {
90  // The first parameter is a correction factor to account for decrease in beam intensity between x
91  // and y scan. It is applied here and not passed to the actual fit function
92  //
93  int nSamples = s_undelayedFitHist->GetNbinsX();
94 
95  if (flag == 3) {
96  s_pullValues.assign(nSamples * 2, 0);
97  }
98 
99  double chiSquare = 0;
100 
101  float delayBaselineAdjust = par[0];
102 
103  // undelayed
104  //
105  for (int isample = 0; isample < nSamples; isample++) {
106  double histValue = s_undelayedFitHist->GetBinContent(isample + 1);
107  double histError = std::max(s_undelayedFitHist->GetBinError(isample + 1), 1.0);
108  double t = s_undelayedFitHist->GetBinCenter(isample + 1);
109 
110  if (t > s_combinedFitTMax) break;
111  if (t < s_combinedFitTMin) continue;
112 
113  double funcVal = s_combinedFitFunc->EvalPar(&t, &par[1]);
114 
115  double pull = (histValue - funcVal) / histError;
116 
117  if (flag == 3) s_pullValues[2 * isample] = pull;
118  chiSquare += pull * pull;
119  }
120 
121  // delayed
122  //
123  for (int isample = 0; isample < nSamples; isample++) {
124  double histValue = s_delayedFitHist->GetBinContent(isample + 1);
125  double histError = std::max(s_delayedFitHist->GetBinError(isample + 1), 1.0);
126  double t = s_delayedFitHist->GetBinCenter(isample + 1);
127 
128  if (t > s_combinedFitTMax) break;
129  if (t < s_combinedFitTMin) continue;
130 
131  double funcVal = s_combinedFitFunc->EvalPar(&t, &par[1]) + delayBaselineAdjust;
132  double pull = (histValue - funcVal) / histError;
133 
134  if (flag == 3) s_pullValues[2 * isample + 1] = pull;
135  chiSquare += pull * pull;
136  }
137 
138  f = chiSquare;
139 }
140 
141 
142 ZDCPulseAnalyzer::ZDCPulseAnalyzer(ZDCMsg::MessageFunctionPtr msgFunc_p, const std::string& tag, int Nsample, float deltaTSample, size_t preSampleIdx, int pedestal,
143  const std::string& fitFunction, int peak2ndDerivMinSample,
144  float peak2ndDerivMinThreshHG, float peak2ndDerivMinThreshLG) :
145  m_msgFunc_p(std::move(msgFunc_p)),
146  m_tag(tag), m_Nsample(Nsample),
147  m_preSampleIdx(preSampleIdx),
148  m_deltaTSample(deltaTSample),
149  m_pedestal(pedestal), m_fitFunction(fitFunction),
150  m_peak2ndDerivMinSample(peak2ndDerivMinSample),
151  m_peak2ndDerivMinThreshLG(peak2ndDerivMinThreshLG),
152  m_peak2ndDerivMinThreshHG(peak2ndDerivMinThreshHG),
153  m_ADCSamplesHGSub(Nsample, 0), m_ADCSamplesLGSub(Nsample, 0),
154  m_ADCSSampSigHG(Nsample, 0), m_ADCSSampSigLG(Nsample, 0),
155  m_samplesSub(Nsample, 0)
156 {
157  // Create the histogram used for fitting
158  //
159  m_tmin = -deltaTSample / 2;
160  m_tmax = m_tmin + ((float) Nsample) * deltaTSample;
162 
163  std::string histName = "ZDCFitHist" + tag;
164  std::string histNameLGRefit = "ZDCFitHist" + tag + "_LGRefit";
165 
166  m_fitHist = std::make_unique<TH1F>(histName.c_str(), "", m_Nsample, m_tmin, m_tmax);
167  m_fitHistLGRefit = std::make_unique<TH1F>(histNameLGRefit.c_str(), "", m_Nsample, m_tmin, m_tmax);
168 
169  m_fitHist->SetDirectory(0);
170  m_fitHistLGRefit->SetDirectory(0);
171 
172  SetDefaults();
173  Reset();
174 }
175 
177  m_msgFunc_p(std::move(msgFunc_p))
178 {
179  SetDefaults();
180 
181  // auto [result, resultString] = ValidateJSONConfig(configJSON);
182  // (*m_msgFunc_p)(ZDCMsg::Debug, "ValidateJSON produced result: " + resultString);
183 
184  auto [result2, resultString2] = ConfigFromJSON(configJSON);
185  (*m_msgFunc_p)(ZDCMsg::Debug, "ConfigFromJSON produced result: "+ resultString2);
186 
188 
189  // Create the histogram used for fitting
190  //
191  if (m_freqMHz >1.e-6) m_deltaTSample = 1000./m_freqMHz;
192  m_tmin = -m_deltaTSample / 2;
195 
196  std::string histName = "ZDCFitHist" + m_tag;
197  std::string histNameLGRefit = "ZDCFitHist" + m_tag + "_LGRefit";
198 
199  m_fitHist = std::make_unique<TH1F>(histName.c_str(), "", m_Nsample, m_tmin, m_tmax);
200  m_fitHistLGRefit = std::make_unique<TH1F>(histNameLGRefit.c_str(), "", m_Nsample, m_tmin, m_tmax);
201 
202  m_fitHist->SetDirectory(0);
203  m_fitHistLGRefit->SetDirectory(0);
204 
205  Reset();
206 }
207 
208 void ZDCPulseAnalyzer::enableDelayed(float deltaT, float pedestalShift, bool fixedBaseline)
209 {
210  m_useDelayed = true;
211  m_useFixedBaseline = fixedBaseline;
212 
213  m_delayedDeltaT = deltaT;
214  m_delayedPedestalDiff = pedestalShift;
215 
216  m_deltaTSample /= 2.;
217 
218  std::string delayedHGName = std::string(m_fitHist->GetName()) + "delayed";
219  std::string delayedLGName = std::string(m_fitHistLGRefit->GetName()) + "delayed";
220 
221  m_delayedHist = std::make_unique<TH1F>(delayedHGName.c_str(), "", m_Nsample, m_tmin + m_delayedDeltaT, m_tmax + m_delayedDeltaT);
222  m_delayedHist->SetDirectory(0);
223 
224  m_delayedHistLGRefit = std::make_unique<TH1F>(delayedLGName.c_str(), "", m_Nsample, m_tmin + m_delayedDeltaT, m_tmax + m_delayedDeltaT);
225  m_delayedHistLGRefit->SetDirectory(0);
226 
227  m_ADCSamplesHGSub.assign(2 * m_Nsample, 0);
228  m_ADCSamplesLGSub.assign(2 * m_Nsample, 0);
229 }
230 
231 void ZDCPulseAnalyzer::enableRepass(float peak2ndDerivMinRepassHG, float peak2ndDerivMinRepassLG)
232 {
233  m_enableRepass = true;
234  m_peak2ndDerivMinRepassHG = peak2ndDerivMinRepassHG;
235  m_peak2ndDerivMinRepassLG = peak2ndDerivMinRepassLG;
236 }
237 
239 {
241 
242  m_nominalTau1 = 1.5;
243  m_nominalTau2 = 5;
244 
245  m_fixTau1 = false;
246  m_fixTau2 = false;
247 
248  m_HGOverflowADC = 3500;
249  m_HGUnderflowADC = 1;
250  m_LGOverflowADC = 3900;
251 
252  // Default values for the gain factors uswed to match low and high gain
253  //
254  m_gainFactorLG = 10;
255  m_gainFactorHG = 1;
256 
257  m_2ndDerivStep = 1;
258 
259  m_noiseSigHG = 1;
260  m_noiseSigLG = 1;
261 
262  m_sigMinLG = 0;
263  m_sigMinHG = 0;
264 
265  m_timeCutMode = 0;
266  m_chisqDivAmpCutLG = 100;
267  m_chisqDivAmpCutHG = 100;
268 
270  m_chisqDivAmpOffsetHG = 150;
271 
272  m_chisqDivAmpPowerLG = 1.2;
273  m_chisqDivAmpPowerHG = 1.2;
274 
275  m_LGT0CorrParams.assign(4, 0);
276  m_HGT0CorrParams.assign(4, 0);
277 
278  // m_defaultFitTMax = m_tmax;
279  // m_defaultFitTMin = m_tmin;
280 
281  m_fitAmpMinHG = 1;
282  m_fitAmpMinLG = 1;
283 
284  m_fitAmpMaxHG = 1500;
285  m_fitAmpMaxLG = 1500;
286 
287  m_haveSignifCuts = false;
288 
289  m_postPulse = false;
290  m_prePulse = false;
291 
292  m_initialPrePulseT0 = -10;
294 
295  m_initialPostPulseT0 = 100;
296 
297  m_initialExpAmp = 0;
298  m_fitPostT0lo = 0;
299 
300  m_useDelayed = false;
301  m_enablePreExcl = false;
302  m_enablePostExcl = false;
303 
305  m_haveNonlinCorr = false;
306  m_quietFits = true;
307  m_saveFitFunc = false;
308  m_fitOptions = "s";
309 }
310 
311 void ZDCPulseAnalyzer::Reset(bool repass)
312 {
313  if (!repass) {
314  m_haveData = false;
315 
316  m_useLowGain = false;
317  m_fail = false;
318  m_HGOverflow = false;
319 
320  m_HGUnderflow = false;
321  m_PSHGOverUnderflow = false;
322  m_LGOverflow = false;
323  m_LGUnderflow = false;
324 
325  m_ExcludeEarly = false;
326  m_ExcludeLate = false;
327 
328  m_adjTimeRangeEvent = false;
329  m_backToHG_pre = false;
330  m_fixPrePulse = false;
331 
332  m_minADCLG = -1;
333  m_maxADCLG = -1;
334  m_minADCHG = -1;
335  m_maxADCHG = -1;
336 
337  m_minADCSampleHG = -1;
338  m_maxADCSampleHG = -1;
339  m_minADCSampleLG = -1;
340  m_maxADCSampleLG = -1;
341 
342  m_ADCPeakHG = -1;
343  m_ADCPeakLG = -1;
344 
345  int sampleVecSize = m_Nsample;
346  if (m_useDelayed) sampleVecSize *= 2;
347 
348  m_ADCSamplesHG.clear();
349  m_ADCSamplesLG.clear();
350 
351  m_ADCSamplesHGSub.clear();
352  m_ADCSamplesLGSub.clear();
353 
354  m_ADCSSampSigHG.assign(sampleVecSize, m_noiseSigHG);
355  m_ADCSSampSigLG.assign(sampleVecSize, m_noiseSigLG);
356 
357  m_minSampleEvt = 0;
358  m_maxSampleEvt = (m_useDelayed ? 2 * m_Nsample - 1 : m_Nsample - 1);
359 
360  m_usedPresampIdx = 0;
361 
364 
365  m_lastHGOverFlowSample = -999;
367  }
368 
369 
372 
373  if (m_initializedFits) {
377  }
378 
379  // -----------------------
380  // Statuses
381  //
382  m_havePulse = false;
383 
384  m_prePulse = false;
385  m_postPulse = false;
386  m_fitFailed = false;
387  m_badChisq = false;
388 
389  m_badT0 = false;
390  m_preExpTail = false;
391  m_repassPulse = false;
392 
393  m_fitMinAmp = false;
394  m_evtLGRefit = false;
395  m_failSigCut = false;
396 
397  // -----------------------
398 
400 
401  m_fitAmplitude = 0;
402  m_ampNoNonLin = 0;
403  m_fitTime = -100;
404  m_fitTimeSub = -100;
405  m_fitTimeCorr = -100;
406  m_fitTCorr2nd = -100;
407 
408  m_fitPreT0 = -100;
409  m_fitPreAmp = -100;
410  m_fitPostT0 = -100;
411  m_fitPostAmp = -100;
412  m_fitExpAmp = -100;
413 
414  m_minDeriv2ndSig = -10;
415  m_preExpSig = -10;
416  m_prePulseSig = -10;
417 
418  m_fitChisq = 0;
419 
420  m_amplitude = 0;
421  m_ampError = 0;
422  m_preSampleAmp = 0;
423  m_preAmplitude = 0;
424  m_postAmplitude = 0;
425  m_expAmplitude = 0;
426  m_bkgdMaxFraction = 0;
427 
428  m_refitLGAmpl = 0;
429  m_refitLGAmpError = 0;
430  m_refitLGChisq = 0;
431  m_refitLGTime = -100;
432  m_refitLGTimeSub = -100;
433 
434  m_initialPrePulseT0 = -10;
436 
437  m_initialPostPulseT0 = 100;
438 
439  m_initialExpAmp = 0;
440  m_fitPostT0lo = 0;
441 
442  m_fitPulls.clear();
443 
444  m_samplesSub.clear();
445  m_samplesDeriv2nd.clear();
446 }
447 
448 
449 void ZDCPulseAnalyzer::setMinimumSignificance(float sigMinHG, float sigMinLG)
450 {
451  m_haveSignifCuts = true;
452  m_sigMinHG = sigMinHG;
453  m_sigMinLG = sigMinLG;
454 }
455 
456 void ZDCPulseAnalyzer::SetGainFactorsHGLG(float gainFactorHG, float gainFactorLG)
457 {
458  m_gainFactorHG = gainFactorHG;
459  m_gainFactorLG = gainFactorLG;
460 }
461 
462 void ZDCPulseAnalyzer::SetFitMinMaxAmp(float minAmpHG, float minAmpLG, float maxAmpHG, float maxAmpLG)
463 {
464  m_fitAmpMinHG = minAmpHG;
465  m_fitAmpMinLG = minAmpLG;
466 
467  m_fitAmpMaxHG = maxAmpHG;
468  m_fitAmpMaxLG = maxAmpLG;
469 }
470 
471 void ZDCPulseAnalyzer::SetTauT0Values(bool fixTau1, bool fixTau2, float tau1, float tau2, float t0HG, float t0LG)
472 {
473  m_fixTau1 = fixTau1;
474  m_fixTau2 = fixTau2;
475  m_nominalTau1 = tau1;
476  m_nominalTau2 = tau2;
477 
478  m_nominalT0HG = t0HG;
479  m_nominalT0LG = t0LG;
480 
481  std::ostringstream ostrm;
482  ostrm << "ZDCPulseAnalyzer::SetTauT0Values:: m_fixTau1=" << m_fixTau1 << " m_fixTau2=" << m_fixTau2 << " m_nominalTau1=" << m_nominalTau1 << " m_nominalTau2=" << m_nominalTau2 << " m_nominalT0HG=" << m_nominalT0HG << " m_nominalT0LG=" << m_nominalT0LG;
483 
484  (*m_msgFunc_p)(ZDCMsg::Info, ostrm.str());
485 
486  m_initializedFits = false;
487 }
488 
490 {
491  if (tmax < m_tmin) {
492  (*m_msgFunc_p)(ZDCMsg::Error, ("ZDCPulseAnalyzer::SetFitTimeMax:: invalid FitTimeMax: " + std::to_string(tmax)));
493  return;
494  }
495 
497 
499 }
500 
501 void ZDCPulseAnalyzer::SetADCOverUnderflowValues(int HGOverflowADC, int HGUnderflowADC, int LGOverflowADC)
502 {
503  m_HGOverflowADC = HGOverflowADC;
504  m_LGOverflowADC = LGOverflowADC;
505  m_HGUnderflowADC = HGUnderflowADC;
506 }
507 
508 void ZDCPulseAnalyzer::SetCutValues(float chisqDivAmpCutHG, float chisqDivAmpCutLG,
509  float deltaT0MinHG, float deltaT0MaxHG,
510  float deltaT0MinLG, float deltaT0MaxLG)
511 {
512  m_chisqDivAmpCutHG = chisqDivAmpCutHG;
513  m_chisqDivAmpCutLG = chisqDivAmpCutLG;
514 
515  m_T0CutLowHG = deltaT0MinHG;
516  m_T0CutLowLG = deltaT0MinLG;
517 
518  m_T0CutHighHG = deltaT0MaxHG;
519  m_T0CutHighLG = deltaT0MaxLG;
520 }
521 
522 void ZDCPulseAnalyzer::enableTimeSigCut(bool AND, float sigCut, const std::string& TF1String,
523  const std::vector<double>& parsHG,
524  const std::vector<double>& parsLG)
525 {
526  m_timeCutMode = AND ? 2 : 1;
527  m_t0CutSig = sigCut;
528 
529  // Make the TF1's that provide the resolution
530  //
531  std::string timeResHGName = "TimeResFuncHG_" + m_tag;
532  std::string timeResLGName = "TimeResFuncLG_" + m_tag;
533 
534  TF1* funcHG_p = new TF1(timeResHGName.c_str(), TF1String.c_str(), 0, m_HGOverflowADC);
535  TF1* funcLG_p = new TF1(timeResLGName.c_str(), TF1String.c_str(), 0, m_LGOverflowADC);
536 
537  if (parsHG.size() != static_cast<unsigned int>(funcHG_p->GetNpar()) ||
538  parsLG.size() != static_cast<unsigned int>(funcLG_p->GetNpar())) {
539  //
540  // generate an error message
541  }
542 
543  funcHG_p->SetParameters(&parsHG[0]);
544  funcLG_p->SetParameters(&parsLG[0]);
545 
546  m_timeResFuncHG_p.reset(funcHG_p);
547  m_timeResFuncLG_p.reset(funcLG_p);
548 
549 }
550 
551 void ZDCPulseAnalyzer::enableFADCCorrections(bool correctPerSample, std::unique_ptr<const TH1>& correHistHG, std::unique_ptr<const TH1>& correHistLG)
552 {
553  m_haveFADCCorrections = true;
554  m_FADCCorrPerSample = correctPerSample;
555 
556  // check for appropriate limits
557  //
558  auto getXmin=[](const TH1 * pH){
559  return pH->GetXaxis()->GetXmin();
560  };
561  auto getXmax=[](const TH1 * pH){
562  return pH->GetXaxis()->GetXmax();
563  };
564  auto xmin= getXmin(correHistHG.get());
565  auto xmax= getXmax(correHistHG.get());
566  if (std::abs(xmin+0.5) > 1e-3 || std::abs(xmax - 4095.5) > 1e-3) {
567  (*m_msgFunc_p)(ZDCMsg::Error, "ZDCPulseAnalyzer::enableFADCCorrections:: invalid high gain correction histogram range: xmin, xmax = " +
568  std::to_string(xmin ) + ", " + std::to_string(xmax) );
569  }
570  else {
571  m_FADCCorrHG = std::move(correHistHG);
572  }
573  xmin= getXmin(correHistLG.get());
574  xmax= getXmax(correHistLG.get());
575  if (std::abs(xmin+0.5) > 1e-3 ||
576  std::abs(xmax - 4095.5) > 1e-3) {
577  (*m_msgFunc_p)(ZDCMsg::Error, "ZDCPulseAnalyzer::enableFADCCorrections:: invalid low gain correction histogram range: xmin, xmax = " +
578  std::to_string(xmin) + ", " + std::to_string(xmax) );
579  }
580  else {
581  m_FADCCorrLG = std::move(correHistLG);
582  }
583 }
584 
585 std::vector<float> ZDCPulseAnalyzer::GetFitPulls(bool refitLG) const
586 {
587  //
588  // If there was no pulse for this event, return an empty vector (see reset() method)
589  //
590  if (!m_havePulse) {
591  return m_fitPulls;
592  }
593 
594  //
595  // The combined (delayed+undelayed) pulse fitting fills out m_fitPulls directly
596  //
597  if (m_useDelayed) {
598  return m_fitPulls;
599  }
600  else {
601  //
602  // When not using combined fitting, We don't have the pre-calculated pulls. Calculate them on the fly.
603  //
604  std::vector<float> pulls(m_Nsample, -100);
605 
606  const TH1* dataHist_p = refitLG ? m_fitHistLGRefit.get() : m_fitHist.get();
607  const TF1* fit_p = static_cast<const TF1*>(dataHist_p->GetListOfFunctions()->Last());
608 
609  for (size_t ibin = 0; ibin < m_Nsample ; ibin++) {
610  float t = dataHist_p->GetBinCenter(ibin + 1);
611  float fitVal = fit_p->Eval(t);
612  float histVal = dataHist_p->GetBinContent(ibin + 1);
613  float histErr = dataHist_p->GetBinError(ibin + 1);
614  float pull = (histVal - fitVal)/histErr;
615  pulls[ibin] = pull;
616  }
617 
618  return pulls;
619  }
620 }
621 
623 {
624  double amplCorrFactor = 1;
625 
626  // If we have FADC correction and we aren't applying it per-sample, do so here
627  //
629  double fadcCorr = highGain ? m_FADCCorrHG->Interpolate(m_fitAmplitude) : m_FADCCorrLG->Interpolate(m_fitAmplitude);
630  amplCorrFactor *= fadcCorr;
631  }
632 
633  double amplCorr = m_fitAmplitude * amplCorrFactor;
634 
635  // If we have a non-linear correction, apply it here
636  // We apply it as an inverse correction - i.e. we divide by a correction
637  // term tha is a sum of coefficients times the ADC minus a reference
638  // to a power. The lowest power is 1, the highest is deteremined by
639  // the number of provided coefficients
640 
641  if (m_haveNonlinCorr) {
642  float invNLCorr = 1.0;
643  float nlPolyArg = (amplCorr - m_nonLinCorrRefADC) / m_nonLinCorrRefScale;
644 
645  if (highGain) {
646  for (size_t power = 1; power <= m_nonLinCorrParamsHG.size(); power++) {
647  invNLCorr += m_nonLinCorrParamsHG[power - 1]*pow(nlPolyArg, power);
648  }
649  }
650  else {
651  for (size_t power = 1; power <= m_nonLinCorrParamsLG.size(); power++) {
652  invNLCorr += m_nonLinCorrParamsLG[power - 1]*pow(nlPolyArg, power);
653  }
654  }
655 
656  amplCorrFactor /= invNLCorr;
657  }
658 
659  return amplCorrFactor;
660 }
661 
663 {
664  float prePulseTMin = 0;
665  float prePulseTMax = prePulseTMin + m_deltaTSample * (m_peak2ndDerivMinSample - m_peak2ndDerivMinTolerance);
666 
667  if (m_fitFunction == "FermiExp") {
668  if (!m_fixTau1 || !m_fixTau2) {
669  //
670  // Use the variable tau version of the expFermiFit
671  //
673  }
674  else {
675  m_defaultFitWrapper = std::unique_ptr<ZDCFitWrapper>(new ZDCFitExpFermiFixedTaus(m_tag, m_tmin, m_tmax, m_nominalTau1, m_nominalTau2));
676  }
677 
678  m_preExpFitWrapper = std::unique_ptr<ZDCFitExpFermiPreExp>(new ZDCFitExpFermiPreExp(m_tag, m_tmin, m_tmax, m_nominalTau1, m_nominalTau2, m_nominalTau2, false));
679  m_prePulseFitWrapper = std::unique_ptr<ZDCPrePulseFitWrapper>(new ZDCFitExpFermiPrePulse(m_tag, m_tmin, m_tmax, m_nominalTau1, m_nominalTau2));
680  }
681  else if (m_fitFunction == "FermiExpRun3") {
682  if (!m_fixTau1 || !m_fixTau2) {
683  //
684  // Use the variable tau version of the expFermiFit
685  //
687  }
688  else {
689  m_defaultFitWrapper = std::unique_ptr<ZDCFitWrapper>(new ZDCFitExpFermiFixedTaus(m_tag, m_tmin, m_tmax, m_nominalTau1, m_nominalTau2));
690  }
691 
692  m_preExpFitWrapper = std::unique_ptr<ZDCFitExpFermiPreExp>(new ZDCFitExpFermiPreExp(m_tag, m_tmin, m_tmax, m_nominalTau1, m_nominalTau2, 6, false));
693  m_prePulseFitWrapper = std::unique_ptr<ZDCPrePulseFitWrapper>(new ZDCFitExpFermiPrePulse(m_tag, m_tmin, m_tmax, m_nominalTau1, m_nominalTau2));
694  }
695  else if (m_fitFunction == "FermiExpLHCf") {
696  //
697  // Use the variable tau version of the expFermiFit
698  //
701 
702  m_preExpFitWrapper = std::unique_ptr<ZDCFitExpFermiLHCfPreExp>(new ZDCFitExpFermiLHCfPreExp(m_tag, m_tmin, m_tmax, m_nominalTau1, m_nominalTau2, 6, false));
703 
704  m_prePulseFitWrapper = std::unique_ptr<ZDCPrePulseFitWrapper>(new ZDCFitExpFermiLHCfPrePulse(m_tag, m_tmin, m_tmax, m_nominalTau1, m_nominalTau2));
705  }
706  else if (m_fitFunction == "FermiExpLinear") {
707  if (!m_fixTau1 || !m_fixTau2) {
708  //
709  // Use the variable tau version of the expFermiFit
710  //
712  }
713  else {
715  }
716 
717  m_preExpFitWrapper = std::unique_ptr<ZDCFitExpFermiPreExp>(new ZDCFitExpFermiPreExp(m_tag, m_tmin, m_tmax, m_nominalTau1, m_nominalTau2, 6, false));
718  m_prePulseFitWrapper = std::unique_ptr<ZDCPrePulseFitWrapper>(new ZDCFitExpFermiLinearPrePulse(m_tag, m_tmin, m_tmax, m_nominalTau1, m_nominalTau2));
719  }
720  else if (m_fitFunction == "ComplexPrePulse") {
721  if (!m_fixTau1 || !m_fixTau2) {
722  //
723  // Use the variable tau version of the expFermiFit
724  //
726  }
727  else {
729  }
730 
731  m_prePulseFitWrapper = std::unique_ptr<ZDCPrePulseFitWrapper>(new ZDCFitComplexPrePulse(m_tag, m_tmin, m_tmax, m_nominalTau1, m_nominalTau2));
732  }
733  else if (m_fitFunction == "GeneralPulse") {
734  if (!m_fixTau1 || !m_fixTau2) {
735  //
736  // Use the variable tau version of the expFermiFit
737  //
739  }
740  else {
742  }
743 
744  m_prePulseFitWrapper = std::unique_ptr<ZDCPrePulseFitWrapper>(new ZDCFitGeneralPulse(m_tag, m_tmin, m_tmax, m_nominalTau1, m_nominalTau2));
745  }
746  else {
747  (*m_msgFunc_p)(ZDCMsg::Fatal, "Wrong fit function type: " + m_fitFunction);
748  }
749 
750  m_prePulseFitWrapper->SetPrePulseT0Range(prePulseTMin, prePulseTMax);
751 
752  m_initializedFits = true;
753 }
754 
755 bool ZDCPulseAnalyzer::LoadAndAnalyzeData(const std::vector<float>& ADCSamplesHG, const std::vector<float>& ADCSamplesLG)
756 {
757  if (m_useDelayed) {
758  (*m_msgFunc_p)(ZDCMsg::Fatal, "ZDCPulseAnalyzer::LoadAndAnalyzeData:: Wrong LoadAndAnalyzeData called -- expecting both delayed and undelayed samples");
759  }
760 
762 
763  // Clear any transient data
764  //
765  Reset(false);
766 
767  // Make sure we have the right number of samples. Should never fail. necessry?
768  //
769  if (ADCSamplesHG.size() != m_Nsample || ADCSamplesLG.size() != m_Nsample) {
770  m_fail = true;
771  return false;
772  }
773 
774  m_ADCSamplesHG = ADCSamplesHG;
775  m_ADCSamplesLG = ADCSamplesLG;
776 
777  return DoAnalysis(false);
778 }
779 
780 bool ZDCPulseAnalyzer::LoadAndAnalyzeData(const std::vector<float>& ADCSamplesHG, const std::vector<float>& ADCSamplesLG,
781  const std::vector<float>& ADCSamplesHGDelayed, const std::vector<float>& ADCSamplesLGDelayed)
782 {
783  if (!m_useDelayed) {
784  (*m_msgFunc_p)(ZDCMsg::Fatal, "ZDCPulseAnalyzer::LoadAndAnalyzeData:: Wrong LoadAndAnalyzeData called -- expecting only undelayed samples");
785  }
786 
788 
789  // Clear any transient data
790  //
791  Reset();
792 
793  // Make sure we have the right number of samples. Should never fail. necessry?
794  //
795  if (ADCSamplesHG.size() != m_Nsample || ADCSamplesLG.size() != m_Nsample ||
796  ADCSamplesHGDelayed.size() != m_Nsample || ADCSamplesLGDelayed.size() != m_Nsample) {
797  m_fail = true;
798  return false;
799  }
800 
801  m_ADCSamplesHG.reserve(m_Nsample*2);
802  m_ADCSamplesLG.reserve(m_Nsample*2);
803 
804  // Now do pedestal subtraction and check for overflows
805  //
806  for (size_t isample = 0; isample < m_Nsample; isample++) {
807  float ADCHG = ADCSamplesHG[isample];
808  float ADCLG = ADCSamplesLG[isample];
809 
810  float ADCHGDelay = ADCSamplesHGDelayed[isample] - m_delayedPedestalDiff;
811  float ADCLGDelay = ADCSamplesLGDelayed[isample] - m_delayedPedestalDiff;
812 
813  if (m_delayedDeltaT > 0) {
814  m_ADCSamplesHG.push_back(ADCHG);
815  m_ADCSamplesHG.push_back(ADCHGDelay);
816 
817  m_ADCSamplesLG.push_back(ADCLG);
818  m_ADCSamplesLG.push_back(ADCLGDelay);
819  }
820  else {
821  m_ADCSamplesHG.push_back(ADCHGDelay);
822  m_ADCSamplesHG.push_back(ADCHG);
823 
824  m_ADCSamplesLG.push_back(ADCLGDelay);
825  m_ADCSamplesLG.push_back(ADCLG);
826  }
827  }
828 
829  return DoAnalysis(false);
830 }
831 
833 {
834  Reset(true);
835 
836  bool result = DoAnalysis(true);
837  if (result && havePulse()) {
838  m_repassPulse = true;
839  }
840 
841  return result;
842 }
843 
845 {
846  //
847  // Dump samples to verbose output
848  //
849  bool doDump = (*m_msgFunc_p)(ZDCMsg::Verbose, "Dumping all samples before subtraction: ");
850 
851  m_NSamplesAna = m_ADCSamplesHG.size();
852 
853  m_ADCSamplesHGSub.assign(m_NSamplesAna, 0);
854  m_ADCSamplesLGSub.assign(m_NSamplesAna, 0);
855 
856  m_useSampleHG.assign(m_NSamplesAna, true);
857  m_useSampleLG.assign(m_NSamplesAna, true);
858 
859  // Now do pedestal subtraction and check for overflows
860  //
863 
864  m_maxADCHG = 0;
865  m_maxADCLG = 0;
866 
867  for (size_t isample = 0; isample < m_NSamplesAna; isample++) {
868  float ADCHG = m_ADCSamplesHG[isample];
869  float ADCLG = m_ADCSamplesLG[isample];
870 
871  //
872  // We always pedestal subtract even if the sample isn't going to be used in analysis
873  // basically for diagnostics
874  //
875  m_ADCSamplesHGSub[isample] = ADCHG - m_pedestal;
876  m_ADCSamplesLGSub[isample] = ADCLG - m_pedestal;
877 
878 
879  // If we have the per-sample FADC correction, we apply it after pedestal correction
880  // since we analyze the FADC response after baseline (~ same as pedestal) subtraction
881  //
883  double fadcCorrHG = m_FADCCorrHG->Interpolate(m_ADCSamplesHGSub[isample]);
884  double fadcCorrLG = m_FADCCorrLG->Interpolate(m_ADCSamplesLGSub[isample]);
885 
886  m_ADCSamplesHGSub[isample] *= fadcCorrHG;
887  m_ADCSamplesLGSub[isample] *= fadcCorrLG;
888 
889  if (doDump) {
890  std::ostringstream dumpString;
891  dumpString << "After FADC correction, sample " << isample << ", HG ADC = " << m_ADCSamplesHGSub[isample]
892  << ", LG ADC = " << m_ADCSamplesLGSub[isample] << std::endl;
893  (*m_msgFunc_p)(ZDCMsg::Verbose, dumpString.str().c_str());
894  }
895  }
896 
897  if (ADCHG > m_maxADCHG) {
898  m_maxADCHG = ADCHG;
899  m_maxADCSampleHG = isample;
900  }
901  else if (ADCHG < m_minADCHG) {
902  m_minADCHG = ADCHG;
903  m_minADCSampleHG = isample;
904  }
905 
906  if (ADCLG > m_maxADCLG) {
907  m_maxADCLG = ADCLG;
908  m_maxADCSampleLG = isample;
909  }
910  else if (ADCLG < m_minADCLG) {
911  m_minADCLG = ADCLG;
912  m_minADCSampleLG = isample;
913  }
914 
915  if (m_useSampleHG[isample]) {
916  if (m_enablePreExcl) {
917  if (ADCHG > m_preExclHGADCThresh && isample < m_maxSamplesPreExcl) {
918  m_useSampleHG[isample] = false;
919  continue;
920  }
921  }
922 
923  if (m_enablePostExcl) {
924  if (ADCHG > m_postExclHGADCThresh && isample >= m_NSamplesAna - m_maxSamplesPostExcl) {
925  m_useSampleHG[isample] = false;
926  continue;
927  }
928  }
929 
930  // If we get here, we've passed the above filters on the sample
931  //
932  if (ADCHG > m_HGOverflowADC) {
933  m_HGOverflow = true;
934 
935  if (isample == m_preSampleIdx) m_PSHGOverUnderflow = true;
936 
937  // Note: the implementation of the explicit pre- and post- sample exclusion should make this
938  // code obselete but before removing it we first need to validate the pre- and post-exclusion
939  //
940  if ((int) isample > m_lastHGOverFlowSample) m_lastHGOverFlowSample = isample;
941  if ((int) isample < m_firstHGOverFlowSample) m_firstHGOverFlowSample = isample;
942  }
943  else if (ADCHG < m_HGUnderflowADC) {
944  m_HGUnderflow = true;
945 
946  if (isample == m_preSampleIdx) m_PSHGOverUnderflow = true;
947  }
948  }
949 
950  // Now the low gain sample
951  //
952  if (m_useSampleLG[isample]) {
953  if (m_enablePreExcl) {
954  if (ADCLG > m_preExclLGADCThresh && isample <= m_maxSamplesPreExcl) {
955  m_useSampleLG[isample] = false;
956  continue;
957  }
958  }
959 
960  if (m_enablePostExcl) {
961  if (ADCLG > m_postExclLGADCThresh && isample >= m_NSamplesAna - m_maxSamplesPostExcl) {
962  m_useSampleLG[isample] = false;
963  m_ExcludeLate = true;
964  continue;
965  }
966  }
967 
968  if (ADCLG > m_LGOverflowADC) {
969  m_LGOverflow = true;
970  m_fail = true;
971  m_amplitude = m_LGOverflowADC * m_gainFactorLG; // Give a vale here even though we know it's wrong because
972  // the user may not check the return value and we know that
973  // amplitude is bigger than this
974  }
975 
976  if (ADCLG == 0) {
977  m_LGUnderflow = true;
978  m_fail = true;
979  }
980  }
981  }
982 
983  // if (doDump) {
984  // (*m_msgFunc_p)(ZDCMsg::Verbose, "Dump of useSamples: ");
985 
986  // std::ostringstream dumpStringUseHG;
987  // dumpStringUseHG << "HG: ";
988  // for (auto val : m_useSampleHG) {
989  // dumpStringUseHG << val << " ";
990  // }
991  // (*m_msgFunc_p)(ZDCMsg::Verbose, dumpStringUseHG.str().c_str());
992 
993  // std::ostringstream dumpStringUseLG;
994  // dumpStringUseLG << "LG: ";
995  // for (auto val : m_useSampleLG) {
996  // dumpStringUseLG << val << " ";
997  // }
998  // (*m_msgFunc_p)(ZDCMsg::Verbose, dumpStringUseLG.str().c_str());
999  // }
1000 
1001  // This ugly code should be obseleted by the introduction of the better, pre- and post-sample exclusion but
1002  // that code still has to be fully validated.
1003  //
1004  // See if we can still use high gain even in the case of HG overflow if the overflow results from
1005  // the first few samples or the last few samples
1006  //
1007  if (m_HGOverflow && !m_HGUnderflow) {
1008  if (m_lastHGOverFlowSample < 2 && m_lastHGOverFlowSample > -1) {
1009  m_HGOverflow = false;
1011  m_adjTimeRangeEvent = true;
1012  m_backToHG_pre = true;
1013  m_ExcludeEarly = true;
1014  }
1015  else if (m_firstHGOverFlowSample < static_cast<int>(m_NSamplesAna) && m_firstHGOverFlowSample >= static_cast<int>(m_NSamplesAna - 2) ) {
1017  m_HGOverflow = false;
1018  m_adjTimeRangeEvent = true;
1019  m_ExcludeLate = true;
1020  }
1021  }
1022 
1023  // (*m_msgFunc_p)(ZDCMsg::Verbose, "ZDCPulseAnalyzer:: " + m_tag + ": ScanAndSubtractSamples done");
1024 
1025  return true;
1026 }
1027 
1029 {
1030  float deriv2ndThreshHG = 0;
1031  float deriv2ndThreshLG = 0;
1032 
1033  if (!repass) {
1035 
1036  deriv2ndThreshHG = m_peak2ndDerivMinThreshHG;
1037  deriv2ndThreshLG = m_peak2ndDerivMinThreshLG;
1038  }
1039  else {
1040  deriv2ndThreshHG = m_peak2ndDerivMinRepassHG;
1041  deriv2ndThreshLG = m_peak2ndDerivMinRepassLG;
1042  }
1043 
1045  if (m_useLowGain) {
1046  // (*m_msgFunc_p)(ZDCMsg::Verbose, "ZDCPulseAnalyzer:: " + m_tag + " using low gain data ");
1047 
1048  auto chisqCutLambda = [cut = m_chisqDivAmpCutLG,
1050  power = m_chisqDivAmpPowerLG]
1051  (float chisq, float amp, unsigned int fitNDoF)->bool
1052  {
1053  double ratio = chisq / (std::pow(amp, power) + offset);
1054  if (chisq/fitNDoF > 2 && ratio > cut) return false;
1055  else return true;
1056  };
1057 
1059  deriv2ndThreshLG, m_noiseSigLG, m_LGT0CorrParams,
1060  chisqCutLambda, m_T0CutLowLG, m_T0CutHighLG);
1061  if (result) {
1062  //
1063  // +++BAC
1064  //
1065  // I have removed some code that attempted a refit in case of failure of the chi-square cut
1066  // Instead, we should implement an analysis of the failure with possible re-fit afterwards
1067  // with possible change of the pre- or post-pulse condition
1068  //
1069  // --BAC
1070  //
1071 
1072  double amplCorrFactor = getAmplitudeCorrection(false);
1073 
1074  //
1075  // Multiply amplitude by gain factor
1076  //
1078  m_amplitude = m_fitAmplitude * amplCorrFactor * m_gainFactorLG;
1079  m_ampError = m_fitAmpError * amplCorrFactor * m_gainFactorLG;
1084 
1085  // BAC: also scale up the 2nd derivative by the gain factor so low and high gain can be treated on the same footing
1086  //
1088  }
1089 
1090  return result;
1091  }
1092  else {
1093  // (*m_msgFunc_p)(ZDCMsg::Verbose, "ZDCPulseAnalyzer:: " + m_tag + " using high gain data ");
1094  auto chisqCutLambda = [cut = m_chisqDivAmpCutHG,
1096  power = m_chisqDivAmpPowerHG, tag = m_tag]
1097  (float chisq, float amp, unsigned int fitNDoF)->bool
1098  {
1099  double ratio = chisq / (std::pow(amp, power) + offset);
1100 
1101  if (chisq/float(fitNDoF) > 2 && ratio > cut) return false;
1102  else return true;
1103  };
1104 
1106  deriv2ndThreshHG, m_noiseSigHG, m_HGT0CorrParams,
1107  chisqCutLambda, m_T0CutLowHG, m_T0CutHighHG);
1108  if (result) {
1109  // +++BAC
1110  //
1111  // I have removed some code that attempted a refit in case of failure of the chi-square cut
1112  // Instead, we should implement an analysis of the failure with possible re-fit afterwards
1113  // with possible change of the pre- or post-pulse condition
1114  //
1115  // --BAC
1116 
1124 
1126 
1127  // If we have a non-linear correction, apply it here
1128  //
1129  // We apply it as an inverse correction - i.e. we divide by a correction
1130  // term tha is a sum of coefficients times the ADC minus a reference
1131  // to a power. The lowest power is 1, the highest is deteremined by
1132  // the number of provided coefficients
1133  //
1134  double amplCorrFactor = getAmplitudeCorrection(true);
1135 
1136  m_amplitude *= amplCorrFactor;
1137  m_ampError *= amplCorrFactor;
1138  }
1139 
1140  // If LG refit has been requested, do it now
1141  //
1142  if (m_LGMode == LGModeRefitLG && m_havePulse) {
1144  DoFit(true);
1145 
1146  double amplCorrFactor = getAmplitudeCorrection(false);
1147  m_refitLGAmplCorr = m_refitLGAmpl*amplCorrFactor;
1148  }
1149 
1150  return result;
1151  }
1152 }
1153 
1154 bool ZDCPulseAnalyzer::AnalyzeData(size_t nSamples, size_t preSampleIdx,
1155  const std::vector<float>& samples, // The samples used for this event
1156  const std::vector<bool>& useSample, // The samples used for this event
1157  float peak2ndDerivMinThresh,
1158  float noiseSig,
1159  const std::vector<float>& t0CorrParams, // The parameters used to correct the t0
1160  ChisqCutLambdatype chisqCutLambda, // Lambda to perform the selection
1161  float minT0Corr, float maxT0Corr // The minimum and maximum corrected T0 values
1162  )
1163 {
1164 
1165  // We keep track of which sample we used to do the subtraction sepaate from m_minSampleEvt
1166  // because the corresponding time is the reference time we provide to the fit function when
1167  // we have pre-pulses and in case we change the m_minSampleEvt after doing the subtraction
1168  //
1169  // e.g. our refit when the chisquare cuts fails
1170  //
1171 
1172  // Find the first used sample in the event
1173  //
1174  bool haveFirst = false;
1175  unsigned int lastUsed = 0;
1176 
1177  for (unsigned int sample = preSampleIdx; sample < nSamples; sample++) {
1178  if (useSample[sample]) {
1179  //
1180  // We're going to use this sample in the analysis, update bookeeping
1181  //
1182  if (!haveFirst) {
1184  haveFirst = true;
1185  }
1186  else {
1187  lastUsed = sample;
1188  }
1189  }
1190  }
1191 
1192  if (lastUsed < m_maxSampleEvt) m_maxSampleEvt = lastUsed;
1193 
1194  // Check to see whether we've changed the range of samples used in the analysis
1195  //
1196  // Should be obseleted with reworking of the fitting using Root::Fit package
1197  //
1198  if (m_minSampleEvt > preSampleIdx) {
1199  m_ExcludeEarly = true;
1200  m_adjTimeRangeEvent = true;
1201  }
1202 
1203  if (m_maxSampleEvt < nSamples - 1) {
1204  m_adjTimeRangeEvent = true;
1205  m_ExcludeLate = true;
1206  }
1207 
1208  // Prepare for subtraction
1209  //
1211  m_preSample = samples[m_usedPresampIdx];
1212 
1213  // std::ostringstream pedMessage;
1214  // pedMessage << "Pedestal index = " << m_usedPresampIdx << ", value = " << m_preSample;
1215  // (*m_msgFunc_p)(ZDCMsg::Verbose, pedMessage.str().c_str());
1216 
1217  m_samplesSub = samples;
1218  m_samplesSig.assign(m_NSamplesAna, noiseSig);
1219 
1220 
1221  //
1222  // When we are combinig delayed and undelayed samples we have to deal with the fact that
1223  // the two readouts can have different noise and thus different baselines. Which is a huge
1224  // headache.
1225  //
1226  if (m_useDelayed) {
1227  if (m_useFixedBaseline) {
1229  }
1230  else {
1231  //
1232  // Use much-improved method to match delayed and undelayed baselines
1233  //
1235  std::ostringstream baselineMsg;
1236  baselineMsg << "Delayed samples baseline correction = " << m_baselineCorr << std::endl;
1237  (*m_msgFunc_p)(ZDCMsg::Debug, baselineMsg.str().c_str());
1238  }
1239 
1240  // Now apply the baseline correction to align ADC values for delayed and undelayed samples
1241  //
1242  for (size_t isample = 0; isample < nSamples; isample++) {
1243  if (isample % 2) m_samplesSub[isample] -= m_baselineCorr;
1244  }
1245 
1246  // If we use one of the delayed samples for presample, we have to adjust the presample value as well
1247  //
1249  }
1250 
1251  // Do the presample subtraction
1252  //
1253  std::for_each(m_samplesSub.begin(), m_samplesSub.end(), [ = ] (float & adcUnsub) {return adcUnsub -= m_preSample;} );
1254 
1255  // Calculate the second derivatives using step size m_2ndDerivStep
1256  //
1258 
1259  // Find the sample which has the lowest 2nd derivative. We loop over the range defined by the
1260  // tolerance on the position of the minimum second derivative. Note: the +1 in the upper iterator is because
1261  // that's where the loop terminates, not the last element.
1262  //
1263  SampleCIter minDeriv2ndIter;
1264 
1266 
1267  minDeriv2ndIter = std::min_element(m_samplesDeriv2nd.begin() + m_peak2ndDerivMinSample - m_peak2ndDerivMinTolerance, m_samplesDeriv2nd.begin() + upperDelta);
1268 
1269  m_minDeriv2nd = *minDeriv2ndIter;
1270  m_minDeriv2ndSig = -m_minDeriv2nd/(std::sqrt(6.0)*noiseSig);
1271 
1272  m_minDeriv2ndIndex = std::distance(m_samplesDeriv2nd.cbegin(), minDeriv2ndIter);
1273 
1274  // BAC 02-04-23 This check turned out to be problematic. Todo: figure out how to replace
1275  //
1276  // // Also check the ADC value for the "peak" sample to make sure it is significant (at least 3 sigma)
1277  // // The factor of sqrt(2) on the noise is because we have done a pre-sample subtraction
1278  // //
1279  if (std::abs(m_minDeriv2nd) >= peak2ndDerivMinThresh) {
1280  m_havePulse = true;
1281  }
1282  else {
1283  m_havePulse = false;
1284  }
1285 
1286  // save the low and high gain ADC values at the peak -- if we have a pulse, at m_minDeriv2ndIndex
1287  // otherwise at m_peak2ndDerivMinSample
1288  //
1289  if (m_havePulse) {
1292  }
1293  else {
1296  }
1297 
1298  // Now decide whether we have a preceeding pulse or not. There are two possible kinds of preceeding pulses:
1299  // 1) exponential tail from a preceeding pulse
1300  // 2) peaked pulse before the main pulse
1301  //
1302  // we can, of course, have both
1303  //
1304 
1305  // To check for exponential tail, test the slope determined by the minimum ADC value (and pre-sample)
1306  // **beware** this can cause trouble in 2015 data where the pulses had overshoot due to the transformers
1307  //
1308 
1309  // Implement a much simpler test for presence of an exponential tail from an OOT pulse
1310  // namely if the slope evaluated at the presample is significantly negative, then
1311  // we have a preceeding pulse.
1312  //
1313  // Note that we don't have to subtract the FADC value corresponding to the presample because
1314  // that subtraction has already been done.
1315  //
1316  if (m_havePulse) {
1317  // If we've alreday excluded early samples, we have almost by construction have negative exponential tail
1318  //
1319  if (m_ExcludeEarly) m_preExpTail = true;
1320 
1321  //
1322  // The subtracted ADC value at m_usedPresampIdx is, by construction, zero
1323  // The next sample has had the pre-sample subtracted, so it represents the initial derivative
1324  //
1325  float derivPresampleSig = m_samplesSub[m_usedPresampIdx+1]/(std::sqrt(2.0)*noiseSig);
1326  if (derivPresampleSig < -5) {
1327  m_preExpTail = true;
1328  m_preExpSig = derivPresampleSig;
1329  }
1330 
1331  for (unsigned int isample = m_usedPresampIdx; isample < m_samplesSub.size(); isample++) {
1332  if (!useSample[isample]) continue;
1333 
1334  float sampleSig = -m_samplesSub[isample]/(std::sqrt(2.0)*noiseSig);
1335 
1336  // Compare the derivative significant to the 2nd derivative significance,
1337  // so we don't waste time dealing with small perturbations on large signals
1338  //
1339  if ((sampleSig > 5 && sampleSig > 0.02*m_minDeriv2ndSig) || sampleSig > 0.5*m_minDeriv2ndSig) {
1340  m_preExpTail = true;
1341  if (sampleSig > m_preExpSig) m_preExpSig = sampleSig;
1342  }
1343  }
1344 
1345  // Now we search for maxima before the main pulse
1346  //
1347  int loopLimit = (m_havePulse ? m_minDeriv2ndIndex - 2 : m_peak2ndDerivMinSample - 2);
1348  int loopStart = m_minSampleEvt == 0 ? 1 : m_minSampleEvt;
1349 
1350  float maxPrepulseSig = 0;
1351  unsigned int maxPrepulseSample = 0;
1352 
1353  for (int isample = loopStart; isample <= loopLimit; isample++) {
1354  if (!useSample[isample]) continue;
1355 
1356  //
1357  // If any of the second derivatives prior to the peak are significantly negative, we have a an extra pulse
1358  // prior to the main one -- as opposed to just an expnential tail
1359  //
1360  float prePulseSig = -m_samplesDeriv2nd[isample]/(std::sqrt(6.0)*noiseSig);
1361 
1362  if ((prePulseSig > 6 && m_samplesDeriv2nd[isample] < 0.05 * m_minDeriv2nd) ||
1363  m_samplesDeriv2nd[isample] < 0.5*m_minDeriv2nd)
1364  {
1365  m_prePulse = true;
1366  if (prePulseSig > maxPrepulseSig) {
1367  maxPrepulseSig = prePulseSig;
1368  maxPrepulseSample = isample;
1369  }
1370  }
1371  }
1372 
1373  if (m_prePulse) {
1374  m_prePulseSig = maxPrepulseSig;
1375 
1376  if (m_preExpTail) {
1377  //
1378  // We have a prepulse. If we already indicated an negative exponential,
1379  // if the prepulse has greater significance, we override the negative exponential
1380  //
1381  if (m_prePulseSig > m_preExpSig) {
1382  m_preExpTail = false;
1383  }
1384  else {
1385  m_prePulse = false;
1386  }
1387  }
1388 
1389  m_initialPrePulseAmp = m_samplesSub[maxPrepulseSample];
1390  m_initialPrePulseT0 = m_deltaTSample * (maxPrepulseSample);
1391  }
1392 
1393  // if (m_preExpTail) m_prePulse = true;
1394 
1395  // -----------------------------------------------------
1396  // Post pulse detection
1397  //
1398  unsigned int postStartIdx = std::max(static_cast<unsigned int>(m_minDeriv2ndIndex + 2),
1399  static_cast<unsigned int>(m_peak2ndDerivMinSample + m_peak2ndDerivMinTolerance + 1));
1400 
1401 
1402  for (int isample = postStartIdx; isample < (int) nSamples - 1; isample++) {
1403  if (!useSample.at(isample)) continue;
1404 
1405  // +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
1406  // BAC 12-01-2024
1407  //
1408  // The following code is commented out as a temporary measure to deal with apparent reflections
1409  // associated with large out-of-time pulses that introduce a "kink" that triggers the derivative
1410  // test. A work-around that doesn't introduce specific code for the 2023 Pb+Pb run but allows
1411  // adaption for this specific issue is going to take some work. For now we leave the 2nd derivative
1412  // test, but will also need to introduce some configurability of the cut -- which we need anyway
1413  //
1414  // Calculate the forward derivative. the pulse should never increase on the tail. If it
1415  // does, we almost certainly have a post-pulse
1416  //
1417  // float deriv = m_samplesSub[isample + 1] - m_samplesSub[isample];
1418  // if (deriv/(std::sqrt(2)*noiseSig) > 6) {
1419  // m_postPulse = true;
1420  // m_maxSampleEvt = isample;
1421  // m_adjTimeRangeEvent = true;
1422  // break;
1423  // }
1424  // else {
1425  //----------------------------------------------------------------------------------------------
1426  //
1427  // Now we check the second derivative which might also indicate a post pulse
1428  // even if the derivative is not sufficiently large
1429  //
1430  // add small 1e-3 in division to avoid floating overflow
1431  //
1432 
1433  float deriv = m_samplesSub.at(isample + 1) - m_samplesSub.at(isample);
1434  float derivSig = deriv/(std::sqrt(2)*noiseSig);
1435  float deriv2ndSig = -m_samplesDeriv2nd.at(isample) / (std::sqrt(6)*noiseSig);
1436 
1437  float deriv2ndTest = m_samplesDeriv2nd.at(isample) / (-m_minDeriv2nd + 1.0e-3);
1438 
1439  if (derivSig > 5) {
1440  //
1441  // Check the 2nd derivative -- we should be at a minimum(?)
1442  //
1443  if (std::abs(deriv2ndTest) > 0.15) {
1444  m_postPulse = true;
1445  m_maxSampleEvt = std::min<int>(isample - (m_2ndDerivStep - 1), m_maxSampleEvt);
1446 
1447  m_adjTimeRangeEvent = true;
1448  break;
1449  }
1450  }
1451 
1452  // The place to apply the cut on samples depends on whether we have found a "minimum" or a "maximum"
1453  // The -m_2ndDerivStep for the minimum accounts for the shift between 2nd derivative and the samples
1454  // if we find a maximum we cut one sample lower
1455  //
1456  if (deriv2ndSig > 5 && deriv2ndTest < -0.1) {
1457  m_postPulse = true;
1458  m_maxSampleEvt = std::min<int>(isample - m_2ndDerivStep, m_maxSampleEvt);
1459 
1460  m_adjTimeRangeEvent = true;
1461  break;
1462  }
1463  }
1464  }
1465 
1466  if (m_postPulse) {
1467  std::ostringstream ostrm;
1468  ostrm << "Post pulse found, m_maxSampleEvt = " << m_maxSampleEvt;
1469  (*m_msgFunc_p)(ZDCMsg::Debug, ostrm.str());
1470  }
1471 
1472 
1473  // -----------------------------------------------------
1474 
1475  // Stop now if we have no pulse or we've detected a failure
1476  //
1477  if (m_fail || !m_havePulse) return false;
1478 
1479  if (!m_useDelayed) DoFit();
1480  else DoFitCombined();
1481 
1482  if (fitFailed()) {
1483  m_fail = true;
1484  }
1485  else {
1486  std::ostringstream ostrm;
1487  // ostrm << "Pulse fit successful with chisquare = " << m_fitChisq;
1488  // (*m_msgFunc_p)(ZDCMsg::Debug, ostrm.str());
1489 
1491 
1492  if (m_timingCorrMode != NoTimingCorr) {
1493  //
1494  // We correct relative to the m_timingCorrRefADC using m_timingCorrScale to scale
1495  //
1496  double t0CorrFact = (m_fitAmplitude - m_timingCorrRefADC) / m_timingCorrScale;
1497  if (m_timingCorrMode == TimingCorrLog) t0CorrFact = std::log(t0CorrFact);
1498 
1499  // Calculate the correction using a polynomial of power determined by the size of the vector.
1500  // For historical reasons we include here a constant term, though it is degenerate with
1501  // the nominal t0 and/or timing calibrations.
1502  //
1503  float correction = 0;
1504 
1505  for (unsigned int ipow = 0; ipow < t0CorrParams.size(); ipow++) {
1506  correction += t0CorrParams[ipow]*std::pow(t0CorrFact, double(ipow));
1507  }
1508 
1509  // The correction represents the offset of the timing value from zero so we subtract the result
1510  //
1512  }
1513 
1514  bool failFixedCut = m_fitTimeCorr < minT0Corr || m_fitTimeCorr > maxT0Corr;
1515 
1516  // Calculate the timing significance. Note this implementation breaks the model for
1517  // how DoAnalysis is currently used by explicitly testing m_useLowGain, but that model
1518  // needs adjustment anyway ... (one step at a time)
1519  //
1520  // Note that to avoid coupling non-linear corrections to the amplitude and the
1521  // timing significance we evaluate the time resolution using the fit amplitude
1522  //
1523  if (m_timeCutMode != 0) {
1524  double timeResolution = 0;
1525  if (m_useLowGain) timeResolution = m_timeResFuncLG_p->Eval(m_fitAmplitude);
1526  else timeResolution = m_timeResFuncHG_p->Eval(m_fitAmplitude);
1527 
1528  m_timeSig = m_fitTimeCorr/timeResolution;
1529  if (std::abs(m_timeSig) > m_t0CutSig) {
1530  //
1531  // We've failed the significance cut. In OR mode (1) we mark the time
1532  // as bad if we've lso failed the fixed cut. In AND mode, we mark it
1533  // as bad regardless of the fixed cut.
1534  //
1535  if (m_timeCutMode == 1) {
1536  if (failFixedCut) m_badT0 = true;
1537  }
1538  else m_badT0 = true;
1539  }
1540  else if (m_timeCutMode == 2 && failFixedCut) m_badT0 = failFixedCut;
1541  }
1542  else {
1543  m_timeSig = -1;
1544  m_badT0 = failFixedCut;
1545  }
1546 
1547  // Now check for valid chisq using lambda function
1548  //
1549  // if (m_fitChisq/m_fitNDoF > 2 && m_fitChisq / (m_fitAmplitude + 1.0e-6) > maxChisqDivAmp) m_badChisq = true;
1550  if (!chisqCutLambda(m_fitChisq, m_fitAmplitude, m_fitNDoF)) m_badChisq = true; }
1551 
1552  return !m_fitFailed;
1553 }
1554 
1555 void ZDCPulseAnalyzer::prepareLGRefit(const std::vector<float>& samplesLG, const std::vector<float>& samplesSig,
1556  const std::vector<bool>& useSamples)
1557 {
1558  m_samplesLGRefit.clear();
1559  m_samplesSigLGRefit.clear();
1560 
1561  float presampleLG = m_ADCSamplesLGSub[m_usedPresampIdx];
1562 
1563  // Do the presample subtraction
1564  //
1565  for (unsigned int idx = 0; idx < samplesLG.size(); idx++) {
1566  m_samplesLGRefit.push_back(samplesLG[idx] - presampleLG);
1567 
1568  if (useSamples[idx]) {
1569  m_samplesSigLGRefit.push_back(samplesSig[idx]);
1570  }
1571  else {
1572  m_samplesSigLGRefit.push_back(0);
1573  }
1574  }
1575 }
1576 
1577 
1578 void ZDCPulseAnalyzer::DoFit(bool refitLG)
1579 {
1580  bool fitLG = m_useLowGain || refitLG;
1581 
1582  TH1* hist_p = nullptr;
1583  FillHistogram(refitLG);
1584 
1585  // Set the initial values
1586  //
1587  float ampInitial, fitAmpMin, fitAmpMax, t0Initial;
1588 
1589  if (fitLG) {
1590  fitAmpMin = m_fitAmpMinLG;
1591  fitAmpMax = m_fitAmpMaxLG;
1592  t0Initial = m_nominalT0LG;
1593  ampInitial = m_ADCPeakLG;
1594  }
1595  else {
1596  ampInitial = m_ADCPeakHG;
1597  fitAmpMin = m_fitAmpMinHG;
1598  fitAmpMax = m_fitAmpMaxHG;
1599  t0Initial = m_nominalT0HG;
1600  }
1601 
1602  if (refitLG) {
1603  hist_p = m_fitHistLGRefit.get();
1604  }
1605  else {
1606  hist_p = m_fitHist.get();
1607  }
1608  if (ampInitial < fitAmpMin) ampInitial = fitAmpMin * 1.5;
1609 
1610 
1611  ZDCFitWrapper* fitWrapper = m_defaultFitWrapper.get();
1612  if (preExpTail()) {
1613  fitWrapper = m_preExpFitWrapper.get();
1614  (static_cast<ZDCPreExpFitWrapper*>(m_preExpFitWrapper.get()))->SetInitialExpPulse(m_initialExpAmp);
1615  }
1616  else if (prePulse()) {
1617  fitWrapper = m_prePulseFitWrapper.get();
1618  (static_cast<ZDCPrePulseFitWrapper*>(fitWrapper))->SetInitialPrePulse(m_initialPrePulseAmp, m_initialPrePulseT0,0,25);
1619  }
1620 
1621  if (m_adjTimeRangeEvent) {
1624 
1625  float fitTReference = m_deltaTSample * m_usedPresampIdx;
1626 
1627  fitWrapper->Initialize(ampInitial, t0Initial, fitAmpMin, fitAmpMax, m_fitTMin, m_fitTMax, fitTReference);
1628  }
1629  else {
1630  fitWrapper->Initialize(ampInitial, t0Initial, fitAmpMin, fitAmpMax);
1631  }
1632 
1633  // Now perform the fit
1634  //
1635  std::string options = m_fitOptions + "Ns";
1636  if (m_quietFits) {
1637  options += "Q";
1638  }
1639 
1640  bool fitFailed = false;
1641 
1642  // dumpTF1(fitWrapper->GetWrapperTF1RawPtr());
1643  checkTF1Limits(fitWrapper->GetWrapperTF1RawPtr());
1644 
1645  //
1646  // Fit the data with the function provided by the fit wrapper
1647  //
1648  TFitResultPtr result_ptr = hist_p->Fit(fitWrapper->GetWrapperTF1RawPtr(), options.c_str(), "", m_fitTMin, m_fitTMax);
1649  int fitStatus = result_ptr;
1650 
1651  //
1652  // If the first fit failed, also check the EDM. If sufficiently small, the failure is almost surely due
1653  // to parameter limits and we just accept the result
1654  //
1655  if (fitStatus != 0 && result_ptr->Edm() > 0.001)
1656  {
1657  //
1658  // We contstrain the fit and try again
1659  //
1660  fitWrapper->ConstrainFit();
1661 
1662  TFitResultPtr constrFitResult_ptr = hist_p->Fit(fitWrapper->GetWrapperTF1RawPtr(), options.c_str(), "", m_fitTMin, m_fitTMax);
1663  fitWrapper->UnconstrainFit();
1664 
1665  if ((int) constrFitResult_ptr != 0) {
1666  //
1667  // Even the constrained fit failed, so we quit.
1668  //
1669  fitFailed = true;
1670  }
1671  else {
1672  // Now we try the fit again with the constraint removed
1673  //
1674  TFitResultPtr unconstrFitResult_ptr = hist_p->Fit(fitWrapper->GetWrapperTF1RawPtr(), options.c_str(), "", m_fitTMin, m_fitTMax);
1675  if ((int) unconstrFitResult_ptr != 0) {
1676  //
1677  // The unconstrained fit failed again, so we redo the constrained fit
1678  //
1679  fitWrapper->ConstrainFit();
1680 
1681  TFitResultPtr constrFit2Result_ptr = hist_p->Fit(fitWrapper->GetWrapperTF1RawPtr(), options.c_str(), "", m_fitTMin, m_fitTMax);
1682  if ((int) constrFit2Result_ptr != 0) {
1683  //
1684  // Even the constrained fit failed the second time, so we quit.
1685  //
1686  fitFailed = true;
1687  }
1688 
1689  result_ptr = constrFit2Result_ptr;
1690  fitWrapper->UnconstrainFit();
1691  }
1692  else {
1693  result_ptr = unconstrFitResult_ptr;
1694  }
1695  }
1696  }
1697 
1698  if (!m_fitFailed && m_saveFitFunc) {
1699  hist_p->GetListOfFunctions()->Clear();
1700 
1701  TF1* func = fitWrapper->GetWrapperTF1RawPtr();
1702  std::string name = func->GetName();
1703 
1704  TF1* copyFunc = static_cast<TF1*>(func->Clone((name + "_copy").c_str()));
1705  hist_p->GetListOfFunctions()->Add(copyFunc);
1706  }
1707 
1708  if (!refitLG) {
1710  m_bkgdMaxFraction = fitWrapper->GetBkgdMaxFraction();
1711  m_fitAmplitude = fitWrapper->GetAmplitude();
1712  m_fitAmpError = fitWrapper->GetAmpError();
1713 
1714  if (preExpTail()) {
1715  m_fitExpAmp = (static_cast<ZDCPreExpFitWrapper*>(m_preExpFitWrapper.get()))->GetExpAmp();
1716  }
1717  else {
1718  m_fitExpAmp = 0;
1719  }
1720 
1721  m_fitTime = fitWrapper->GetTime();
1722  m_fitTimeSub = m_fitTime - t0Initial;
1723 
1724  m_fitChisq = result_ptr->Chi2();
1725  m_fitNDoF = result_ptr->Ndf();
1726 
1727  m_fitTau1 = fitWrapper->GetTau1();
1728  m_fitTau2 = fitWrapper->GetTau2();
1729 
1730  // Here we need to check if the fit amplitude is small (close) enough to fitAmpMin.
1731  // with "< 1+epsilon" where epsilon ~ 1%
1732  if (m_fitAmplitude < fitAmpMin * 1.01) {
1733  m_fitMinAmp = true;
1734  }
1735 
1736  if (m_haveSignifCuts) {
1737  float sigMinCut = fitLG ? m_sigMinLG : m_sigMinHG;
1738 
1739  if (m_fitAmpError > 1e-6) {
1740  if (m_fitAmplitude/m_fitAmpError < sigMinCut) m_failSigCut = true;
1741  }
1742  }
1743  }
1744  else {
1745  m_evtLGRefit = true;
1746  m_refitLGFitAmpl = fitWrapper->GetAmplitude();
1747  m_refitLGAmpl = fitWrapper->GetAmplitude() * m_gainFactorLG;
1748  m_refitLGAmpError = fitWrapper->GetAmpError() * m_gainFactorLG;
1749  m_refitLGChisq = result_ptr->Chi2();
1750  m_refitLGTime = fitWrapper->GetTime();
1751  m_refitLGTimeSub = m_refitLGTime - t0Initial;
1752  }
1753 }
1754 
1756 {
1757  bool fitLG = refitLG || m_useLowGain;
1758  TH1* hist_p = nullptr, *delayedHist_p = nullptr;
1759 
1760  FillHistogram(refitLG);
1761  if (refitLG) {
1762  hist_p = m_fitHistLGRefit.get();
1763  delayedHist_p = m_delayedHistLGRefit.get();
1764  }
1765  else {
1766  hist_p = m_fitHist.get();
1767  delayedHist_p = m_delayedHist.get();
1768  }
1769 
1770  float fitAmpMin, fitAmpMax, t0Initial, ampInitial;
1771  if (fitLG) {
1772  ampInitial = m_ADCPeakLG;
1773  fitAmpMin = m_fitAmpMinLG;
1774  fitAmpMax = m_fitAmpMaxLG;
1775  t0Initial = m_nominalT0LG;
1776  }
1777  else {
1778  ampInitial = m_ADCPeakHG;
1779  fitAmpMin = m_fitAmpMinHG;
1780  fitAmpMax = m_fitAmpMaxHG;
1781  t0Initial = m_nominalT0HG;
1782  }
1783 
1784  // Set the initial values
1785  //
1786  if (ampInitial < fitAmpMin) ampInitial = fitAmpMin * 1.5;
1787 
1788  ZDCFitWrapper* fitWrapper = m_defaultFitWrapper.get();
1789  //if (prePulse()) fitWrapper = fitWrapper = m_preExpFitWrapper.get();
1790 
1791  if (preExpTail()) {
1792  fitWrapper = m_preExpFitWrapper.get();
1793  (static_cast<ZDCPreExpFitWrapper*>(m_preExpFitWrapper.get()))->SetInitialExpPulse(m_initialExpAmp);
1794  }
1795  else if (prePulse()) {
1796  fitWrapper = m_prePulseFitWrapper.get();
1797  (static_cast<ZDCPrePulseFitWrapper*>(fitWrapper))->SetInitialPrePulse(m_initialPrePulseAmp, m_initialPrePulseT0,0,25);
1798  }
1799 
1800  // Initialize the fit wrapper for this event, specifying a
1801  // per-event fit range if necessary
1802  //
1803  if (m_adjTimeRangeEvent) {
1806 
1807  float fitTReference = m_deltaTSample * m_usedPresampIdx;
1808 
1809  fitWrapper->Initialize(ampInitial, t0Initial, fitAmpMin, fitAmpMax, m_fitTMin, m_fitTMax, fitTReference);
1810  }
1811  else {
1812  fitWrapper->Initialize(ampInitial, t0Initial, fitAmpMin, fitAmpMax);
1813  }
1814 
1815  // Set up the virtual fitter
1816  //
1817  TFitter* theFitter = nullptr;
1818 
1819  if (prePulse()) {
1821 
1822  theFitter = m_prePulseCombinedFitter.get();
1823  }
1824  else {
1826 
1827  theFitter = m_defaultCombinedFitter.get();
1828  }
1829 
1830  // dumpTF1(fitWrapper->GetWrapperTF1RawPtr());
1831 
1832  // Set the static pointers to histograms and function for use in FCN
1833  //
1834  s_undelayedFitHist = hist_p;
1835  s_delayedFitHist = delayedHist_p;
1836  s_combinedFitFunc = fitWrapper->GetWrapperTF1RawPtr();
1839 
1840 
1841  size_t numFitPar = theFitter->GetNumberTotalParameters();
1842 
1843  theFitter->GetMinuit()->fISW[4] = -1;
1844 
1845  // Now perform the fit
1846  //
1847  if (m_quietFits) {
1848  theFitter->GetMinuit()->fISW[4] = -1;
1849 
1850  int ierr= 0;
1851  theFitter->GetMinuit()->mnexcm("SET NOWarnings",nullptr,0,ierr);
1852  }
1853  else theFitter->GetMinuit()->fISW[4] = 0;
1854 
1855  // Only include baseline shift in fit for pre-pulses. Otherwise baseline matching should work
1856  //
1857  if (prePulse()) {
1858  theFitter->SetParameter(0, "delayBaselineAdjust", 0, 0.01, -100, 100);
1859  theFitter->ReleaseParameter(0);
1860  }
1861  else {
1862  theFitter->SetParameter(0, "delayBaselineAdjust", 0, 0.01, -100, 100);
1863  theFitter->FixParameter(0);
1864  }
1865 
1866 
1867  double arglist[100];
1868  arglist[0] = 5000; // number of function calls
1869  arglist[1] = 0.01; // tolerance
1870  int status = theFitter->ExecuteCommand("MIGRAD", arglist, 2);
1871 
1872  double fitAmp = theFitter->GetParameter(1);
1873 
1874  // Capture the chi-square etc.
1875  //
1876  double chi2, edm, errdef;
1877  int nvpar, nparx;
1878 
1879  theFitter->GetStats(chi2, edm, errdef, nvpar, nparx);
1880 
1881  // Here we need to check if fitAmp is small (close) enough to fitAmpMin.
1882  // with "< 1+epsilon" where epsilon ~ 1%
1883  //
1884  if (status || fitAmp < fitAmpMin * 1.01 || edm > 0.01){
1885 
1886  //
1887  // We first retry the fit with no baseline adjust
1888  //
1889  theFitter->SetParameter(0, "delayBaselineAdjust", 0, 0.01, -100, 100);
1890  theFitter->FixParameter(0);
1891 
1892  // Here we need to check if fitAmp is small (close) enough to fitAmpMin.
1893  // with "< 1+epsilon" where epsilon ~ 1%
1894  if (fitAmp < fitAmpMin * 1.01) {
1895  if (m_adjTimeRangeEvent) {
1896  float fitTReference = m_deltaTSample * m_usedPresampIdx;
1897  fitWrapper->Initialize(ampInitial, t0Initial, fitAmpMin, fitAmpMax, m_fitTMin, m_fitTMax, fitTReference);
1898  }
1899  else {
1900  fitWrapper->Initialize(ampInitial, t0Initial, fitAmpMin, fitAmpMax);
1901  }
1902  }
1903 
1904  status = theFitter->ExecuteCommand("MIGRAD", arglist, 2);
1905  if (status != 0) {
1906  //
1907  // Fit failed event with no baseline adust, so be it
1908  //
1909  theFitter->ReleaseParameter(0);
1910  m_fitFailed = true;
1911  }
1912  else {
1913  //
1914  // The fit succeeded with no baseline adjust, re-fit allowing for baseline adjust using
1915  // the parameters from previous fit as the starting point
1916  //
1917  theFitter->ReleaseParameter(0);
1918  status = theFitter->ExecuteCommand("MIGRAD", arglist, 2);
1919 
1920  if (status) {
1921  // Since we know the fit can succeed without the baseline adjust, go back to it
1922  //
1923  theFitter->SetParameter(0, "delayBaselineAdjust", 0, 0.01, -100, 100);
1924  theFitter->FixParameter(0);
1925  status = theFitter->ExecuteCommand("MIGRAD", arglist, 2);
1926  }
1927  }
1928  }
1929  else m_fitFailed = false;
1930 
1931  // Check to see if the fit forced the amplitude to the minimum, if so set the corresponding status bit
1932  //
1933  fitAmp = theFitter->GetParameter(1);
1934 
1935  // Here we need to check if fitAmp is small (close) enough to fitAmpMin.
1936  // with "< 1+epsilon" where epsilon ~ 1%
1937  if (fitAmp < fitAmpMin * 1.01) {
1938  m_fitMinAmp = true;
1939  }
1940 
1941  // Check that the amplitude passes minimum significance requirement
1942  //
1943  float sigMinCut = fitLG ? m_sigMinLG : m_sigMinHG;
1944  if (m_fitAmpError > 1e-6) {
1945  if (m_fitAmplitude/m_fitAmpError < sigMinCut) m_failSigCut = true;
1946  }
1947 
1948  if (!m_quietFits) theFitter->GetMinuit()->fISW[4] = -1;
1949 
1950  std::vector<double> funcParams(numFitPar - 1);
1951  std::vector<double> funcParamErrs(numFitPar - 1);
1952 
1953  // Save the baseline shift between delayed and undelayed samples
1954  //
1955  m_delayedBaselineShift = theFitter->GetParameter(0);
1956 
1957  // Capture and store the fit function parameteds and errors
1958  //
1959  for (size_t ipar = 1; ipar < numFitPar; ipar++) {
1960  funcParams[ipar - 1] = theFitter->GetParameter(ipar);
1961  funcParamErrs[ipar - 1] = theFitter->GetParError(ipar);
1962  }
1963 
1964  s_combinedFitFunc->SetParameters(&funcParams[0]);
1965  s_combinedFitFunc->SetParErrors(&funcParamErrs[0]);
1966 
1967  // Capture the chi-square etc.
1968  //
1969  theFitter->GetStats(chi2, edm, errdef, nvpar, nparx);
1970 
1971  int ndf = 2 * m_Nsample - nvpar;
1972 
1973  s_combinedFitFunc->SetChisquare(chi2);
1974  s_combinedFitFunc->SetNDF(ndf);
1975 
1976  // add to list of functions
1977  if (m_saveFitFunc) {
1978  s_undelayedFitHist->GetListOfFunctions()->Clear();
1979  s_undelayedFitHist->GetListOfFunctions()->Add(s_combinedFitFunc);
1980 
1981  s_delayedFitHist->GetListOfFunctions()->Clear();
1982  s_delayedFitHist->GetListOfFunctions()->Add(s_combinedFitFunc);
1983  }
1984 
1985  if (!refitLG) {
1986  // Save the pull values from the last call to FCN
1987  //
1988  arglist[0] = 3; // number of function calls
1989  theFitter->ExecuteCommand("Cal1fcn", arglist, 1);
1991 
1992  m_fitAmplitude = fitWrapper->GetAmplitude();
1993  m_fitTime = fitWrapper->GetTime();
1994  if (prePulse()) {
1995  m_fitPreT0 = (static_cast<ZDCPrePulseFitWrapper*>(m_prePulseFitWrapper.get()))->GetPreT0();
1996  m_fitPreAmp = (static_cast<ZDCPrePulseFitWrapper*>(m_prePulseFitWrapper.get()))->GetPreAmp();
1997  m_fitPostT0 = (static_cast<ZDCPrePulseFitWrapper*>(m_prePulseFitWrapper.get()))->GetPostT0();
1998  m_fitPostAmp = (static_cast<ZDCPrePulseFitWrapper*>(m_prePulseFitWrapper.get()))->GetPostAmp();
1999  }
2000 
2001  if (preExpTail()) {
2002  m_fitExpAmp = (static_cast<ZDCPreExpFitWrapper*>(m_preExpFitWrapper.get()))->GetExpAmp();
2003  }
2004  else {
2005  m_fitExpAmp = 0;
2006  }
2007 
2008  m_fitTimeSub = m_fitTime - t0Initial;
2009  m_fitChisq = chi2;
2010  m_fitNDoF = ndf;
2011 
2012  m_fitTau1 = fitWrapper->GetTau1();
2013  m_fitTau2 = fitWrapper->GetTau2();
2014 
2015  m_fitAmpError = fitWrapper->GetAmpError();
2016  m_bkgdMaxFraction = fitWrapper->GetBkgdMaxFraction();
2017  }
2018  else {
2019  m_evtLGRefit = true;
2020  m_refitLGFitAmpl = fitWrapper->GetAmplitude();
2021  m_refitLGAmpl = fitWrapper->GetAmplitude() * m_gainFactorLG;
2022  m_refitLGAmpError = fitWrapper->GetAmpError() * m_gainFactorLG;
2023  m_refitLGChisq = chi2;
2024  m_refitLGTime = fitWrapper->GetTime();
2025  m_refitLGTimeSub = m_refitLGTime - t0Initial;
2026  }
2027 }
2028 
2030 {
2031  for (int ipar = 0; ipar < func->GetNpar(); ipar++) {
2032  double parLimitLow, parLimitHigh;
2033 
2034  func->GetParLimits(ipar, parLimitLow, parLimitHigh);
2035 
2036  (*m_msgFunc_p)(ZDCMsg::Debug, (
2037  "ZDCPulseAnalyzer name=" + std::string(func->GetName())
2038  + " ipar=" + std::to_string(ipar)
2039  + " parLimitLow=" + std::to_string(parLimitLow)
2040  + " parLimitHigh="+ std::to_string(parLimitHigh)
2041  )
2042  );
2043 
2044  //if (std::abs(parLimitHigh / parLimitLow - 1) > 1e-6) {
2045  if (std::abs(parLimitHigh - parLimitLow) > (1e-6)*std::abs(parLimitLow)) {
2046  double value = func->GetParameter(ipar);
2047  if (value >= parLimitHigh) {
2048  value = parLimitHigh * 0.9;
2049  }
2050  else if (value <= parLimitLow) {
2051  value = parLimitLow + 0.1*std::abs(parLimitLow);
2052  }
2053  func->SetParameter(ipar, value);
2054  }
2055  }
2056 }
2057 
2058 std::unique_ptr<TFitter> ZDCPulseAnalyzer::MakeCombinedFitter(TF1* func)
2059 {
2060  TVirtualFitter::SetDefaultFitter("Minuit");
2061 
2062  size_t nFitParams = func->GetNpar() + 1;
2063  std::unique_ptr<TFitter> fitter = std::make_unique<TFitter>(nFitParams);
2064 
2065  fitter->GetMinuit()->fISW[4] = -1;
2066  fitter->SetParameter(0, "delayBaselineAdjust", 0, 0.01, -100, 100);
2067 
2068  for (size_t ipar = 0; ipar < nFitParams - 1; ipar++) {
2069  double parLimitLow, parLimitHigh;
2070 
2071  func->GetParLimits(ipar, parLimitLow, parLimitHigh);
2072  if (std::abs(parLimitHigh / parLimitLow - 1) < 1e-6) {
2073  double value = func->GetParameter(ipar);
2074  double lowLim = std::min(value * 0.99, value * 1.01);
2075  double highLim = std::max(value * 0.99, value * 1.01);
2076 
2077  fitter->SetParameter(ipar + 1, func->GetParName(ipar), func->GetParameter(ipar), 0.01, lowLim, highLim);
2078  fitter->FixParameter(ipar + 1);
2079  }
2080  else {
2081  double value = func->GetParameter(ipar);
2082  if (value >= parLimitHigh) value = parLimitHigh * 0.99;
2083  else if (value <= parLimitLow) value = parLimitLow * 1.01;
2084 
2085  double step = std::min((parLimitHigh - parLimitLow)/100., (value - parLimitLow)/100.);
2086 
2087  fitter->SetParameter(ipar + 1, func->GetParName(ipar), value, step, parLimitLow, parLimitHigh);
2088  }
2089  }
2090 
2092 
2093  return fitter;
2094 }
2095 
2096 void ZDCPulseAnalyzer::UpdateFitterTimeLimits(TFitter* fitter, ZDCFitWrapper* wrapper, bool prePulse)
2097 {
2098  double parLimitLow, parLimitHigh;
2099 
2100  auto func_p = wrapper->GetWrapperTF1();
2101  func_p->GetParLimits(1, parLimitLow, parLimitHigh);
2102 
2103  fitter->SetParameter(2, func_p->GetParName(1), func_p->GetParameter(1), 0.01, parLimitLow, parLimitHigh);
2104 
2105  if (prePulse) {
2106  unsigned int parIndex = static_cast<ZDCPrePulseFitWrapper*>(wrapper)->GetPreT0ParIndex();
2107 
2108  func_p->GetParLimits(parIndex, parLimitLow, parLimitHigh);
2109  fitter->SetParameter(parIndex + 1, func_p->GetParName(parIndex), func_p->GetParameter(parIndex), 0.01, parLimitLow, parLimitHigh);
2110  }
2111 }
2112 
2114 {
2115  (*m_msgFunc_p)(ZDCMsg::Info, ("ZDCPulseAnalyzer dump for tag = " + m_tag));
2116 
2117  (*m_msgFunc_p)(ZDCMsg::Info, ("Presample index, value = " + std::to_string(m_preSampleIdx) + ", " + std::to_string(m_preSample)));
2118 
2119  if (m_useDelayed) {
2120  (*m_msgFunc_p)(ZDCMsg::Info, ("using delayed samples with delta T = " + std::to_string(m_delayedDeltaT) + ", and pedestalDiff == " +
2122  }
2123 
2124  std::ostringstream message1;
2125  message1 << "samplesSub ";
2126  for (size_t sample = 0; sample < m_samplesSub.size(); sample++) {
2127  message1 << ", [" << sample << "] = " << m_samplesSub[sample];
2128  }
2129  (*m_msgFunc_p)(ZDCMsg::Info, message1.str());
2130 
2131  std::ostringstream message3;
2132  message3 << "samplesDeriv2nd ";
2133  for (size_t sample = 0; sample < m_samplesDeriv2nd.size(); sample++) {
2134  message3 << ", [" << sample << "] = " << m_samplesDeriv2nd[sample];
2135  }
2136  (*m_msgFunc_p)(ZDCMsg::Info, message3.str());
2137 
2138  (*m_msgFunc_p)(ZDCMsg::Info, ("minimum 2nd deriv sample, value = " + std::to_string(m_minDeriv2ndIndex) + ", " + std::to_string(m_minDeriv2nd)));
2139 }
2140 
2141 void ZDCPulseAnalyzer::dumpTF1(const TF1* func) const
2142 {
2143  std::string message = "Dump of TF1: " + std::string(func->GetName());
2144  bool continueDump = (*m_msgFunc_p)(ZDCMsg::Verbose, std::move(message));
2145  if (!continueDump) return;
2146 
2147  unsigned int npar = func->GetNpar();
2148  for (unsigned int ipar = 0; ipar < npar; ipar++) {
2149  std::ostringstream msgstr;
2150 
2151  double parMin = 0, parMax = 0;
2152  func->GetParLimits(ipar, parMin, parMax);
2153 
2154  msgstr << "Parameter " << ipar << ", value = " << func->GetParameter(ipar) << ", error = "
2155  << func->GetParError(ipar) << ", min = " << parMin << ", max = " << parMax;
2156  (*m_msgFunc_p)(ZDCMsg::Verbose, msgstr.str());
2157  }
2158 }
2159 
2161 {
2162  std::ostringstream ostrStream;
2163 
2164  (*m_msgFunc_p)(ZDCMsg::Info, ("ZDCPulserAnalyzer:: settings for instance: " + m_tag));
2165 
2166  ostrStream << "Nsample = " << m_Nsample << " at frequency " << m_freqMHz << " MHz, preSample index = "
2167  << m_preSampleIdx << ", nominal pedestal = " << m_pedestal;
2168 
2169  (*m_msgFunc_p)(ZDCMsg::Info, ostrStream.str()); ostrStream.str(""); ostrStream.clear();
2170 
2171  ostrStream << "LG mode = " << m_LGMode << ", gainFactor HG = " << m_gainFactorHG << ", gainFactor LG = " << m_gainFactorLG << ", noise sigma HG = " << m_noiseSigHG << ", noiseSigLG = " << m_noiseSigLG;
2172  (*m_msgFunc_p)(ZDCMsg::Info, ostrStream.str()); ostrStream.str(""); ostrStream.clear();
2173 
2174  ostrStream << "peak sample = " << m_peak2ndDerivMinSample << ", tolerance = " << m_peak2ndDerivMinTolerance
2175  << ", 2ndDerivThresh HG, LG = " << m_peak2ndDerivMinThreshHG << ", " << m_peak2ndDerivMinThreshLG
2176  << ", 2nd deriv step = " << m_2ndDerivStep;
2177  (*m_msgFunc_p)(ZDCMsg::Info, ostrStream.str()); ostrStream.str(""); ostrStream.clear();
2178 
2179  if (m_useDelayed) {
2180  ostrStream << "using delayed samples with delta T = " << m_delayedDeltaT << ", and default pedestalDiff == "
2182  (*m_msgFunc_p)(ZDCMsg::Info, ostrStream.str()); ostrStream.str(""); ostrStream.clear();
2183  }
2184 
2185  ostrStream <<"Fit function = " << m_fitFunction << "fixTau1 = " << m_fixTau1 << ", fixTau2 = " << m_fixTau2
2186  << ", nominalTau1 = " << m_nominalTau1 << ", nominalTau2 = " << m_nominalTau2 << "\n"
2187  << ", nominalT0HG = " << m_nominalT0HG << ", nominalT0LG = " << m_nominalT0LG
2188  << ", t0Cuts HG = [" << m_T0CutLowHG << ", " << m_T0CutHighHG << "], t0Cuts LG = ["
2189  << m_T0CutLowLG << ", " << m_T0CutHighLG << "]";
2190  (*m_msgFunc_p)(ZDCMsg::Info, ostrStream.str()); ostrStream.str(""); ostrStream.clear();
2191 
2192  ostrStream << "HGOverflowADC = " << m_HGOverflowADC << ", HGUnderflowADC = " << m_HGUnderflowADC
2193  << ", LGOverflowADC = "<< m_LGOverflowADC;
2194  (*m_msgFunc_p)(ZDCMsg::Info, ostrStream.str()); ostrStream.str(""); ostrStream.clear();
2195 
2196  ostrStream << "chisqDivAmpCutHG = " << m_chisqDivAmpCutHG << ", chisqDivAmpCutLG=" << m_chisqDivAmpCutLG
2197  << ", chisqDivAmpOffsetHG = " << m_chisqDivAmpOffsetHG << ", chisqDivAmpOffsetLG = " << m_chisqDivAmpOffsetLG
2198  << ", chisqDivAmpPowerHG = " << m_chisqDivAmpPowerHG << ", chisqDivAmpPowerLG = " << m_chisqDivAmpPowerLG;
2199  (*m_msgFunc_p)(ZDCMsg::Info, ostrStream.str()); ostrStream.str(""); ostrStream.clear();
2200 
2201  if ( m_enableRepass) {
2202  ostrStream << "Repass enabled with peak2ndDerivMinRepassHG = " << m_peak2ndDerivMinRepassHG << ", peak2ndDerivMinRepassLG = " << m_peak2ndDerivMinRepassLG;
2203  (*m_msgFunc_p)(ZDCMsg::Info, ostrStream.str()); ostrStream.str(""); ostrStream.clear();
2204  }
2205 
2206  if (m_enablePreExcl) {
2207  ostrStream << "Pre-exclusion enabled for up to " << m_maxSamplesPreExcl << ", samples with ADC threshold HG = "
2208  << m_preExclHGADCThresh << ", LG = " << m_preExclLGADCThresh;
2209  (*m_msgFunc_p)(ZDCMsg::Info, ostrStream.str()); ostrStream.str(""); ostrStream.clear();
2210  }
2211  if (m_enablePostExcl) {
2212  ostrStream << "Post-exclusion enabled for up to " << m_maxSamplesPostExcl << ", samples with ADC threshold HG = "
2213  << m_postExclHGADCThresh << ", LG = " << m_postExclLGADCThresh;
2214  (*m_msgFunc_p)(ZDCMsg::Info, ostrStream.str()); ostrStream.str(""); ostrStream.clear();
2215  }
2216  if (m_haveSignifCuts) {
2217  ostrStream << "Minimum significance cuts applied: HG min. sig. = " << m_sigMinHG << ", LG min. sig. " << m_sigMinLG;
2218  (*m_msgFunc_p)(ZDCMsg::Info, ostrStream.str()); ostrStream.str(""); ostrStream.clear();
2219  }
2220 
2221 }
2222 
2224 {
2225  unsigned int statusMask = 0;
2226 
2227  if (havePulse()) statusMask |= 1 << PulseBit;
2228  if (useLowGain()) statusMask |= 1 << LowGainBit;
2229  if (failed()) statusMask |= 1 << FailBit;
2230  if (HGOverflow()) statusMask |= 1 << HGOverflowBit;
2231 
2232  if (HGUnderflow()) statusMask |= 1 << HGUnderflowBit;
2233  if (PSHGOverUnderflow()) statusMask |= 1 << PSHGOverUnderflowBit;
2234  if (LGOverflow()) statusMask |= 1 << LGOverflowBit;
2235  if (LGUnderflow()) statusMask |= 1 << LGUnderflowBit;
2236 
2237  if (prePulse()) statusMask |= 1 << PrePulseBit;
2238  if (postPulse()) statusMask |= 1 << PostPulseBit;
2239  if (fitFailed()) statusMask |= 1 << FitFailedBit;
2240  if (badChisq()) statusMask |= 1 << BadChisqBit;
2241 
2242  if (badT0()) statusMask |= 1 << BadT0Bit;
2243  if (excludeEarlyLG()) statusMask |= 1 << ExcludeEarlyLGBit;
2244  if (excludeLateLG()) statusMask |= 1 << ExcludeLateLGBit;
2245  if (preExpTail()) statusMask |= 1 << preExpTailBit;
2246  if (fitMinimumAmplitude()) statusMask |= 1 << FitMinAmpBit;
2247  if (repassPulse()) statusMask |= 1 << RepassPulseBit;
2248  if (armSumInclude()) statusMask |= 1 << ArmSumIncludeBit;
2249  if (failSigCut()) statusMask |= 1 << FailSigCutBit;
2250 
2251  return statusMask;
2252 }
2253 
2254 std::shared_ptr<TGraphErrors> ZDCPulseAnalyzer::GetCombinedGraph(bool LGRefit)
2255 {
2256  //
2257  // We defer filling the histogram if we don't have a pulse until the histogram is requested
2258  //
2259  GetHistogramPtr(LGRefit);
2260 
2261  TH1* hist_p = nullptr, *delayedHist_p = nullptr;
2262  if (LGRefit) {
2263  hist_p = m_fitHistLGRefit.get();
2264  delayedHist_p = m_delayedHistLGRefit.get();
2265  }
2266  else {
2267  hist_p = m_fitHist.get();
2268  delayedHist_p = m_delayedHist.get();
2269  }
2270 
2271  std::shared_ptr<TGraphErrors> theGraph = std::make_shared<TGraphErrors>(TGraphErrors(2 * m_Nsample));
2272  size_t npts = 0;
2273 
2274  for (int ipt = 0; ipt < hist_p->GetNbinsX(); ipt++) {
2275  theGraph->SetPoint(npts, hist_p->GetBinCenter(ipt + 1), hist_p->GetBinContent(ipt + 1));
2276  theGraph->SetPointError(npts++, 0, hist_p->GetBinError(ipt + 1));
2277  }
2278 
2279  for (int iDelayPt = 0; iDelayPt < delayedHist_p->GetNbinsX(); iDelayPt++) {
2280  theGraph->SetPoint(npts, delayedHist_p->GetBinCenter(iDelayPt + 1), delayedHist_p->GetBinContent(iDelayPt + 1) - m_delayedBaselineShift);
2281  theGraph->SetPointError(npts++, 0, delayedHist_p->GetBinError(iDelayPt + 1));
2282  }
2283  if (m_havePulse) {
2284  TF1* func_p = static_cast<TF1*>(hist_p->GetListOfFunctions()->Last());
2285  if (func_p) {
2286  theGraph->GetListOfFunctions()->Add(new TF1(*func_p));
2287  hist_p->GetListOfFunctions()->SetOwner (false);
2288  }
2289  }
2290  theGraph->SetName(( std::string(hist_p->GetName()) + "combinaed").c_str());
2291 
2292  theGraph->SetMarkerStyle(20);
2293  theGraph->SetMarkerColor(1);
2294 
2295  return theGraph;
2296 }
2297 
2298 
2299 std::shared_ptr<TGraphErrors> ZDCPulseAnalyzer::GetGraph(bool forceLG)
2300 {
2301  //
2302  // We defer filling the histogram if we don't have a pulse until the histogram is requested
2303  //
2304  const TH1* hist_p = GetHistogramPtr(forceLG);
2305 
2306  std::shared_ptr<TGraphErrors> theGraph = std::make_shared<TGraphErrors>(TGraphErrors(m_Nsample));
2307  size_t npts = 0;
2308 
2309  for (int ipt = 0; ipt < hist_p->GetNbinsX(); ipt++) {
2310  theGraph->SetPoint(npts, hist_p->GetBinCenter(ipt + 1), hist_p->GetBinContent(ipt + 1));
2311  theGraph->SetPointError(npts++, 0, hist_p->GetBinError(ipt + 1));
2312  }
2313 
2314  TF1* func_p = static_cast<TF1*>(hist_p->GetListOfFunctions()->Last());
2315  theGraph->GetListOfFunctions()->Add(func_p);
2316  theGraph->SetName(( std::string(hist_p->GetName()) + "not_combinaed").c_str());
2317 
2318  theGraph->SetMarkerStyle(20);
2319  theGraph->SetMarkerColor(1);
2320 
2321  return theGraph;
2322 }
2323 
2324 
2325 std::vector<float> ZDCPulseAnalyzer::CalculateDerivative(const std::vector <float>& inputData, unsigned int step)
2326 {
2327  unsigned int nSamples = inputData.size();
2328 
2329  // So we pad at the beginning and end based on step (i.e. with step - 1 zeros). Fill out the fill vector with zeros initially
2330  //
2331  unsigned int vecSize = 2*(step - 1) + nSamples - step - 1;
2332  std::vector<float> results(vecSize, 0);
2333 
2334  // Now fill out the values
2335  //
2336  unsigned int fillIdx = step - 1;
2337 
2338  for (unsigned int sample = 0; sample < nSamples - step; sample++) {
2339  int deriv = inputData[sample + step] - inputData[sample];
2340  results.at(fillIdx++) = deriv;
2341  }
2342 
2343  return results;
2344 }
2345 
2346 std::vector<float> ZDCPulseAnalyzer::Calculate2ndDerivative(const std::vector <float>& inputData, unsigned int step)
2347 {
2348  unsigned int nSamples = inputData.size();
2349 
2350  // We start with two zero entries for which we can't calculate the double-step derivative
2351  // and would pad with two zero entries at the end. Start by initializing
2352  //
2353  unsigned int vecSize = 2*step + nSamples - step - 1;
2354  std::vector<float> results(vecSize, 0);
2355 
2356  unsigned int fillIndex = step;
2357  for (unsigned int sample = step; sample < nSamples - step; sample++) {
2358  int deriv2nd = inputData[sample + step] + inputData[sample - step] - 2*inputData[sample];
2359  results.at(fillIndex++) = deriv2nd;
2360  }
2361 
2362  return results;
2363 }
2364 
2365 // Implement a more general method for the nasty problem (Runs 1 & 2 only) of matching
2366 // the baselines on the delayed and undelayed data. Instead of specifically looking
2367 // for early samples or late samples, find the region of the waveform with the smallest
2368 // combination of slope and second derivative in both sets of samples and match there.
2369 //
2370 // Depending on the second derivative
2371 // we match using linear or (negative) exponential interpolation
2372 //
2373 // Note: the samples have already been combined so we distinguish by even or odd index
2374 // we use step = 2 for the derivative and 2nd derivative calculation to handle
2375 // the delayed and undelayed separately.
2376 //
2377 //
2378 float ZDCPulseAnalyzer::obtainDelayedBaselineCorr(const std::vector<float>& samples)
2379 {
2380  const unsigned int nsamples = samples.size();
2381 
2382  std::vector<float> derivVec = CalculateDerivative(samples, 2);
2383  std::vector<float> deriv2ndVec = Calculate2ndDerivative(samples, 2);
2384 
2385  // Now step through and check even and odd samples values for 2nd derivative and derivative
2386  // we start with index 2 since the 2nd derivative calculation has 2 initial zeros with nstep = 2
2387  //
2388  float minScore = 1.0e9;
2389  unsigned int minIndex = 0;
2390 
2391  for (unsigned int idx = 2; idx < nsamples - 1; idx++) {
2392  float deriv = derivVec[idx];
2393  float prevDeriv = derivVec[idx - 1];
2394 
2395  float derivDiff = deriv - prevDeriv;
2396 
2397  float deriv2nd = deriv2ndVec[idx];
2398  if (idx > nsamples - 2) deriv2nd = deriv2ndVec[idx - 1];
2399 
2400  // Calculate a score based on the actual derivatives (squared) and 2nd derivatives (squared)
2401  // and the slope differences (squared). The relative weights are not adjustable for now
2402  //
2403  float score = (deriv*deriv + 2*derivDiff*derivDiff +
2404  0.5*deriv2nd*deriv2nd);
2405 
2406  if (score < minScore) {
2407  minScore = score;
2408  minIndex = idx;
2409  }
2410  }
2411 
2412  // We use four samples, two each of "even" and "odd".
2413  // Because of the way the above analysis is done, we can always
2414  // Go back one even and one odd sample and forward one odd sample.
2415  //
2416  //if minIndex is < 2 or >samples.size() the result is undefined; prevent this:
2417  if (minIndex<2 or (minIndex+1) >=nsamples){
2418  throw std::out_of_range("minIndex out of range in ZDCPulseAnalyzer::obtainDelayedBaselineCorr");
2419  }
2420  float sample0 = samples[minIndex - 2];
2421  float sample1 = samples[minIndex - 1];
2422  float sample2 = samples[minIndex];
2423  float sample3 = samples[minIndex + 1];
2424 
2425  // Possibility -- implement logarithmic interpolation for large 2nd derivative?
2426  //
2427  float baselineCorr = (0.5 * (sample1 - sample0 + sample3 - sample2) -
2428  0.25 * (sample3 - sample1 + sample2 - sample0));
2429 
2430  if (minIndex % 2 != 0) baselineCorr =-baselineCorr;
2431 
2432  return baselineCorr;
2433 }
2434 
2435 std::pair<bool, std::string> ZDCPulseAnalyzer::ValidateJSONConfig(const JSON& config)
2436 {
2437  bool result = true;
2438  std::string resultString = "success";
2439 
2440  for (auto [key, descr] : JSONConfigParams) {
2441  auto iter = config.find(key);
2442  if (iter != config.end()) {
2443  //
2444  // Check type consistency
2445  //
2446  auto jsonType = iter.value().type();
2447  auto paramType = std::get<0>(descr);
2448  if (jsonType != paramType) {
2449  result = false;
2450  resultString = "Bad type for parameter " + key + ", type in JSON = " + std::to_string((unsigned int) jsonType) ;
2451  break;
2452  }
2453 
2454  size_t paramSize = std::get<1>(descr);
2455  size_t jsonSize = iter.value().size();
2456  if (jsonSize != paramSize) {
2457  result = false;
2458  resultString = "Bad length for parameter " + key + ", length in JSON = " + std::to_string(jsonSize) ;
2459  break;
2460  }
2461  }
2462  else {
2463  bool required = std::get<2>(descr);
2464  if (required) {
2465  result = false;
2466  resultString = "Missing required parameter " + key;
2467  break;
2468  }
2469  }
2470  }
2471 
2472  if (result) {
2473  //
2474  // Now check that the parameters in the JSON object are in the master list of parameters
2475  //
2476  for (auto [key, value] : config.items()) {
2477  //
2478  // Look for this key in the list of allowed parameters. Yes, it's a slow
2479  // search but we only do it once at configuration time.
2480  //
2481  // bool found = false;
2482  auto iter = JSONConfigParams.find(key);
2483  if (iter == JSONConfigParams.end()) {
2484  result = false;
2485  resultString = "Unknown parameter, key = " + key;
2486  break;
2487  }
2488  }
2489  }
2490 
2491  return {result, resultString};
2492 }
2493 
2494 std::pair<bool, std::string> ZDCPulseAnalyzer::ConfigFromJSON(const JSON& config)
2495 {
2496  bool result = true;
2497  std::string resultString = "success";
2498 
2499  for (auto [key, value] : config.items()) {
2500  //
2501  // Big if statement to process configuration parameters
2502  //
2503  if (key == "Nsample") m_Nsample = value;
2504  else if (key == "tag") m_tag = value;
2505  else if (key == "LGMode") m_LGMode = value;
2506  else if (key == "Nsample") m_Nsample = value;
2507  else if (key == "useDelayed") m_useDelayed = value;
2508  else if (key == "preSampleIdx") m_preSampleIdx = value;
2509  else if (key == "FADCFreqMHz") m_freqMHz = value;
2510  else if (key == "nominalPedestal") m_pedestal = value;
2511  else if (key == "fitFunction") m_fitFunction = value;
2512  else if (key == "peakSample") m_peak2ndDerivMinSample = value;
2513  else if (key == "peakTolerance") m_peak2ndDerivMinTolerance = value;
2514  else if (key == "2ndDerivThreshHG") m_peak2ndDerivMinThreshHG = value;
2515  else if (key == "2ndDerivThreshLG") m_peak2ndDerivMinThreshLG = value;
2516  else if (key == "2ndDerivStep") m_2ndDerivStep = value;
2517  else if (key == "HGOverflowADC") m_HGOverflowADC = value;
2518  else if (key == "HGUnderflowADC") m_HGUnderflowADC = value;
2519  else if (key == "LGOverflowADC") m_LGOverflowADC = value;
2520  else if (key == "nominalT0HG") m_nominalT0HG = value;
2521  else if (key == "nominalT0LG") m_nominalT0LG = value;
2522  else if (key == "nominalTau1") m_nominalTau1 = value;
2523  else if (key == "nominalTau2") m_nominalTau2 = value;
2524  else if (key == "fixTau1") m_fixTau1 = value;
2525  else if (key == "fixTau2") m_fixTau2 = value;
2526  else if (key == "T0CutsHG") {
2527  m_T0CutLowHG = value[0];
2528  m_T0CutHighHG = value[1];
2529  }
2530  else if (key == "T0CutsLG") {
2531  m_T0CutLowLG = value[0];
2532  m_T0CutHighLG = value[1];
2533  }
2534  else if (key == "chisqDivAmpCutHG") m_chisqDivAmpCutHG = value;
2535  else if (key == "chisqDivAmpCutLG") m_chisqDivAmpCutLG = value;
2536  else if (key == "chisqDivAmpOffsetHG") m_chisqDivAmpOffsetHG = value;
2537  else if (key == "chisqDivAmpOffsetLG") m_chisqDivAmpOffsetLG = value;
2538  else if (key == "chisqDivAmpPowerHG") m_chisqDivAmpPowerHG = value;
2539  else if (key == "chisqDivAmpPowerLG") m_chisqDivAmpPowerLG = value;
2540  else if (key == "gainFactorHG") m_gainFactorHG = value;
2541  else if (key == "gainFactorLG") m_gainFactorLG = value;
2542  else if (key == "noiseSigmaHG") m_noiseSigHG = value;
2543  else if (key == "noiseSigmaLG") m_noiseSigLG = value;
2544  else if (key == "enableRepass") m_enableRepass = value;
2545  else if (key == "Repass2ndDerivThreshHG") m_peak2ndDerivMinRepassHG = value;
2546  else if (key == "Repass2ndDerivThreshLG") m_peak2ndDerivMinRepassLG = value;
2547  else if (key == "fitAmpMinMaxHG") {
2548  m_fitAmpMinHG = value[0];
2549  m_fitAmpMaxHG = value[1];
2550  }
2551  else if (key == "fitAmpMinMaxLG") {
2552  m_fitAmpMinLG = value[0];
2553  m_fitAmpMaxLG = value[1];
2554  }
2555  else if (key == "quietFits") {
2556  m_quietFits = value[0];
2557  }
2558  else if (key == "enablePreExclusion") {
2559  m_enablePreExcl = true;
2563  }
2564  else if (key == "enablePostExclusion") {
2565  m_enablePostExcl = true;
2569  }
2570  else if (key == "ampMinSignifHGLG") {
2571  m_haveSignifCuts = true;
2572  m_sigMinHG = value[0];
2573  m_sigMinLG = value[1];
2574  }
2575  else if (key == "enableFADCCorrections") {
2576  auto fileNameJson = value["filename"];
2577  auto doPerSampleCorrJson = value["doPerSampleCorr"];
2578 
2579  if (fileNameJson.is_null() || doPerSampleCorrJson.is_null()) {
2580  result = false;
2581  std::string resultString = "failure processing enableFADCCorrections object";
2582  break;
2583  }
2584 
2585  m_haveFADCCorrections = true;
2586  m_FADCCorrPerSample = doPerSampleCorrJson;
2587  m_fadcCorrFileName = fileNameJson;
2588  }
2589  else {
2590  result = false;
2591  std::string resultString = "unprocessed parameter";
2592  break;
2593  }
2594  }
2595 
2596  return {result, resultString};
2597 }
ZDCFitWrapper::GetAmplitude
virtual float GetAmplitude() const =0
ZDCPulseAnalyzer::SetupFitFunctions
void SetupFitFunctions()
Definition: ZDCPulseAnalyzer.cxx:662
ZDCPulseAnalyzer::m_LGOverflow
bool m_LGOverflow
Definition: ZDCPulseAnalyzer.h:248
ZDCPulseAnalyzer::m_fitTau2
float m_fitTau2
Definition: ZDCPulseAnalyzer.h:318
ZDCPulseAnalyzer::m_fitPostT0
float m_fitPostT0
Definition: ZDCPulseAnalyzer.h:323
ZDCPulseAnalyzer::RepassPulseBit
@ RepassPulseBit
Definition: ZDCPulseAnalyzer.h:49
ZDCPulseAnalyzer::m_defaultFitWrapper
std::unique_ptr< ZDCFitWrapper > m_defaultFitWrapper
Definition: ZDCPulseAnalyzer.h:211
covarianceTool.ndf
ndf
Definition: covarianceTool.py:678
ZDCPulseAnalyzer::m_failSigCut
bool m_failSigCut
Definition: ZDCPulseAnalyzer.h:264
AllowedVariables::e
e
Definition: AsgElectronSelectorTool.cxx:37
ZDCPulseAnalyzer::m_defaultT0Max
float m_defaultT0Max
Definition: ZDCPulseAnalyzer.h:161
ZDCPulseAnalyzer::m_expAmplitude
float m_expAmplitude
Definition: ZDCPulseAnalyzer.h:332
ZDCPulseAnalyzer::m_minADCSampleLG
int m_minADCSampleLG
Definition: ZDCPulseAnalyzer.h:283
LArSamples::FitterData::fitter
const ShapeFitter * fitter
Definition: ShapeFitter.cxx:23
ZDCPulseAnalyzer::preExpTail
bool preExpTail() const
Definition: ZDCPulseAnalyzer.h:609
ZDCPulseAnalyzer::m_prePulseFitWrapper
std::unique_ptr< ZDCPrePulseFitWrapper > m_prePulseFitWrapper
Definition: ZDCPulseAnalyzer.h:212
ZDCPulseAnalyzer::m_useSampleHG
std::vector< bool > m_useSampleHG
Definition: ZDCPulseAnalyzer.h:355
ZDCPulseAnalyzer::m_enableRepass
bool m_enableRepass
Definition: ZDCPulseAnalyzer.h:110
ZDCPulseAnalyzer::SetGainFactorsHGLG
void SetGainFactorsHGLG(float gainFactorHG, float gainFactorLG)
Definition: ZDCPulseAnalyzer.cxx:456
ZDCPulseAnalyzer::m_ADCSamplesLG
std::vector< float > m_ADCSamplesLG
Definition: ZDCPulseAnalyzer.h:350
ZDCPulseAnalyzer::FitMinAmpBit
@ FitMinAmpBit
Definition: ZDCPulseAnalyzer.h:48
ZDCPulseAnalyzer::m_samplesSub
std::vector< float > m_samplesSub
Definition: ZDCPulseAnalyzer.h:360
ZDCPulseAnalyzer::m_fitPostAmp
float m_fitPostAmp
Definition: ZDCPulseAnalyzer.h:324
ZDCPulseAnalyzer::m_chisqDivAmpCutLG
float m_chisqDivAmpCutLG
Definition: ZDCPulseAnalyzer.h:143
ZDCPulseAnalyzer::enableRepass
void enableRepass(float peak2ndDerivMinRepassHG, float peak2ndDerivMinRepassLG)
Definition: ZDCPulseAnalyzer.cxx:231
ZDCPulseAnalyzer::m_NSamplesAna
unsigned int m_NSamplesAna
Definition: ZDCPulseAnalyzer.h:348
ZDCFitExpFermiVariableTausLHCf
Definition: ZDCFitWrapper.h:231
createLinkingScheme.iter
iter
Definition: createLinkingScheme.py:62
ZDCPulseAnalyzer::m_HGOverflowADC
int m_HGOverflowADC
Definition: ZDCPulseAnalyzer.h:127
ZDCPulseAnalyzer::m_enablePreExcl
bool m_enablePreExcl
Definition: ZDCPulseAnalyzer.h:176
ZDCPulseAnalyzer::m_LGT0CorrParams
std::vector< float > m_LGT0CorrParams
Definition: ZDCPulseAnalyzer.h:190
get_generator_info.result
result
Definition: get_generator_info.py:21
ZDCPulseAnalyzer::m_fitAmplitude
float m_fitAmplitude
Definition: ZDCPulseAnalyzer.h:310
ZDCPulseAnalyzer::m_HGOverflow
bool m_HGOverflow
Definition: ZDCPulseAnalyzer.h:244
AddEmptyComponent.histName
string histName
Definition: AddEmptyComponent.py:64
ZDCPulseAnalyzer::m_fitTau1
float m_fitTau1
Definition: ZDCPulseAnalyzer.h:317
ZDCPulseAnalyzer::m_useFixedBaseline
bool m_useFixedBaseline
Definition: ZDCPulseAnalyzer.h:224
ZDCPulseAnalyzer::m_postAmplitude
float m_postAmplitude
Definition: ZDCPulseAnalyzer.h:331
ZDCPulseAnalyzer::m_adjTimeRangeEvent
bool m_adjTimeRangeEvent
Definition: ZDCPulseAnalyzer.h:217
ZDCPulseAnalyzer::FillHistogram
void FillHistogram(bool refitLG)
Definition: ZDCPulseAnalyzer.h:410
ZDCPulseAnalyzer::m_haveFADCCorrections
bool m_haveFADCCorrections
Definition: ZDCPulseAnalyzer.h:199
ZDCPulseAnalyzer::m_minDeriv2nd
float m_minDeriv2nd
Definition: ZDCPulseAnalyzer.h:293
ZDCPulseAnalyzer::ReanalyzeData
bool ReanalyzeData()
Definition: ZDCPulseAnalyzer.cxx:832
ZDCPrePulseFitWrapper
Definition: ZDCFitWrapper.h:116
ZDCPulseAnalyzer::m_FADCCorrHG
std::unique_ptr< const TH1 > m_FADCCorrHG
Definition: ZDCPulseAnalyzer.h:202
ZDCFitExpFermiFixedTaus
Definition: ZDCFitWrapper.h:298
ZDCPulseAnalyzer::m_peak2ndDerivMinThreshLG
float m_peak2ndDerivMinThreshLG
Definition: ZDCPulseAnalyzer.h:105
ZDCPulseAnalyzer::m_2ndDerivStep
size_t m_2ndDerivStep
Definition: ZDCPulseAnalyzer.h:102
ZDCPulseAnalyzer::m_quietFits
bool m_quietFits
Definition: ZDCPulseAnalyzer.h:98
ZDCPulseAnalyzer::m_nonLinCorrRefScale
float m_nonLinCorrRefScale
Definition: ZDCPulseAnalyzer.h:195
ZDCPulseAnalyzer::m_timeCutMode
unsigned int m_timeCutMode
Definition: ZDCPulseAnalyzer.h:159
ZDCPulseAnalyzer::GetHistogramPtr
const TH1 * GetHistogramPtr(bool refitLG=false)
Definition: ZDCPulseAnalyzer.h:721
ZDCPulseAnalyzer::MakeCombinedFitter
static std::unique_ptr< TFitter > MakeCombinedFitter(TF1 *func)
Definition: ZDCPulseAnalyzer.cxx:2058
ZDCPulseAnalyzer::m_fitAmpMaxHG
float m_fitAmpMaxHG
Definition: ZDCPulseAnalyzer.h:167
ZDCPulseAnalyzer::m_backToHG_pre
bool m_backToHG_pre
Definition: ZDCPulseAnalyzer.h:268
ZDCFitWrapper
Definition: ZDCFitWrapper.h:19
TRT::Track::chiSquare
@ chiSquare
Definition: InnerDetector/InDetCalibEvent/TRT_CalibData/TRT_CalibData/TrackInfo.h:67
ZDCPulseAnalyzer::ArmSumIncludeBit
@ ArmSumIncludeBit
Definition: ZDCPulseAnalyzer.h:50
ZDCPulseAnalyzer::FailSigCutBit
@ FailSigCutBit
Definition: ZDCPulseAnalyzer.h:51
ZDCPulseAnalyzer.h
ZDCPulseAnalyzer::m_refitLGAmpl
float m_refitLGAmpl
Definition: ZDCPulseAnalyzer.h:337
ZDCPulseAnalyzer::m_fitTimeSub
float m_fitTimeSub
Definition: ZDCPulseAnalyzer.h:313
ZDCPulseAnalyzer::FailBit
@ FailBit
Definition: ZDCPulseAnalyzer.h:30
ZDCPulseAnalyzer::m_delayedPedestalDiff
float m_delayedPedestalDiff
Definition: ZDCPulseAnalyzer.h:226
ZDCPulseAnalyzer::enableTimeSigCut
void enableTimeSigCut(bool AND, float sigCut, const std::string &TF1String, const std::vector< double > &parsHG, const std::vector< double > &parsLG)
Definition: ZDCPulseAnalyzer.cxx:522
max
constexpr double max()
Definition: ap_fixedTest.cxx:33
ZDCPulseAnalyzer::DoFitCombined
void DoFitCombined(bool refitLG=false)
Definition: ZDCPulseAnalyzer.cxx:1755
ZDCFitWrapper::ConstrainFit
virtual void ConstrainFit()=0
ZDCPulseAnalyzer::m_sigMinHG
float m_sigMinHG
Definition: ZDCPulseAnalyzer.h:171
ZDCPulseAnalyzer::LGOverflow
bool LGOverflow() const
Definition: ZDCPulseAnalyzer.h:598
ZDCFitWrapper::GetTime
virtual float GetTime() const =0
min
constexpr double min()
Definition: ap_fixedTest.cxx:26
ZDCPulseAnalyzer::m_timingCorrRefADC
float m_timingCorrRefADC
Definition: ZDCPulseAnalyzer.h:188
ZDCPulseAnalyzer::UpdateFitterTimeLimits
void UpdateFitterTimeLimits(TFitter *fitter, ZDCFitWrapper *wrapper, bool prePulse)
Definition: ZDCPulseAnalyzer.cxx:2096
ZDCPulseAnalyzer::GetStatusMask
unsigned int GetStatusMask() const
Definition: ZDCPulseAnalyzer.cxx:2223
ZDCPulseAnalyzer::m_postPulse
bool m_postPulse
Definition: ZDCPulseAnalyzer.h:252
ZDCPulseAnalyzer::m_defaultCombinedFitter
std::unique_ptr< TFitter > m_defaultCombinedFitter
Definition: ZDCPulseAnalyzer.h:231
ZDCPulseAnalyzer::m_nominalT0HG
float m_nominalT0HG
Definition: ZDCPulseAnalyzer.h:131
ZDCPulseAnalyzer::PostPulseBit
@ PostPulseBit
Definition: ZDCPulseAnalyzer.h:39
ZDCPulseAnalyzer::m_minADCLG
float m_minADCLG
Definition: ZDCPulseAnalyzer.h:282
ZDCPulseAnalyzer::m_delayedHist
std::unique_ptr< TH1 > m_delayedHist
Definition: ZDCPulseAnalyzer.h:227
ZDCPulseAnalyzer::SetCutValues
void SetCutValues(float chisqDivAmpCutHG, float chisqDivAmpCutLG, float deltaT0MinHG, float deltaT0MaxHG, float deltaT0MinLG, float deltaT0MaxLG)
Definition: ZDCPulseAnalyzer.cxx:508
ZDCPulseAnalyzer::m_fitHistLGRefit
std::unique_ptr< TH1 > m_fitHistLGRefit
Definition: ZDCPulseAnalyzer.h:208
ZDCPulseAnalyzer::m_haveNonlinCorr
bool m_haveNonlinCorr
Definition: ZDCPulseAnalyzer.h:193
ZDCPulseAnalyzer::m_preExpTail
bool m_preExpTail
Definition: ZDCPulseAnalyzer.h:259
ZDCPulseAnalyzer::m_bkgdMaxFraction
float m_bkgdMaxFraction
Definition: ZDCPulseAnalyzer.h:333
ZDCPulseAnalyzer::m_tag
std::string m_tag
Definition: ZDCPulseAnalyzer.h:88
ZDCPulseAnalyzer::ScanAndSubtractSamples
bool ScanAndSubtractSamples()
Definition: ZDCPulseAnalyzer.cxx:844
ZDCFitExpFermiLinearPrePulse
Definition: ZDCFitWrapper.h:850
Root::AND
@ AND
Definition: TGRLCollection.h:32
ZDCPulseAnalyzer::m_sigMinLG
float m_sigMinLG
Definition: ZDCPulseAnalyzer.h:172
ZDCPulseAnalyzer::LowGainBit
@ LowGainBit
Definition: ZDCPulseAnalyzer.h:29
athena.value
value
Definition: athena.py:124
ZDCPulseAnalyzer::ZDCPulseAnalyzer
ZDCPulseAnalyzer(ZDCMsg::MessageFunctionPtr msgFunc_p, const std::string &tag, int Nsample, float deltaTSample, size_t preSampleIdx, int pedestal, const std::string &fitFunction, int peak2ndDerivMinSample, float peak2DerivMinThreshHG, float peak2DerivMinThreshLG)
Definition: ZDCPulseAnalyzer.cxx:142
ZDCPreExpFitWrapper
Definition: ZDCFitWrapper.h:144
ZDCPulseAnalyzer::m_prePulseSig
float m_prePulseSig
Definition: ZDCPulseAnalyzer.h:303
ZDCPulseAnalyzer::m_HGUnderflow
bool m_HGUnderflow
Definition: ZDCPulseAnalyzer.h:246
ZDCPulseAnalyzer::useLowGain
bool useLowGain() const
Definition: ZDCPulseAnalyzer.h:592
ZDCPulseAnalyzer::SetDefaults
void SetDefaults()
Definition: ZDCPulseAnalyzer.cxx:238
ZDCPulseAnalyzer::armSumInclude
bool armSumInclude() const
Definition: ZDCPulseAnalyzer.h:612
ZDCPulseAnalyzer::HGUnderflowBit
@ HGUnderflowBit
Definition: ZDCPulseAnalyzer.h:33
ZDCFitWrapper::GetWrapperTF1
virtual std::shared_ptr< TF1 > GetWrapperTF1()
Definition: ZDCFitWrapper.h:110
ZDCPulseAnalyzer::m_fixTau2
bool m_fixTau2
Definition: ZDCPulseAnalyzer.h:138
read_hist_ntuple.t
t
Definition: read_hist_ntuple.py:5
const
bool const RAWDATA *ch2 const
Definition: LArRodBlockPhysicsV0.cxx:560
ZDCJSONConfig::JSONParamList
std::map< std::string, JSONParamDescr > JSONParamList
Definition: ZDCJSONConfig.h:20
ZDCPulseAnalyzer::m_defaultFitTMax
float m_defaultFitTMax
Definition: ZDCPulseAnalyzer.h:140
MuonR4::to_string
std::string to_string(const SectorProjector proj)
Definition: MsTrackSeeder.cxx:66
ZDCPulseAnalyzer::m_delayedBaselineShift
float m_delayedBaselineShift
Definition: ZDCPulseAnalyzer.h:334
ZDCPulseAnalyzer::m_maxADCSampleHG
int m_maxADCSampleHG
Definition: ZDCPulseAnalyzer.h:279
ZDCPulseAnalyzer::m_badT0
bool m_badT0
Definition: ZDCPulseAnalyzer.h:256
ZDCPulseAnalyzer::m_fitMinAmp
bool m_fitMinAmp
Definition: ZDCPulseAnalyzer.h:262
ZDCMsg::Info
@ Info
Definition: ZDCMsg.h:20
ZDCPulseAnalyzer::s_delayedFitHist
static TH1 * s_delayedFitHist
Definition: ZDCPulseAnalyzer.h:79
ZDCPulseAnalyzer::BadChisqBit
@ BadChisqBit
Definition: ZDCPulseAnalyzer.h:41
ReweightUtils.message
message
Definition: ReweightUtils.py:15
ZDCPulseAnalyzer::m_samplesDeriv2nd
std::vector< float > m_samplesDeriv2nd
Definition: ZDCPulseAnalyzer.h:366
ZDCFitExpFermiLHCfPreExp
Definition: ZDCFitWrapper.h:672
ZDCPulseAnalyzer::ChisqCutLambdatype
std::function< bool(float, float, float)> ChisqCutLambdatype
Definition: ZDCPulseAnalyzer.h:389
ZDCPulseAnalyzer::repassPulse
bool repassPulse() const
Definition: ZDCPulseAnalyzer.h:611
ZDCFitWrapper::GetTau2
virtual float GetTau2() const =0
ZDCPulseAnalyzer::m_haveData
bool m_haveData
Definition: ZDCPulseAnalyzer.h:239
config
Definition: PhysicsAnalysis/AnalysisCommon/AssociationUtils/python/config.py:1
ZDCPulseAnalyzer::badChisq
bool badChisq() const
Definition: ZDCPulseAnalyzer.h:604
ZDCPulseAnalyzer::PSHGOverUnderflowBit
@ PSHGOverUnderflowBit
Definition: ZDCPulseAnalyzer.h:34
ZDCPulseAnalyzer::HGOverflow
bool HGOverflow() const
Definition: ZDCPulseAnalyzer.h:594
ZDCPulseAnalyzer::m_maxADCSampleLG
int m_maxADCSampleLG
Definition: ZDCPulseAnalyzer.h:284
ZDCPulseAnalyzer::s_combinedFitFunc
static TF1 * s_combinedFitFunc
Definition: ZDCPulseAnalyzer.h:80
ZDCPulseAnalyzer::m_ADCPeakHG
float m_ADCPeakHG
Definition: ZDCPulseAnalyzer.h:286
ZDCPulseAnalyzer::LGUnderflowBit
@ LGUnderflowBit
Definition: ZDCPulseAnalyzer.h:36
ZDCPulseAnalyzer::m_nominalTau2
float m_nominalTau2
Definition: ZDCPulseAnalyzer.h:135
ZDCPulseAnalyzer::m_prePulse
bool m_prePulse
Definition: ZDCPulseAnalyzer.h:251
ZDCPulseAnalyzer::s_combinedFitTMin
static float s_combinedFitTMin
Definition: ZDCPulseAnalyzer.h:82
ZDCPulseAnalyzer::dumpConfiguration
void dumpConfiguration() const
Definition: ZDCPulseAnalyzer.cxx:2160
ZDCPulseAnalyzer::PSHGOverUnderflow
bool PSHGOverUnderflow() const
Definition: ZDCPulseAnalyzer.h:597
ZDCPulseAnalyzer::m_useSampleLG
std::vector< bool > m_useSampleLG
Definition: ZDCPulseAnalyzer.h:354
ZDCPulseAnalyzer::m_ampNoNonLin
float m_ampNoNonLin
Definition: ZDCPulseAnalyzer.h:327
ZDCPulseAnalyzer::m_fitFunction
std::string m_fitFunction
Definition: ZDCPulseAnalyzer.h:101
ZDCPulseAnalyzer::LoadAndAnalyzeData
bool LoadAndAnalyzeData(const std::vector< float > &ADCSamplesHG, const std::vector< float > &ADCSamplesLG)
Definition: ZDCPulseAnalyzer.cxx:755
ZDCPulseAnalyzer::m_preExclHGADCThresh
unsigned int m_preExclHGADCThresh
Definition: ZDCPulseAnalyzer.h:178
ZDCPulseAnalyzer::m_pedestal
int m_pedestal
Definition: ZDCPulseAnalyzer.h:93
ZDCPulseAnalyzer::SampleCIter
std::vector< float >::const_iterator SampleCIter
Definition: ZDCPulseAnalyzer.h:73
ZDCPulseAnalyzer::m_fixPrePulse
bool m_fixPrePulse
Definition: ZDCPulseAnalyzer.h:261
ZDCPulseAnalyzer::m_fitPulls
std::vector< float > m_fitPulls
Definition: ZDCPulseAnalyzer.h:373
tools.zlumi_mc_cf.correction
def correction(mu, runmode, campaign, run=None)
Definition: zlumi_mc_cf.py:4
ZDCPulseAnalyzer::m_preSample
float m_preSample
Definition: ZDCPulseAnalyzer.h:274
ZDCPulseAnalyzer::TimingCorrLog
@ TimingCorrLog
Definition: ZDCPulseAnalyzer.h:61
ZDCPulseAnalyzer::m_fitOptions
std::string m_fitOptions
Definition: ZDCPulseAnalyzer.h:126
ZDCMsg::Verbose
@ Verbose
Definition: ZDCMsg.h:18
ZDCPulseAnalyzer::ValidateJSONConfig
std::pair< bool, std::string > ValidateJSONConfig(const JSON &config)
Definition: ZDCPulseAnalyzer.cxx:2435
ZDCFitWrapper::UnconstrainFit
virtual void UnconstrainFit()=0
ZDCPulseAnalyzer::m_postExclHGADCThresh
unsigned int m_postExclHGADCThresh
Definition: ZDCPulseAnalyzer.h:182
ZDCPulseAnalyzer::LGOverflowBit
@ LGOverflowBit
Definition: ZDCPulseAnalyzer.h:35
ZDCPulseAnalyzer::m_nominalT0LG
float m_nominalT0LG
Definition: ZDCPulseAnalyzer.h:132
ZDCPulseAnalyzer::m_minADCHG
float m_minADCHG
Definition: ZDCPulseAnalyzer.h:276
ZDCPulseAnalyzer::m_fitAmpMinLG
float m_fitAmpMinLG
Definition: ZDCPulseAnalyzer.h:165
ZDCPulseAnalyzer::m_badChisq
bool m_badChisq
Definition: ZDCPulseAnalyzer.h:254
ZDCPulseAnalyzer::m_t0CutSig
float m_t0CutSig
Definition: ZDCPulseAnalyzer.h:158
ZDCPulseAnalyzer::m_saveFitFunc
bool m_saveFitFunc
Definition: ZDCPulseAnalyzer.h:99
ZDCFitWrapper::GetAmpError
virtual float GetAmpError() const =0
ZDCPulseAnalyzer::m_peak2ndDerivMinRepassHG
float m_peak2ndDerivMinRepassHG
Definition: ZDCPulseAnalyzer.h:112
ZDCPulseAnalyzer::m_initialPrePulseAmp
float m_initialPrePulseAmp
Definition: ZDCPulseAnalyzer.h:306
ZDCPulseAnalyzer::m_evtLGRefit
bool m_evtLGRefit
Definition: ZDCPulseAnalyzer.h:336
ZDCPulseAnalyzer::PrePulseBit
@ PrePulseBit
Definition: ZDCPulseAnalyzer.h:38
ZDCPulseAnalyzer::m_peak2ndDerivMinRepassLG
float m_peak2ndDerivMinRepassLG
Definition: ZDCPulseAnalyzer.h:111
ZDCPulseAnalyzer::m_fitPostT0lo
float m_fitPostT0lo
Definition: ZDCPulseAnalyzer.h:299
FullCPAlgorithmsTest_eljob.sample
sample
Definition: FullCPAlgorithmsTest_eljob.py:116
ZDCPulseAnalyzer::m_fixTau1
bool m_fixTau1
Definition: ZDCPulseAnalyzer.h:137
ZDCFitWrapper::GetWrapperTF1RawPtr
virtual TF1 * GetWrapperTF1RawPtr() const
Definition: ZDCFitWrapper.h:112
ZDCJSONConfig::JSON
nlohmann::json JSON
Definition: ZDCJSONConfig.h:18
xmin
double xmin
Definition: listroot.cxx:60
ZDCPulseAnalyzer::m_fitTMax
float m_fitTMax
Definition: ZDCPulseAnalyzer.h:296
ZDCPulseAnalyzer::m_fitChisq
float m_fitChisq
Definition: ZDCPulseAnalyzer.h:319
ZDCPulseAnalyzer::CombinedPulsesFCN
static void CombinedPulsesFCN(int &numParam, double *, double &f, double *par, int flag)
Definition: ZDCPulseAnalyzer.cxx:88
ZDCPulseAnalyzer::m_noiseSigHG
float m_noiseSigHG
Definition: ZDCPulseAnalyzer.h:121
ZDCPulseAnalyzer::m_maxADCLG
float m_maxADCLG
Definition: ZDCPulseAnalyzer.h:281
ZDCPulseAnalyzer::excludeLateLG
bool excludeLateLG() const
Definition: ZDCPulseAnalyzer.h:608
ZDCPulseAnalyzer::m_fail
bool m_fail
Definition: ZDCPulseAnalyzer.h:243
ZDCPulseAnalyzer::ExcludeLateLGBit
@ ExcludeLateLGBit
Definition: ZDCPulseAnalyzer.h:45
ZDCPulseAnalyzer::m_fitPreAmp
float m_fitPreAmp
Definition: ZDCPulseAnalyzer.h:322
ZDCPulseAnalyzer::LGUnderflow
bool LGUnderflow() const
Definition: ZDCPulseAnalyzer.h:599
ZDCPulseAnalyzer::m_refitLGChisq
float m_refitLGChisq
Definition: ZDCPulseAnalyzer.h:341
master.flag
bool flag
Definition: master.py:29
ZDCPulseAnalyzer::m_deltaTSample
float m_deltaTSample
Definition: ZDCPulseAnalyzer.h:92
ZDCPulseAnalyzer::m_nonLinCorrParamsLG
std::vector< float > m_nonLinCorrParamsLG
Definition: ZDCPulseAnalyzer.h:197
ZDCPulseAnalyzer::m_gainFactorLG
float m_gainFactorLG
Definition: ZDCPulseAnalyzer.h:117
ZDCPulseAnalyzer::m_minSampleEvt
unsigned int m_minSampleEvt
Definition: ZDCPulseAnalyzer.h:219
ZDCPulseAnalyzer::m_noiseSigLG
float m_noiseSigLG
Definition: ZDCPulseAnalyzer.h:122
chi2
double chi2(TH1 *h0, TH1 *h1)
Definition: comparitor.cxx:525
ZDCPulseAnalyzer::m_gainFactorHG
float m_gainFactorHG
Definition: ZDCPulseAnalyzer.h:116
ZDCPulseAnalyzer::m_fitFailed
bool m_fitFailed
Definition: ZDCPulseAnalyzer.h:253
ZDCPulseAnalyzer::fitFailed
bool fitFailed() const
Definition: ZDCPulseAnalyzer.h:603
ZDCPulseAnalyzer::m_maxSampleEvt
unsigned int m_maxSampleEvt
Definition: ZDCPulseAnalyzer.h:220
BindingsTest.cut
cut
This script demonstrates how to call a C++ class from Python Also how to use PyROOT is shown.
Definition: BindingsTest.py:13
ZDCPulseAnalyzer::m_timingCorrScale
float m_timingCorrScale
Definition: ZDCPulseAnalyzer.h:189
ZDCPulseAnalyzer::m_initialExpAmp
float m_initialExpAmp
Definition: ZDCPulseAnalyzer.h:292
ZDCPulseAnalyzer::m_maxADCHG
float m_maxADCHG
Definition: ZDCPulseAnalyzer.h:277
ZDCPulseAnalyzer::m_freqMHz
float m_freqMHz
Definition: ZDCPulseAnalyzer.h:91
ZDCPulseAnalyzer::m_ADCSamplesHG
std::vector< float > m_ADCSamplesHG
Definition: ZDCPulseAnalyzer.h:349
ZDCPulseAnalyzer::LGModeForceLG
@ LGModeForceLG
Definition: ZDCPulseAnalyzer.h:57
python.StandardJetMods.pull
pull
Definition: StandardJetMods.py:309
ZDCPulseAnalyzer::m_T0CutLowHG
float m_T0CutLowHG
Definition: ZDCPulseAnalyzer.h:153
ZDCPulseAnalyzer::m_fitNDoF
float m_fitNDoF
Definition: ZDCPulseAnalyzer.h:320
add-xsec-uncert-quadrature-N.results
dictionary results
Definition: add-xsec-uncert-quadrature-N.py:39
ZDCFitWrapper::GetBkgdMaxFraction
virtual float GetBkgdMaxFraction() const =0
hist_file_dump.f
f
Definition: hist_file_dump.py:140
ZDCPulseAnalyzer::GetCombinedGraph
std::shared_ptr< TGraphErrors > GetCombinedGraph(bool forceLG=false)
Definition: ZDCPulseAnalyzer.cxx:2254
ZDCPulseAnalyzer::m_havePulse
bool m_havePulse
Definition: ZDCPulseAnalyzer.h:241
python.TransformConfig.descr
descr
print "%s.properties()" % self.__name__
Definition: TransformConfig.py:359
ZDCPulseAnalyzer::AnalyzeData
bool AnalyzeData(size_t nSamples, size_t preSample, const std::vector< float > &samples, const std::vector< bool > &useSamples, float peak2ndDerivMinThresh, float noiseSig, const std::vector< float > &toCorrParams, ChisqCutLambdatype chisqCutLambda, float minT0Corr, float maxT0Corr)
Definition: ZDCPulseAnalyzer.cxx:1154
python.AtlRunQueryLib.options
options
Definition: AtlRunQueryLib.py:378
ZDCPulseAnalyzer::m_refitLGAmplCorr
float m_refitLGAmplCorr
Definition: ZDCPulseAnalyzer.h:339
ZDCPulseAnalyzer::m_nonLinCorrRefADC
float m_nonLinCorrRefADC
Definition: ZDCPulseAnalyzer.h:194
ZDCPulseAnalyzer::m_T0CutHighLG
float m_T0CutHighLG
Definition: ZDCPulseAnalyzer.h:151
ZDCPulseAnalyzer::ConfigFromJSON
std::pair< bool, std::string > ConfigFromJSON(const JSON &config)
Definition: ZDCPulseAnalyzer.cxx:2494
ZDCFitExpFermiPrePulse
Definition: ZDCFitWrapper.h:361
ZDCPulseAnalyzer::m_ampError
float m_ampError
Definition: ZDCPulseAnalyzer.h:328
ZDCFitExpFermiVariableTaus
Definition: ZDCFitWrapper.h:164
ZDCPulseAnalyzer::HGOverflowBit
@ HGOverflowBit
Definition: ZDCPulseAnalyzer.h:31
ZDCPulseAnalyzer::SetADCOverUnderflowValues
void SetADCOverUnderflowValues(int HGOverflowADC, int HGUnderflowADC, int LGOverflowADC)
Definition: ZDCPulseAnalyzer.cxx:501
ZDCPulseAnalyzer::m_ExcludeLate
bool m_ExcludeLate
Definition: ZDCPulseAnalyzer.h:258
ZDCPulseAnalyzer::m_preExclLGADCThresh
unsigned int m_preExclLGADCThresh
Definition: ZDCPulseAnalyzer.h:179
ZDCPulseAnalyzer::m_ADCPeakLG
float m_ADCPeakLG
Definition: ZDCPulseAnalyzer.h:287
ZDCPulseAnalyzer::m_refitLGFitAmpl
float m_refitLGFitAmpl
Definition: ZDCPulseAnalyzer.h:338
ZDCPulseAnalyzer::m_enablePostExcl
bool m_enablePostExcl
Definition: ZDCPulseAnalyzer.h:181
perfmonmt-printer.required
required
Definition: perfmonmt-printer.py:184
ZDCPulseAnalyzer::JSONConfigParams
static const ZDCJSONConfig::JSONParamList JSONConfigParams
Definition: ZDCPulseAnalyzer.h:71
ZDCPulseAnalyzer::obtainDelayedBaselineCorr
static float obtainDelayedBaselineCorr(const std::vector< float > &samples)
Definition: ZDCPulseAnalyzer.cxx:2378
ZDCPulseAnalyzer::m_ADCSamplesHGSub
std::vector< float > m_ADCSamplesHGSub
Definition: ZDCPulseAnalyzer.h:351
ZDCPulseAnalyzer::m_fitPreT0
float m_fitPreT0
Definition: ZDCPulseAnalyzer.h:321
ZDCPulseAnalyzer::m_useDelayed
bool m_useDelayed
Definition: ZDCPulseAnalyzer.h:108
lumiFormat.array
array
Definition: lumiFormat.py:91
ZDCMsg::Error
@ Error
Definition: ZDCMsg.h:22
ZDCPulseAnalyzer::SetTauT0Values
void SetTauT0Values(bool fixTau1, bool fixTau2, float tau1, float tau2, float t0HG, float t0LG)
Definition: ZDCPulseAnalyzer.cxx:471
ZDCPulseAnalyzer::HGUnderflow
bool HGUnderflow() const
Definition: ZDCPulseAnalyzer.h:596
ZDCPulseAnalyzer::m_chisqDivAmpPowerLG
float m_chisqDivAmpPowerLG
Definition: ZDCPulseAnalyzer.h:147
ZDCPulseAnalyzer::SetFitMinMaxAmp
void SetFitMinMaxAmp(float minAmpHG, float minAmpLG, float maxAmpHG, float maxAmpLG)
Definition: ZDCPulseAnalyzer.cxx:462
ZDCPulseAnalyzer::m_LGMode
unsigned int m_LGMode
Definition: ZDCPulseAnalyzer.h:94
ZDCPulseAnalyzer::m_HGUnderflowADC
int m_HGUnderflowADC
Definition: ZDCPulseAnalyzer.h:128
ZDCPulseAnalyzer::checkTF1Limits
void checkTF1Limits(TF1 *func)
Definition: ZDCPulseAnalyzer.cxx:2029
ZDCPulseAnalyzer::m_useLowGain
bool m_useLowGain
Definition: ZDCPulseAnalyzer.h:242
ZDCPulseAnalyzer::m_ADCSSampSigLG
std::vector< float > m_ADCSSampSigLG
Definition: ZDCPulseAnalyzer.h:358
ZDCPulseAnalyzer::Calculate2ndDerivative
static std::vector< float > Calculate2ndDerivative(const std::vector< float > &inputData, unsigned int step)
Definition: ZDCPulseAnalyzer.cxx:2346
ZDCPulseAnalyzer::m_repassPulse
bool m_repassPulse
Definition: ZDCPulseAnalyzer.h:263
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:240
createCoolChannelIdFile.par
par
Definition: createCoolChannelIdFile.py:28
ZDCPulseAnalyzer::m_ExcludeEarly
bool m_ExcludeEarly
Definition: ZDCPulseAnalyzer.h:257
ZDCPulseAnalyzer::m_preSampleIdx
unsigned int m_preSampleIdx
Definition: ZDCPulseAnalyzer.h:90
ZDCPulseAnalyzer::m_delayedDeltaT
float m_delayedDeltaT
Definition: ZDCPulseAnalyzer.h:225
ZDCFitWrapper::Initialize
void Initialize(float initialAmp, float initialT0, float ampMin, float ampMax)
Definition: ZDCFitWrapper.cxx:9
ZDCPulseAnalyzer::GetGraph
std::shared_ptr< TGraphErrors > GetGraph(bool forceLG=false)
Definition: ZDCPulseAnalyzer.cxx:2299
ZDCPulseAnalyzer::m_preExpSig
float m_preExpSig
Definition: ZDCPulseAnalyzer.h:302
ZDCPulseAnalyzer::JSON
ZDCJSONConfig::JSON JSON
Definition: ZDCPulseAnalyzer.h:26
ZDCFitWrapper::GetTau1
virtual float GetTau1() const =0
ZDCPulseAnalyzer::m_timeResFuncLG_p
std::unique_ptr< const TF1 > m_timeResFuncLG_p
Definition: ZDCPulseAnalyzer.h:157
ZDCPulseAnalyzer::ExcludeEarlyLGBit
@ ExcludeEarlyLGBit
Definition: ZDCPulseAnalyzer.h:44
ZDCPulseAnalyzer::m_minDeriv2ndSig
float m_minDeriv2ndSig
Definition: ZDCPulseAnalyzer.h:301
ZDCPulseAnalyzer::m_fitTime
float m_fitTime
Definition: ZDCPulseAnalyzer.h:312
ZDCPulseAnalyzer::m_preExpFitWrapper
std::unique_ptr< ZDCPreExpFitWrapper > m_preExpFitWrapper
Definition: ZDCPulseAnalyzer.h:213
ZDCMsg::MessageFunctionPtr
std::shared_ptr< MessageFunction > MessageFunctionPtr
Definition: ZDCMsg.h:14
ZDCPulseAnalyzer::m_peak2ndDerivMinSample
size_t m_peak2ndDerivMinSample
Definition: ZDCPulseAnalyzer.h:103
ZDCPulseAnalyzer::SetFitTimeMax
void SetFitTimeMax(float tmax)
Definition: ZDCPulseAnalyzer.cxx:489
ZDCPulseAnalyzer::m_delayedHistLGRefit
std::unique_ptr< TH1 > m_delayedHistLGRefit
Definition: ZDCPulseAnalyzer.h:228
ZDCPulseAnalyzer::failed
bool failed() const
Definition: ZDCPulseAnalyzer.h:593
ZDCFitComplexPrePulse
Definition: ZDCFitWrapper.h:963
ZDCPulseAnalyzer::m_usedPresampIdx
int m_usedPresampIdx
Definition: ZDCPulseAnalyzer.h:273
ZDCPulseAnalyzer::NoTimingCorr
@ NoTimingCorr
Definition: ZDCPulseAnalyzer.h:61
ZDCPulseAnalyzer::badT0
bool badT0() const
Definition: ZDCPulseAnalyzer.h:606
ZDCPulseAnalyzer::m_lastHGOverFlowSample
int m_lastHGOverFlowSample
Definition: ZDCPulseAnalyzer.h:345
ZDCPulseAnalyzer::m_chisqDivAmpOffsetHG
float m_chisqDivAmpOffsetHG
Definition: ZDCPulseAnalyzer.h:146
xAOD::score
@ score
Definition: TrackingPrimitives.h:514
ZDCPulseAnalyzer::m_samplesLGRefit
std::vector< float > m_samplesLGRefit
Definition: ZDCPulseAnalyzer.h:363
ReadOfcFromCool.nsamples
nsamples
Definition: ReadOfcFromCool.py:115
JSON
ZDCJSONConfig::JSON JSON
Definition: ZDCPulseAnalyzer.cxx:20
ZDCPulseAnalyzer::m_chisqDivAmpCutHG
float m_chisqDivAmpCutHG
Definition: ZDCPulseAnalyzer.h:144
ZDCPulseAnalyzer::m_peak2ndDerivMinTolerance
size_t m_peak2ndDerivMinTolerance
Definition: ZDCPulseAnalyzer.h:104
ZDCPulseAnalyzer::m_baselineCorr
float m_baselineCorr
Definition: ZDCPulseAnalyzer.h:269
ZDCPulseAnalyzer::m_PSHGOverUnderflow
bool m_PSHGOverUnderflow
Definition: ZDCPulseAnalyzer.h:247
ZDCFitExpFermiVariableTausRun3
Definition: ZDCFitWrapper.h:292
ZDCPulseAnalyzer::m_chisqDivAmpPowerHG
float m_chisqDivAmpPowerHG
Definition: ZDCPulseAnalyzer.h:148
ZDCPulseAnalyzer::DoFit
void DoFit(bool refitLG=false)
Definition: ZDCPulseAnalyzer.cxx:1578
ZDCPulseAnalyzer::m_preAmplitude
float m_preAmplitude
Definition: ZDCPulseAnalyzer.h:330
ZDCPulseAnalyzer::m_tmin
float m_tmin
Definition: ZDCPulseAnalyzer.h:95
python.CaloAddPedShiftConfig.int
int
Definition: CaloAddPedShiftConfig.py:45
python.compareTCTs.ratio
ratio
Definition: compareTCTs.py:294
ZDCPulseAnalyzer::prepareLGRefit
void prepareLGRefit(const std::vector< float > &samplesLG, const std::vector< float > &samplesSig, const std::vector< bool > &useSamples)
Definition: ZDCPulseAnalyzer.cxx:1555
ZDCPulseAnalyzer::failSigCut
bool failSigCut() const
Definition: ZDCPulseAnalyzer.h:613
ZDCPulseAnalyzer::m_ADCSSampSigHG
std::vector< float > m_ADCSSampSigHG
Definition: ZDCPulseAnalyzer.h:357
ZDCPulseAnalyzer::m_fitAmpMaxLG
float m_fitAmpMaxLG
Definition: ZDCPulseAnalyzer.h:168
ZDCPulseAnalyzer::DoAnalysis
bool DoAnalysis(bool repass)
Definition: ZDCPulseAnalyzer.cxx:1028
ZDCPulseAnalyzer::Reset
void Reset(bool reanalyze=false)
Definition: ZDCPulseAnalyzer.cxx:311
ZDCPulseAnalyzer::m_minADCSampleHG
int m_minADCSampleHG
Definition: ZDCPulseAnalyzer.h:278
ZDCPulseAnalyzer::LGModeNormal
@ LGModeNormal
Definition: ZDCPulseAnalyzer.h:56
ZDCPulseAnalyzer::m_fitAmpError
float m_fitAmpError
Definition: ZDCPulseAnalyzer.h:311
ZDCPulseAnalyzer::m_fitTimeCorr
float m_fitTimeCorr
Definition: ZDCPulseAnalyzer.h:314
ZDCPulseAnalyzer::m_samplesSigLGRefit
std::vector< float > m_samplesSigLGRefit
Definition: ZDCPulseAnalyzer.h:364
ZDCPulseAnalyzer::m_peak2ndDerivMinThreshHG
float m_peak2ndDerivMinThreshHG
Definition: ZDCPulseAnalyzer.h:106
ZDCPulseAnalyzer::m_timeSig
float m_timeSig
Definition: ZDCPulseAnalyzer.h:315
ZDCPulseAnalyzer::dump
void dump() const
Definition: ZDCPulseAnalyzer.cxx:2113
ZDCPulseAnalyzer::prePulse
bool prePulse() const
Definition: ZDCPulseAnalyzer.h:601
ZDCPulseAnalyzer::m_T0CutLowLG
float m_T0CutLowLG
Definition: ZDCPulseAnalyzer.h:150
ZDCPulseAnalyzer::setMinimumSignificance
void setMinimumSignificance(float sigMinHG, float sigMinLG)
Definition: ZDCPulseAnalyzer.cxx:449
python.CaloCondTools.log
log
Definition: CaloCondTools.py:20
ZDCPulseAnalyzer::m_defaultT0Min
float m_defaultT0Min
Definition: ZDCPulseAnalyzer.h:162
xmax
double xmax
Definition: listroot.cxx:61
convertTimingResiduals.offset
offset
Definition: convertTimingResiduals.py:71
ZDCPulseAnalyzer::BadT0Bit
@ BadT0Bit
Definition: ZDCPulseAnalyzer.h:43
LArNewCalib_DelayDump_OFC_Cali.idx
idx
Definition: LArNewCalib_DelayDump_OFC_Cali.py:69
ZDCPulseAnalyzer::m_maxSamplesPostExcl
unsigned int m_maxSamplesPostExcl
Definition: ZDCPulseAnalyzer.h:184
ZDCPulseAnalyzer::m_amplitude
float m_amplitude
Definition: ZDCPulseAnalyzer.h:326
ZDCPulseAnalyzer::m_initializedFits
bool m_initializedFits
Definition: ZDCPulseAnalyzer.h:210
ZDCPulseAnalyzer::excludeEarlyLG
bool excludeEarlyLG() const
Definition: ZDCPulseAnalyzer.h:607
ZDCPulseAnalyzer::m_defaultFitTMin
float m_defaultFitTMin
Definition: ZDCPulseAnalyzer.h:141
ZDCPulseAnalyzer::m_fitTMin
float m_fitTMin
Definition: ZDCPulseAnalyzer.h:297
ZDCPulseAnalyzer::dumpTF1
void dumpTF1(const TF1 *) const
Definition: ZDCPulseAnalyzer.cxx:2141
ZDCPulseAnalyzer::m_minDeriv2ndIndex
int m_minDeriv2ndIndex
Definition: ZDCPulseAnalyzer.h:294
ZDCPulseAnalyzer::m_LGUnderflow
bool m_LGUnderflow
Definition: ZDCPulseAnalyzer.h:249
ZDCPulseAnalyzer::m_LGOverflowADC
int m_LGOverflowADC
Definition: ZDCPulseAnalyzer.h:129
LArCellBinning.step
step
Definition: LArCellBinning.py:158
ZDCPulseAnalyzer::m_ADCSamplesLGSub
std::vector< float > m_ADCSamplesLGSub
Definition: ZDCPulseAnalyzer.h:352
ZDCPulseAnalyzer::LGModeRefitLG
@ LGModeRefitLG
Definition: ZDCPulseAnalyzer.h:58
ZDCPulseAnalyzer::m_T0CutHighHG
float m_T0CutHighHG
Definition: ZDCPulseAnalyzer.h:154
ZDCPulseAnalyzer::m_firstHGOverFlowSample
int m_firstHGOverFlowSample
Definition: ZDCPulseAnalyzer.h:346
ZDCPulseAnalyzer::GetFitPulls
std::vector< float > GetFitPulls(bool forceLG=false) const
Definition: ZDCPulseAnalyzer.cxx:585
ZDCPulseAnalyzer::s_undelayedFitHist
static TH1 * s_undelayedFitHist
Definition: ZDCPulseAnalyzer.h:78
ZDCMsg::Debug
@ Debug
Definition: ZDCMsg.h:19
ZDCPulseAnalyzer::m_fitAmpMinHG
float m_fitAmpMinHG
Definition: ZDCPulseAnalyzer.h:164
ZDCMsg::Fatal
@ Fatal
Definition: ZDCMsg.h:23
ZDCPulseAnalyzer::m_FADCCorrPerSample
bool m_FADCCorrPerSample
Definition: ZDCPulseAnalyzer.h:201
ZDCPulseAnalyzer::fitMinimumAmplitude
bool fitMinimumAmplitude() const
Definition: ZDCPulseAnalyzer.h:610
ZDCPulseAnalyzer::m_timingCorrMode
unsigned int m_timingCorrMode
Definition: ZDCPulseAnalyzer.h:187
ZDCPulseAnalyzer::m_maxSamplesPreExcl
unsigned int m_maxSamplesPreExcl
Definition: ZDCPulseAnalyzer.h:177
ZDCPulseAnalyzer::enableFADCCorrections
void enableFADCCorrections(bool correctPerSample, std::unique_ptr< const TH1 > &correHistHG, std::unique_ptr< const TH1 > &correHistLG)
Definition: ZDCPulseAnalyzer.cxx:551
LArDigits2NtupleDumper.nSamples
nSamples
Definition: LArDigits2NtupleDumper.py:85
merge.status
status
Definition: merge.py:16
ZDCPulseAnalyzer::m_FADCCorrLG
std::unique_ptr< const TH1 > m_FADCCorrLG
Definition: ZDCPulseAnalyzer.h:203
ZDCPulseAnalyzer::FitFailedBit
@ FitFailedBit
Definition: ZDCPulseAnalyzer.h:40
CaloCondBlobAlgs_fillNoiseFromASCII.tag
string tag
Definition: CaloCondBlobAlgs_fillNoiseFromASCII.py:23
ZDCPulseAnalyzer::m_initialPostPulseT0
float m_initialPostPulseT0
Definition: ZDCPulseAnalyzer.h:308
ZDCPulseAnalyzer::enableDelayed
void enableDelayed(float deltaT, float pedestalShift, bool fixedBaseline=false)
Definition: ZDCPulseAnalyzer.cxx:208
ZDCPulseAnalyzer::m_fitExpAmp
float m_fitExpAmp
Definition: ZDCPulseAnalyzer.h:325
ZDCPulseAnalyzer::m_prePulseCombinedFitter
std::unique_ptr< TFitter > m_prePulseCombinedFitter
Definition: ZDCPulseAnalyzer.h:230
ZDCPulseAnalyzer::m_fitHist
std::unique_ptr< TH1 > m_fitHist
Definition: ZDCPulseAnalyzer.h:207
ZDCPulseAnalyzer::m_timeResFuncHG_p
std::unique_ptr< const TF1 > m_timeResFuncHG_p
Definition: ZDCPulseAnalyzer.h:156
ZDCPulseAnalyzer::m_preSampleAmp
float m_preSampleAmp
Definition: ZDCPulseAnalyzer.h:329
ZDCPulseAnalyzer::m_chisqDivAmpOffsetLG
float m_chisqDivAmpOffsetLG
Definition: ZDCPulseAnalyzer.h:145
ZDCPulseAnalyzer::m_nonLinCorrParamsHG
std::vector< float > m_nonLinCorrParamsHG
Definition: ZDCPulseAnalyzer.h:196
ZDCPulseAnalyzer::getAmplitudeCorrection
double getAmplitudeCorrection(bool highGain)
Definition: ZDCPulseAnalyzer.cxx:622
ZDCPulseAnalyzer::m_nominalTau1
float m_nominalTau1
Definition: ZDCPulseAnalyzer.h:134
ZDCPulseAnalyzer::m_fitTCorr2nd
float m_fitTCorr2nd
Definition: ZDCPulseAnalyzer.h:316
pow
constexpr int pow(int base, int exp) noexcept
Definition: ap_fixedTest.cxx:15
ZDCPulseAnalyzer::s_combinedFitTMax
static float s_combinedFitTMax
Definition: ZDCPulseAnalyzer.h:81
ZDCPulseAnalyzer::PulseBit
@ PulseBit
Definition: ZDCPulseAnalyzer.h:28
ZDCPulseAnalyzer::m_postExclLGADCThresh
unsigned int m_postExclLGADCThresh
Definition: ZDCPulseAnalyzer.h:183
Amg::distance
float distance(const Amg::Vector3D &p1, const Amg::Vector3D &p2)
calculates the distance between two point in 3D space
Definition: GeoPrimitivesHelpers.h:54
ZDCPulseAnalyzer::m_refitLGAmpError
float m_refitLGAmpError
Definition: ZDCPulseAnalyzer.h:340
ZDCPulseAnalyzer::m_samplesSig
std::vector< float > m_samplesSig
Definition: ZDCPulseAnalyzer.h:361
ZDCPulseAnalyzer::m_initialPrePulseT0
float m_initialPrePulseT0
Definition: ZDCPulseAnalyzer.h:305
ZDCPulseAnalyzer::postPulse
bool postPulse() const
Definition: ZDCPulseAnalyzer.h:602
ZDCPulseAnalyzer::m_refitLGTimeSub
float m_refitLGTimeSub
Definition: ZDCPulseAnalyzer.h:343
ZDCPulseAnalyzer::preExpTailBit
@ preExpTailBit
Definition: ZDCPulseAnalyzer.h:46
ZDCPulseAnalyzer::m_fadcCorrFileName
std::string m_fadcCorrFileName
Definition: ZDCPulseAnalyzer.h:200
ZDCFitExpFermiLHCfPrePulse
Definition: ZDCFitWrapper.h:465
ZDCPulseAnalyzer::havePulse
bool havePulse() const
Definition: ZDCPulseAnalyzer.h:591
ZDCPulseAnalyzer::m_HGT0CorrParams
std::vector< float > m_HGT0CorrParams
Definition: ZDCPulseAnalyzer.h:191
ZDCPulseAnalyzer::s_pullValues
static std::vector< float > s_pullValues
Definition: ZDCPulseAnalyzer.h:83
ZDCFitExpFermiPreExp
Definition: ZDCFitWrapper.h:581
ZDCFitExpFermiLinearFixedTaus
Definition: ZDCFitWrapper.h:785
ZDCPulseAnalyzer::m_tmax
float m_tmax
Definition: ZDCPulseAnalyzer.h:96
ZDCPulseAnalyzer::m_haveSignifCuts
bool m_haveSignifCuts
Definition: ZDCPulseAnalyzer.h:170
ZDCFitGeneralPulse
Definition: ZDCFitWrapper.h:1078
ZDCPulseAnalyzer::m_refitLGTime
float m_refitLGTime
Definition: ZDCPulseAnalyzer.h:342
ZDCPulseAnalyzer::m_Nsample
unsigned int m_Nsample
Definition: ZDCPulseAnalyzer.h:89
ZDCPulseAnalyzer::CalculateDerivative
static std::vector< float > CalculateDerivative(const std::vector< float > &inputData, unsigned int step)
Definition: ZDCPulseAnalyzer.cxx:2325
python.LArMinBiasAlgConfig.float
float
Definition: LArMinBiasAlgConfig.py:65
mapkey::key
key
Definition: TElectronEfficiencyCorrectionTool.cxx:37