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