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