10 #include <top/modules/TOPTimeBaseCalibrator/TOPTimeBaseCalibratorModule.h>
11 #include <top/geometry/TOPGeometryPar.h>
14 #include <framework/datastore/StoreArray.h>
17 #include <framework/logging/Logger.h>
20 #include <top/dataobjects/TOPDigit.h>
21 #include <top/dataobjects/TOPRawDigit.h>
22 #include <top/dbobjects/TOPSampleTimes.h>
28 #include <TMatrixDSym.h>
56 setDescription(
"Sample time calibrator");
60 addParam(
"moduleID", m_moduleID,
"slot ID to calibrate.");
61 addParam(
"minTimeDiff", m_minTimeDiff,
62 "lower bound on time difference [samples].", 0.0);
63 addParam(
"maxTimeDiff", m_maxTimeDiff,
64 "upper bound on time difference [samples].", 100.0);
65 addParam(
"directoryName", m_directoryName,
66 "name (with path) of the directory for the output root files.",
68 addParam(
"minHits", m_minHits,
69 "minimal required hits per channel.", (
unsigned) 1000);
70 addParam(
"numIterations", m_numIterations,
71 "number of iTBC iterations", (
unsigned) 100);
72 addParam(
"numConvIter", m_conv_iter,
73 "Max number of iTBC iterations for conversion.", (
unsigned) 50);
74 addParam(
"minChi2Change", m_dchi2_min,
75 "Minimal of chi2 change in iterations.", 0.2);
76 addParam(
"dtMin", m_dt_min,
77 "minimum Delta T of raw calpulse in iTBC", 20.0);
78 addParam(
"dtMax", m_dt_max,
79 "maximum Delta T of raw calpulse in iTBC", 24.0);
80 addParam(
"xStep", m_xstep,
81 "unit for an interation of delta(X_s)", 0.020);
82 addParam(
"devStep", m_dev_step,
83 "a step size to calculate the value of d(chisq)/dxval", 0.001);
84 addParam(
"chgStep", m_change_xstep,
85 "update m_xstep if m_dchi2dxv < m_change_step", 0.015);
86 addParam(
"newStep", m_new_xstep,
87 "a new step for delta(X_s) if m_dchi2dxv < m_change_step", 2.0 * m_xstep);
88 addParam(
"sigm2_exp", m_sigm2_exp,
89 "(sigma_0(dT))**2 for nomarlization of chisq = sum{dT^2/sigma^2}", 0.0424 * 0.0424);
91 addParam(
"method", m_method,
"method: 0 - profile histograms only, "
92 "1 - matrix inversion, 2 - iterative, "
93 "3 - matrix inversion w/ singular value decomposition.", (
unsigned) 1);
94 addParam(
"useFallingEdge", m_useFallingEdge,
95 "if true, use cal pulse falling edge instead of rising edge",
false);
99 TOPTimeBaseCalibratorModule::~TOPTimeBaseCalibratorModule()
103 void TOPTimeBaseCalibratorModule::initialize()
107 m_digits.isRequired();
110 const auto* geo = TOPGeometryPar::Instance()->getGeometry();
111 if (!geo->isModuleIDValid(m_moduleID))
112 B2ERROR(
"Invalid module ID: " << m_moduleID);
115 if (m_directoryName.empty()) m_directoryName =
"./";
116 if (m_directoryName !=
"./") gSystem->mkdir(m_directoryName.c_str(), kTRUE);
119 m_syncTimeBase = geo->getNominalTDC().getSyncTimeBase();
122 m_goodHits = TH2F(
"goodHits",
"pulse height vs. pulse width of all good hits",
123 100, 0, 10, 100, 0, 2000);
124 m_goodHits.SetXTitle(
"pulse width (FWHM) [ns]");
125 m_goodHits.SetYTitle(
"pulse height [ADC counts]");
126 m_calPulseFirst = TH2F(
"calPulseFirst",
127 "pulse height vs. pulse width of the first calibration pulse",
128 100, 0, 10, 100, 0, 2000);
129 m_calPulseFirst.SetXTitle(
"pulse width (FWHM) [ns]");
130 m_calPulseFirst.SetYTitle(
"pulse height [ADC counts]");
131 m_calPulseSecond = TH2F(
"calPulseSecond",
132 "pulse height vs. pulse width of the second calibration pulse",
133 100, 0, 10, 100, 0, 2000);
134 m_calPulseSecond.SetXTitle(
"pulse width (FWHM) [ns]");
135 m_calPulseSecond.SetYTitle(
"pulse height [ADC counts]");
140 void TOPTimeBaseCalibratorModule::beginRun()
144 void TOPTimeBaseCalibratorModule::event()
146 const auto* geo = TOPGeometryPar::Instance()->getGeometry();
147 double sampleWidth = geo->getNominalTDC().getSampleWidth();
149 vector<Hit> hits[c_NumChannels];
151 for (
const auto& digit : m_digits) {
152 if (digit.getModuleID() != m_moduleID)
continue;
153 if (digit.getHitQuality() != TOPDigit::c_Junk) {
154 m_goodHits.Fill(digit.getPulseWidth(), digit.getPulseHeight());
156 if (digit.getHitQuality() != TOPDigit::c_CalPulse)
continue;
157 double rawTime = digit.getRawTime();
158 double errScaleFactor = 1;
159 if (m_useFallingEdge) {
160 const auto* rawDigit = digit.getRelated<
TOPRawDigit>();
162 B2ERROR(
"No relation to TOPRawDigit - can't determine falling edge time error");
167 rawTime += rawDigit->getFWHM();
168 errScaleFactor = rawDigit->getCFDFallingTimeError(1.0) / rawDigit->getCFDLeadingTimeError(1.0);
170 double t = rawTime + digit.getFirstWindow() * c_WindowSize;
172 B2ERROR(
"Got negative sample number - digit ignored");
175 double et = digit.getTimeError() / sampleWidth * errScaleFactor;
177 B2ERROR(
"Time error is not given - digit ignored");
180 unsigned channel = digit.getChannel();
181 if (channel < c_NumChannels) {
182 hits[channel].push_back(
Hit(t, et, digit.getPulseHeight(), digit.getPulseWidth()));
186 for (
unsigned channel = 0; channel < c_NumChannels; channel++) {
187 const auto& channelHits = hits[channel];
188 if (channelHits.size() == 2) {
189 double t0 = channelHits[0].time;
190 double t1 = channelHits[1].time;
191 auto diff = fabs(t0 - t1);
192 if (diff < m_minTimeDiff)
continue;
193 if (diff > m_maxTimeDiff)
continue;
194 double sig0 = channelHits[0].timeErr;
195 double sig1 = channelHits[1].timeErr;
196 double sigma = sqrt(sig0 * sig0 + sig1 * sig1);
197 m_ntuples[channel].push_back(
TwoTimes(t0, t1, sigma));
199 m_calPulseFirst.Fill(channelHits[0].pulseWidth, channelHits[0].pulseHeight);
200 m_calPulseSecond.Fill(channelHits[1].pulseWidth, channelHits[1].pulseHeight);
202 m_calPulseSecond.Fill(channelHits[0].pulseWidth, channelHits[0].pulseHeight);
203 m_calPulseFirst.Fill(channelHits[1].pulseWidth, channelHits[1].pulseHeight);
205 }
else if (channelHits.size() > 2) {
206 B2WARNING(
"More than two cal pulses per channel found - ignored");
213 void TOPTimeBaseCalibratorModule::endRun()
217 void TOPTimeBaseCalibratorModule::terminate()
220 const auto& feMapper = TOPGeometryPar::Instance()->getFrontEndMapper();
222 for (
unsigned bs = 0; bs < c_NumBoardstacks; bs ++) {
226 const auto* feMap = feMapper.getMap(m_moduleID, bs);
228 B2ERROR(
"No front-end map available for slot " << m_moduleID <<
", BS" << bs);
231 unsigned scrodID = feMap->getScrodID();
235 string fileName = m_directoryName +
"/tbcSlot";
236 if (m_moduleID < 10) {
237 fileName +=
"0" + to_string(m_moduleID);
239 fileName += to_string(m_moduleID);
241 fileName +=
"_" + to_string(bs) +
"-scrod" +
242 to_string(scrodID) +
".root";
243 TFile* fout = TFile::Open(fileName.c_str(),
"recreate");
245 B2ERROR(
"Can't open the output file " << fileName);
252 m_calPulseFirst.Write();
253 m_calPulseSecond.Write();
255 B2INFO(
"Fitting time base corrections for SCROD " << scrodID
256 <<
", output file: " << fileName);
260 TH1F Hchan(
"channels",
"channels with cal pulse",
261 c_NumScrodChannels, 0, c_NumScrodChannels);
262 Hchan.SetXTitle(
"channel number");
263 Hchan.SetYTitle(
"double cal pulse counts");
264 TH1F Hsuccess(
"success",
"successfuly fitted channels",
265 c_NumScrodChannels, 0, c_NumScrodChannels);
266 Hsuccess.SetXTitle(
"channel number");
267 TH1F Hchi2(
"chi2",
"normalized chi2",
268 c_NumScrodChannels, 0, c_NumScrodChannels);
269 Hchi2.SetXTitle(
"channel number");
270 Hchi2.SetYTitle(
"chi2/ndf");
271 TH1F Hndf(
"ndf",
"degrees of freedom",
272 c_NumScrodChannels, 0, c_NumScrodChannels);
273 Hndf.SetXTitle(
"channel number");
274 Hndf.SetYTitle(
"ndf");
275 TH1F HDeltaT(
"DeltaT",
"Fitted double pulse delay time",
276 c_NumScrodChannels, 0, c_NumScrodChannels);
277 HDeltaT.SetXTitle(
"channel number");
278 HDeltaT.SetYTitle(
"#Delta T [ns]");
282 for (
unsigned chan = 0; chan < c_NumScrodChannels; chan++) {
283 unsigned channel = chan + bs * c_NumScrodChannels;
284 auto& ntuple = m_ntuples[channel];
285 Hchan.SetBinContent(chan + 1, ntuple.size());
286 if (ntuple.size() > m_minHits) {
287 bool ok = calibrateChannel(ntuple, scrodID, chan, Hchi2, Hndf, HDeltaT);
288 if (ok) Hsuccess.Fill(chan);
290 B2INFO(
"... channel " << chan <<
" statistics too low ("
291 << ntuple.size() <<
" double cal pulses)");
308 bool TOPTimeBaseCalibratorModule::calibrateChannel(vector<TwoTimes>& ntuple,
309 unsigned scrodID,
unsigned chan,
310 TH1F& Hchi2, TH1F& Hndf,
316 string name =
"timeDiff_ch" + to_string(chan);
317 string forWhat =
"scrod " + to_string(scrodID) +
" channel " + to_string(chan);
318 string title =
"Cal pulse time difference vs. sample for " + forWhat;
319 TProfile prof(name.c_str(), title.c_str(), c_TimeAxisSize, 0, c_TimeAxisSize,
320 m_minTimeDiff, m_maxTimeDiff,
"S");
321 prof.SetXTitle(
"sample number");
322 prof.SetYTitle(
"time difference [samples]");
324 name =
"sampleOccup_ch" + to_string(chan);
325 title =
"Occupancy for " + forWhat;
326 TH1F hist(name.c_str(), title.c_str(), c_TimeAxisSize, 0, c_TimeAxisSize);
327 hist.SetXTitle(
"sample number");
328 hist.SetYTitle(
"entries per sample");
330 vector<double> profMeans(c_TimeAxisSize, (m_maxTimeDiff + m_minTimeDiff) / 2);
331 double sigma = (m_maxTimeDiff - m_minTimeDiff) / 6;
333 for (
int iter = 0; iter < 5; iter++) {
338 for (
auto& twoTimes : ntuple) {
339 int sample = int(twoTimes.t1) % c_TimeAxisSize;
340 double diff = twoTimes.t2 - twoTimes.t1;
341 double Ddiff = diff - profMeans[sample];
342 if (fabs(Ddiff) < 3 * sigma) {
343 prof.Fill(sample, diff);
348 twoTimes.good =
true;
350 twoTimes.good =
false;
353 for (
int i = 0; i < c_TimeAxisSize; i++) {
354 profMeans[i] = prof.GetBinContent(i + 1);
356 if (n == 0)
return false;
359 sigma = sqrt(x2 - x * x);
366 double meanTimeDifference = 0;
367 for (
auto& x : profMeans) meanTimeDifference += x;
368 meanTimeDifference /= profMeans.size();
369 meanTimeDifference -= int(meanTimeDifference / c_TimeAxisSize) * c_TimeAxisSize;
375 B2INFO(
"... channel " << chan <<
": "
376 << ntuple.size() <<
" double cal pulses)");
379 return matrixInversion(ntuple, scrodID, chan, meanTimeDifference,
380 Hchi2, Hndf, HDeltaT);
382 return iterativeTBC(ntuple, scrodID, chan, meanTimeDifference,
383 Hchi2, Hndf, HDeltaT);
385 B2ERROR(
"Singuler value decomposition not implemented yet");
388 B2ERROR(
"Unknown method " << m_method);
395 bool TOPTimeBaseCalibratorModule::matrixInversion(
const vector<TwoTimes>& ntuple,
396 unsigned scrodID,
unsigned chan,
397 double meanTimeDifference,
398 TH1F& Hchi2, TH1F& Hndf,
404 TMatrixDSym A(c_TimeAxisSize);
405 vector<double> b(c_TimeAxisSize, 0.0);
407 for (
const auto& twoTimes : ntuple) {
408 if (!twoTimes.good)
continue;
410 vector<double> m(c_TimeAxisSize, 0.0);
411 int i1 = int(twoTimes.t1);
412 m[i1 % c_TimeAxisSize] = 1.0 - (twoTimes.t1 - i1);
413 int i2 = int(twoTimes.t2);
414 m[i2 % c_TimeAxisSize] = twoTimes.t2 - i2;
415 i2 = i1 + (i2 - i1) % c_TimeAxisSize;
416 for (
int k = i1 + 1; k < i2; k++) m[k % c_TimeAxisSize] = 1;
418 double relSigma = twoTimes.sigma / meanTimeDifference;
419 double sig2 = relSigma * relSigma;
421 for (
int jj = i1; jj < i2 + 1; jj++) {
422 int j = jj % c_TimeAxisSize;
423 for (
int kk = i1; kk < i2 + 1; kk++) {
424 int k = kk % c_TimeAxisSize;
425 A(j, k) += m[j] * m[k] / sig2;
428 for (
int k = 0; k < c_TimeAxisSize; k++) b[k] += m[k] / sig2;
433 string forWhat =
"scrod " + to_string(scrodID) +
" channel " + to_string(chan);
434 saveAsHistogram(A,
"matA_ch" + to_string(chan),
"Matrix for " + forWhat);
435 saveAsHistogram(b,
"vecB_ch" + to_string(chan),
"Right side for " + forWhat);
442 B2INFO(
"... channel " << chan <<
" failed");
446 vector<double> x(c_TimeAxisSize, 0.0);
447 for (
int k = 0; k < c_TimeAxisSize; k++) {
448 for (
int j = 0; j < c_TimeAxisSize; j++) {
449 x[k] += A(k, j) * b[j];
456 int ndf = -c_TimeAxisSize;
458 for (
const auto& twoTimes : ntuple) {
459 if (!twoTimes.good)
continue;
461 vector<double> m(c_TimeAxisSize, 0.0);
462 int i1 = int(twoTimes.t1);
463 m[i1 % c_TimeAxisSize] = 1.0 - (twoTimes.t1 - i1);
464 int i2 = int(twoTimes.t2);
465 m[i2 % c_TimeAxisSize] = twoTimes.t2 - i2;
466 i2 = i1 + (i2 - i1) % c_TimeAxisSize;
467 for (
int k = i1 + 1; k < i2; k++) m[k % c_TimeAxisSize] = 1;
469 for (
int k = 0; k < c_TimeAxisSize; k++) s += m[k] * x[k];
470 double relSigma = twoTimes.sigma / meanTimeDifference;
471 double sig2 = relSigma * relSigma;
472 chi2 += s * s / sig2;
475 Hchi2.SetBinContent(chan + 1, chi2 / ndf);
476 Hndf.SetBinContent(chan + 1, ndf);
481 for (
auto xi : x) sum += xi;
486 double DeltaT = 2 * m_syncTimeBase / sum;
487 for (
auto& xi : x) xi *= DeltaT;
488 HDeltaT.SetBinContent(chan + 1, DeltaT);
491 for (
int k = 0; k < c_TimeAxisSize; k++) err.push_back(sqrt(A(k, k)) * DeltaT);
493 vector<double> sampleTimes;
494 sampleTimes.push_back(0);
495 for (
auto xi : x) sampleTimes.push_back(xi + sampleTimes.back());
499 saveAsHistogram(A,
"invA_ch" + to_string(chan),
"Inverted matrix for " + forWhat);
500 saveAsHistogram(x, err,
"dt_ch" + to_string(chan),
"Sample time bins for " + forWhat,
501 "sample number",
"#Delta t [ns]");
502 saveAsHistogram(sampleTimes,
"sampleTimes_ch" + to_string(chan),
503 "Time base corrections for " + forWhat,
"sample number",
"t [ns]");
507 string name =
"timeDiffcal_ch" + to_string(chan);
508 string title =
"Calibrated cal pulse time difference vs. sample for " + forWhat;
509 TH2F Hcor(name.c_str(), title.c_str(), c_TimeAxisSize, 0, c_TimeAxisSize,
510 100, DeltaT - 0.5, DeltaT + 0.5);
511 Hcor.SetXTitle(
"sample number");
512 Hcor.SetYTitle(
"time difference [ns]");
513 Hcor.SetStats(kTRUE);
516 timeBase.
setTimeAxis(sampleTimes, sampleTimes.back() / 2);
518 for (
const auto& twoTimes : ntuple) {
519 if (!twoTimes.good)
continue;
520 double dt = timeBase.
getDeltaTime(0, twoTimes.t2, twoTimes.t1);
521 int sample = int(twoTimes.t1) % c_TimeAxisSize;
522 Hcor.Fill(sample, dt);
526 B2INFO(
"... channel " << chan <<
" OK (chi^2/ndf = " << chi2 / ndf
527 <<
", ndf = " << ndf <<
")");
533 bool TOPTimeBaseCalibratorModule::iterativeTBC(
const std::vector<TwoTimes>& ntuple,
534 unsigned scrodID,
unsigned chan,
535 double meanTimeDifference,
536 TH1F& Hchi2, TH1F& Hndf,
539 std::vector<double> xval(c_TimeAxisSize + 1, 0.0);
540 double wx = 2 * m_syncTimeBase / c_TimeAxisSize;
541 for (
int i = 0; i < c_TimeAxisSize + 1; i++) xval[i] = i * wx;
543 B2INFO(
"TimeBaseCalibration starts for channel#" << chan);
545 double pre_chi2 = 10000000.0;
546 unsigned num_small_dev = 0;
548 for (
unsigned j = 0; j < m_numIterations; j++) {
550 Iteration(ntuple, xval);
551 double this_chi2 = Chisq(ntuple, xval);
552 if (this_chi2 < 0)
continue;
553 double deltaChi2 = pre_chi2 - this_chi2;
554 if (deltaChi2 < -m_dchi2_min)
break;
555 if (fabs(deltaChi2) < m_deltamin) num_small_dev++;
556 if (num_small_dev > m_conv_iter)
break;
557 pre_chi2 = this_chi2;
562 double chi2 = Chisq(ntuple, xval);
563 Hchi2.SetBinContent(chan + 1, chi2);
564 Hndf.SetBinContent(chan + 1, m_good);
569 for (
auto xi : xval) sum += xi;
575 double DeltaT = meanTimeDifference * (2 * m_syncTimeBase / c_TimeAxisSize);
576 HDeltaT.SetBinContent(chan + 1, DeltaT);
578 std::vector<double> timeInterval;
579 for (
int i = 0; i < c_TimeAxisSize; i++)timeInterval.push_back(xval[i + 1] - xval[i]);
582 std::vector<double> sampleTimes;
583 for (
auto xi : xval) sampleTimes.push_back(xi);
586 std::string forWhat =
"scrod " + to_string(scrodID) +
" channel " + to_string(chan);
587 saveAsHistogram(timeInterval,
"dt_ch" + to_string(chan),
"Sample time bins for " + forWhat,
588 "sample number",
"#Delta t [ns]");
589 saveAsHistogram(sampleTimes,
"sampleTimes_ch" + to_string(chan),
590 "Time base corrections for " + forWhat,
"sample number",
"t [ns]");
593 std::string name =
"timeDiffcal_ch" + to_string(chan);
594 std::string title =
"Calibrated cal pulse time difference vs. sample for " + forWhat;
595 TH2F Hcor(name.c_str(), title.c_str(), c_TimeAxisSize, 0, c_TimeAxisSize,
596 100, DeltaT - 0.5, DeltaT + 0.5);
597 Hcor.SetXTitle(
"sample number");
598 Hcor.SetYTitle(
"time difference [ns]");
599 Hcor.SetStats(kTRUE);
602 timeBase.
setTimeAxis(sampleTimes, sampleTimes.back() / 2);
604 for (
const auto& twoTimes : ntuple) {
605 if (!twoTimes.good)
continue;
606 double dt = timeBase.
getDeltaTime(0, twoTimes.t2, twoTimes.t1);
607 int sample = int(twoTimes.t1) % c_TimeAxisSize;
608 Hcor.Fill(sample, dt);
612 B2INFO(
"... channel " << chan <<
" OK (chi^2/ndf = " << chi2
613 <<
", ndf = " << m_good <<
")");
618 void TOPTimeBaseCalibratorModule::Iteration(
const std::vector<TwoTimes>& ntuple, std::vector<double>& xval)
620 for (
int i = 0; i < c_TimeAxisSize; i++) {
621 double wdth = xval[i + 1] - xval[i];
622 if (wdth < m_min_binwidth) {
623 xval[i] = xval[i] - 0.5 * fabs(wdth) - 0.5 * m_min_binwidth;
624 xval[i + 1] = xval[i + 1] + 0.5 * fabs(wdth) + 0.5 * m_min_binwidth;
626 if (wdth > m_max_binwidth) {
627 xval[i] = xval[i] - 0.5 * fabs(wdth) - 0.5 * m_max_binwidth;
628 xval[i + 1] = xval[i + 1] + 0.5 * fabs(wdth) + 0.5 * m_max_binwidth;
633 for (
int i = 0; i < c_TimeAxisSize; i++) xval[i] = xval[i] - xval[0];
635 std::vector<double> xxval(c_TimeAxisSize + 1, 0.0);
636 for (
int i = 0; i < c_TimeAxisSize + 1; i++) xxval[i] = xval[i];
638 double chi2_0 = Chisq(ntuple, xxval);
639 if (chi2_0 < 0) B2ERROR(
"iTBC chisq_0<0! xval has problem.");
641 std::vector<double> dr_chi2(c_TimeAxisSize + 1, 0.0);
642 TH1D hdrsamp_try(
"hdrsamp_try",
"dchi2/dx distribution", 100, -0.01, 0.01);
644 for (
int smp = 1; smp < c_TimeAxisSize; smp++) {
645 xxval[smp] = xval[smp] + m_dev_step;
646 double chi2_ch = Chisq(ntuple, xxval);
647 if (chi2_ch < 0)
continue;
648 dr_chi2[smp] = (chi2_ch - chi2_0) / m_dev_step;
649 hdrsamp_try.Fill(dr_chi2[smp]);
650 xxval[smp] = xval[smp];
653 for (
int smp = 1; smp < c_TimeAxisSize; smp++) {
654 double vx_it_step = dr_chi2[smp] * m_xstep;
655 xval[smp] = xval[smp] - vx_it_step;
659 m_dchi2dxv = hdrsamp_try.GetRMS();
661 if (fabs(m_dchi2dxv) < m_change_xstep) m_xstep = m_new_xstep;
664 double TOPTimeBaseCalibratorModule::Chisq(
const std::vector<TwoTimes>& ntuple,
const std::vector<double>& xxval)
671 for (
const auto& twoTimes : ntuple) {
672 if (!twoTimes.good)
continue;
674 std::vector<double> m(c_TimeAxisSize, 0.0);
676 int i1 = int(twoTimes.t1);
677 double fr = twoTimes.t1 - i1;
678 double samp0 = i1 % 256;
679 double ctdc1 = xxval[samp0] + fr * (xxval[samp0 + 1] - xxval[samp0]);
680 int i2 = int(twoTimes.t2);
681 fr = twoTimes.t2 - i2;
682 double samp1 = i2 % 256;
683 double ctdc2 = xxval[samp1] + fr * (xxval[samp1 + 1] - xxval[samp1]);
685 if (samp1 > samp0) cdt = ctdc2 - ctdc1;
686 else cdt = ctdc2 - ctdc1 + m_syncTimeBase * 2;
688 if (cdt < m_dt_max && cdt > m_dt_min) {
698 double mean = sum1 / m_good;
699 chi2 = (sum2 - m_good * mean * mean) / m_good / m_sigm2_exp;
706 void TOPTimeBaseCalibratorModule::saveAsHistogram(
const std::vector<double>& vec,
707 const std::string& name,
708 const std::string& title,
709 const std::string& xTitle,
710 const std::string& yTitle)
const
712 if (vec.empty())
return;
714 TH1F h(name.c_str(), title.c_str(), vec.size(), 0, vec.size());
715 h.SetXTitle(xTitle.c_str());
716 h.SetYTitle(yTitle.c_str());
717 if (name.find(
"Fit") != string::npos) h.SetLineColor(2);
719 for (
unsigned i = 0; i < vec.size(); i++) h.SetBinContent(i + 1, vec[i]);
725 void TOPTimeBaseCalibratorModule::saveAsHistogram(
const vector<double>& vec,
726 const vector<double>& err,
729 const string& xTitle,
730 const string& yTitle)
const
732 if (vec.empty())
return;
734 TH1F h(name.c_str(), title.c_str(), vec.size(), 0, vec.size());
735 h.SetXTitle(xTitle.c_str());
736 h.SetYTitle(yTitle.c_str());
738 for (
unsigned i = 0; i < vec.size(); i++) h.SetBinContent(i + 1, vec[i]);
739 for (
unsigned i = 0; i < err.size(); i++) h.SetBinError(i + 1, err[i]);
745 void TOPTimeBaseCalibratorModule::saveAsHistogram(
const TMatrixDSym& M,
747 const string& title)
const
749 int n = M.GetNrows();
750 TH2F h(name.c_str(), title.c_str(), n, 0, n, n, 0, n);
751 h.SetXTitle(
"columns");
754 for (
int j = 0; j < n; j++) {
755 for (
int k = 0; k < n; k++) {
756 h.SetBinContent(j + 1, n - k, M(j, k));
Class to store unpacked raw data (hits in feature-extraction format) It provides also calculation of ...
Calibration constants of a singe ASIC channel: time axis (sample times)
double getDeltaTime(int window, double sample2, double sample1) const
Returns time difference between sample2 and sample1.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
void setTimeAxis(double syncTimeBase)
Sets equidistant time axis (uncalibrated).
Abstract base class for different kinds of events.
Structure to hold some of the calpulse data.
Structure to hold calpulse raw times expressed in samples since sample 0 of window 0.