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