Belle II Software development
TOPTimeBaseCalibratorModule.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9// Own header.
10#include <top/modules/TOPTimeBaseCalibrator/TOPTimeBaseCalibratorModule.h>
11
12// TOP headers.
13#include <top/geometry/TOPGeometryPar.h>
14
15// framework - DataStore
16#include <framework/datastore/StoreArray.h>
17
18// framework aux
19#include <framework/logging/Logger.h>
20
21// DataStore classes
22#include <top/dataobjects/TOPDigit.h>
23#include <top/dataobjects/TOPRawDigit.h>
24#include <top/dbobjects/TOPSampleTimes.h>
25
26// Root
27#include <TFile.h>
28#include <TProfile.h>
29#include <TSystem.h>
30#include <TMatrixDSym.h>
31
32#include <iostream>
33
34using namespace std;
35
36namespace Belle2 {
41
42 using namespace TOP;
43
44 //-----------------------------------------------------------------
46 //-----------------------------------------------------------------
47
48 REG_MODULE(TOPTimeBaseCalibrator);
49
50 //-----------------------------------------------------------------
51 // Implementation
52 //-----------------------------------------------------------------
53
55
56 {
57 // set module description
58 setDescription("Sample time calibrator");
59 // setPropertyFlags(c_ParallelProcessingCertified | c_InitializeInProcess);
60
61 // Add parameters
62 addParam("moduleID", m_moduleID, "slot ID to calibrate.");
63 addParam("minTimeDiff", m_minTimeDiff,
64 "lower bound on time difference [samples].", 0.0);
65 addParam("maxTimeDiff", m_maxTimeDiff,
66 "upper bound on time difference [samples].", 100.0);
67 addParam("directoryName", m_directoryName,
68 "name (with path) of the directory for the output root files.",
69 string(""));
70 addParam("minHits", m_minHits,
71 "minimal required hits per channel.", (unsigned) 1000);
72 addParam("numIterations", m_numIterations,
73 "number of iTBC iterations", (unsigned) 100);
74 addParam("numConvIter", m_conv_iter,
75 "Max number of iTBC iterations for conversion.", (unsigned) 50);
76 addParam("minChi2Change", m_dchi2_min,
77 "Minimal of chi2 change in iterations.", 0.2);
78 addParam("dtMin", m_dt_min,
79 "minimum Delta T of raw calpulse in iTBC", 20.0);
80 addParam("dtMax", m_dt_max,
81 "maximum Delta T of raw calpulse in iTBC", 24.0);
82 addParam("xStep", m_xstep,
83 "unit for an interaction of delta(X_s)", 0.020);
84 addParam("devStep", m_dev_step,
85 "a step size to calculate the value of d(chisq)/dxval", 0.001);
86 addParam("chgStep", m_change_xstep,
87 "update m_xstep if m_dchi2dxv < m_change_step", 0.015);
88 addParam("newStep", m_new_xstep,
89 "a new step for delta(X_s) if m_dchi2dxv < m_change_step", 2.0 * m_xstep);
90 addParam("sigm2_exp", m_sigm2_exp,
91 "(sigma_0(dT))**2 for nomarlization of chisq = sum{dT^2/sigma^2}", 0.0424 * 0.0424);
92
93 addParam("method", m_method, "method: 0 - profile histograms only, "
94 "1 - matrix inversion, 2 - iterative, "
95 "3 - matrix inversion w/ singular value decomposition.", (unsigned) 1);
96 addParam("useFallingEdge", m_useFallingEdge,
97 "if true, use cal pulse falling edge instead of rising edge", false);
98 addParam("saveMatrix", m_saveMatrix, "if true, save also matrix and its inverse in a root file", false);
99
100 }
101
103 {
104
105 // input
106 m_digits.isRequired();
107
108 // checks
109 const auto* geo = TOPGeometryPar::Instance()->getGeometry();
110 if (!geo->isModuleIDValid(m_moduleID))
111 B2ERROR("Invalid module ID: " << m_moduleID);
112
113 // check for existence and mkdir if not
114 if (m_directoryName.empty()) m_directoryName = "./";
115 if (m_directoryName != "./") gSystem->mkdir(m_directoryName.c_str(), kTRUE);
116
117 // synchronization time corresponding to two ASIC windows
118 m_syncTimeBase = geo->getNominalTDC().getSyncTimeBase();
119
120 // control histograms
121 m_goodHits = TH2F("goodHits", "pulse height vs. pulse width of all good hits",
122 100, 0, 10, 100, 0, 2000);
123 m_goodHits.SetXTitle("pulse width (FWHM) [ns]");
124 m_goodHits.SetYTitle("pulse height [ADC counts]");
125 m_calPulseFirst = TH2F("calPulseFirst",
126 "pulse height vs. pulse width of the first calibration pulse",
127 100, 0, 10, 100, 0, 2000);
128 m_calPulseFirst.SetXTitle("pulse width (FWHM) [ns]");
129 m_calPulseFirst.SetYTitle("pulse height [ADC counts]");
130 m_calPulseSecond = TH2F("calPulseSecond",
131 "pulse height vs. pulse width of the second calibration pulse",
132 100, 0, 10, 100, 0, 2000);
133 m_calPulseSecond.SetXTitle("pulse width (FWHM) [ns]");
134 m_calPulseSecond.SetYTitle("pulse height [ADC counts]");
135
136 }
137
138
140 {
141 const auto* geo = TOPGeometryPar::Instance()->getGeometry();
142 double sampleWidth = geo->getNominalTDC().getSampleWidth();
143
144 vector<Hit> hits[c_NumChannels];
145
146 for (const auto& digit : m_digits) {
147 if (digit.getModuleID() != m_moduleID) continue;
148 if (digit.getHitQuality() != TOPDigit::c_Junk) {
149 m_goodHits.Fill(digit.getPulseWidth(), digit.getPulseHeight());
150 }
151 if (digit.getHitQuality() != TOPDigit::c_CalPulse) continue;
152 double rawTime = digit.getRawTime();
153 double errScaleFactor = 1;
154 if (m_useFallingEdge) {
155 const auto* rawDigit = digit.getRelated<TOPRawDigit>();
156 if (!rawDigit) {
157 B2ERROR("No relation to TOPRawDigit - can't determine falling edge time error");
158 continue;
159 }
160 // rawTime may include corrections due to window number discontinuity,
161 // therefore one must add the width and not just use getCFDFallingTime()
162 rawTime += rawDigit->getFWHM();
163 errScaleFactor = rawDigit->getCFDFallingTimeError(1.0) / rawDigit->getCFDLeadingTimeError(1.0);
164 }
165 double t = rawTime + digit.getFirstWindow() * c_WindowSize;
166 if (t < 0) {
167 B2ERROR("Got negative sample number - digit ignored");
168 continue;
169 }
170 double et = digit.getTimeError() / sampleWidth * errScaleFactor;
171 if (et <= 0) {
172 B2ERROR("Time error is not given - digit ignored");
173 continue;
174 }
175 unsigned channel = digit.getChannel();
176 if (channel < c_NumChannels) {
177 hits[channel].push_back(Hit(t, et, digit.getPulseHeight(), digit.getPulseWidth()));
178 }
179 }
180
181 for (unsigned channel = 0; channel < c_NumChannels; channel++) {
182 const auto& channelHits = hits[channel];
183 if (channelHits.size() == 2) {
184 double t0 = channelHits[0].time;
185 double t1 = channelHits[1].time;
186 auto diff = fabs(t0 - t1); // since not sorted yet
187 if (diff < m_minTimeDiff) continue;
188 if (diff > m_maxTimeDiff) continue;
189 double sig0 = channelHits[0].timeErr;
190 double sig1 = channelHits[1].timeErr;
191 double sigma = sqrt(sig0 * sig0 + sig1 * sig1);
192 m_ntuples[channel].push_back(TwoTimes(t0, t1, sigma));
193 if (t0 < t1) { // check, since not sorted yet
194 m_calPulseFirst.Fill(channelHits[0].pulseWidth, channelHits[0].pulseHeight);
195 m_calPulseSecond.Fill(channelHits[1].pulseWidth, channelHits[1].pulseHeight);
196 } else {
197 m_calPulseSecond.Fill(channelHits[0].pulseWidth, channelHits[0].pulseHeight);
198 m_calPulseFirst.Fill(channelHits[1].pulseWidth, channelHits[1].pulseHeight);
199 }
200 } else if (channelHits.size() > 2) {
201 B2WARNING("More than two cal pulses per channel found - ignored");
202 }
203 }
204
205 }
206
207
209 {
210
211 const auto& feMapper = TOPGeometryPar::Instance()->getFrontEndMapper();
212
213 for (unsigned bs = 0; bs < c_NumBoardstacks; bs ++) {
214
215 // determine scrod ID
216
217 const auto* feMap = feMapper.getMap(m_moduleID, bs);
218 if (!feMap) {
219 B2ERROR("No front-end map available for slot " << m_moduleID << ", BS" << bs);
220 continue;
221 }
222 unsigned scrodID = feMap->getScrodID();
223
224 // open output root file
225
226 string fileName = m_directoryName + "/tbcSlot";
227 if (m_moduleID < 10) {
228 fileName += "0" + to_string(m_moduleID);
229 } else {
230 fileName += to_string(m_moduleID);
231 }
232 fileName += "_" + to_string(bs) + "-scrod" +
233 to_string(scrodID) + ".root";
234 TFile* fout = TFile::Open(fileName.c_str(), "recreate");
235 if (!fout) {
236 B2ERROR("Can't open the output file " << fileName);
237 continue;
238 }
239
240 // write control histograms
241
242 m_goodHits.Write();
243 m_calPulseFirst.Write();
244 m_calPulseSecond.Write();
245
246 B2INFO("Fitting time base corrections for SCROD " << scrodID
247 << ", output file: " << fileName);
248
249 // book histograms
250
251 TH1F Hchan("channels", "channels with cal pulse",
252 c_NumScrodChannels, 0, c_NumScrodChannels);
253 Hchan.SetXTitle("channel number");
254 Hchan.SetYTitle("double cal pulse counts");
255 TH1F Hsuccess("success", "successfully fitted channels",
256 c_NumScrodChannels, 0, c_NumScrodChannels);
257 Hsuccess.SetXTitle("channel number");
258 TH1F Hchi2("chi2", "normalized chi2",
259 c_NumScrodChannels, 0, c_NumScrodChannels);
260 Hchi2.SetXTitle("channel number");
261 Hchi2.SetYTitle("chi2/ndf");
262 TH1F Hndf("ndf", "degrees of freedom",
263 c_NumScrodChannels, 0, c_NumScrodChannels);
264 Hndf.SetXTitle("channel number");
265 Hndf.SetYTitle("ndf");
266 TH1F HDeltaT("DeltaT", "Fitted double pulse delay time",
267 c_NumScrodChannels, 0, c_NumScrodChannels);
268 HDeltaT.SetXTitle("channel number");
269 HDeltaT.SetYTitle("#Delta T [ns]");
270
271 // loop over channels of a single SCROD
272
273 for (unsigned chan = 0; chan < c_NumScrodChannels; chan++) {
274 unsigned channel = chan + bs * c_NumScrodChannels;
275 auto& ntuple = m_ntuples[channel];
276 Hchan.SetBinContent(chan + 1, ntuple.size());
277 if (ntuple.size() > m_minHits) {
278 bool ok = calibrateChannel(ntuple, scrodID, chan, Hchi2, Hndf, HDeltaT);
279 if (ok) Hsuccess.Fill(chan);
280 } else {
281 B2INFO("... channel " << chan << " statistics too low ("
282 << ntuple.size() << " double cal pulses)");
283 }
284 }
285
286 // write histograms and close the file
287
288 Hchan.Write();
289 Hsuccess.Write();
290 Hchi2.Write();
291 Hndf.Write();
292 HDeltaT.Write();
293 fout->Close();
294
295 }
296 }
297
298
300 unsigned scrodID, unsigned chan,
301 TH1F& Hchi2, TH1F& Hndf,
302 TH1F& HDeltaT)
303 {
304
305 // determine outlayer removal cuts
306
307 string name = "timeDiff_ch" + to_string(chan);
308 string forWhat = "scrod " + to_string(scrodID) + " channel " + to_string(chan);
309 string title = "Cal pulse time difference vs. sample for " + forWhat;
310 TProfile prof(name.c_str(), title.c_str(), c_TimeAxisSize, 0, c_TimeAxisSize,
312 prof.SetXTitle("sample number");
313 prof.SetYTitle("time difference [samples]");
314
315 name = "sampleOccup_ch" + to_string(chan);
316 title = "Occupancy for " + forWhat;
317 TH1F hist(name.c_str(), title.c_str(), c_TimeAxisSize, 0, c_TimeAxisSize);
318 hist.SetXTitle("sample number");
319 hist.SetYTitle("entries per sample");
320
321 vector<double> profMeans(c_TimeAxisSize, (m_maxTimeDiff + m_minTimeDiff) / 2);
322 double sigma = (m_maxTimeDiff - m_minTimeDiff) / 6;
323
324 for (int iter = 0; iter < 5; iter++) {
325 prof.Reset();
326 double x = 0;
327 double x2 = 0;
328 int n = 0;
329 for (auto& twoTimes : ntuple) {
330 int sample = int(twoTimes.t1) % c_TimeAxisSize;
331 double diff = twoTimes.t2 - twoTimes.t1;
332 double Ddiff = diff - profMeans[sample];
333 if (fabs(Ddiff) < 3 * sigma) {
334 prof.Fill(sample, diff);
335 hist.Fill(sample);
336 x += Ddiff;
337 x2 += Ddiff * Ddiff;
338 n++;
339 twoTimes.good = true;
340 } else {
341 twoTimes.good = false;
342 }
343 }
344 for (int i = 0; i < c_TimeAxisSize; i++) {
345 profMeans[i] = prof.GetBinContent(i + 1);
346 }
347 if (n == 0) return false;
348 x2 /= n;
349 x /= n;
350 sigma = sqrt(x2 - x * x);
351 }
352 prof.Write();
353 hist.Write();
354
355 // calculate average time difference
356
357 double meanTimeDifference = 0;
358 for (const auto& x : profMeans) meanTimeDifference += x;
359 meanTimeDifference /= profMeans.size();
360 meanTimeDifference -= int(meanTimeDifference / static_cast<double>(c_TimeAxisSize)) * c_TimeAxisSize;
361
362 // perform the fit
363
364 switch (m_method) {
365 case 0:
366 B2INFO("... channel " << chan << ": "
367 << ntuple.size() << " double cal pulses)");
368 return false;
369 case 1:
370 return matrixInversion(ntuple, scrodID, chan, meanTimeDifference,
371 Hchi2, Hndf, HDeltaT);
372 case 2:
373 return iterativeTBC(ntuple, scrodID, chan, meanTimeDifference,
374 Hchi2, Hndf, HDeltaT);
375 case 3:
376 B2ERROR("Singuler value decomposition not implemented yet");
377 return false;
378 default:
379 B2ERROR("Unknown method " << m_method);
380 return false;
381 }
382
383 }
384
385
386 bool TOPTimeBaseCalibratorModule::matrixInversion(const vector<TwoTimes>& ntuple,
387 unsigned scrodID, unsigned chan,
388 double meanTimeDifference,
389 TH1F& Hchi2, TH1F& Hndf,
390 TH1F& HDeltaT)
391 {
392
393 // Ax = b: construct matrix A and right side vector b
394
395 TMatrixDSym A(c_TimeAxisSize);
396 vector<double> b(c_TimeAxisSize, 0.0);
397
398 for (const auto& twoTimes : ntuple) {
399 if (!twoTimes.good) continue;
400
401 vector<double> m(c_TimeAxisSize, 0.0);
402 int i1 = int(twoTimes.t1);
403 m[i1 % c_TimeAxisSize] = 1.0 - (twoTimes.t1 - i1);
404 int i2 = int(twoTimes.t2);
405 m[i2 % c_TimeAxisSize] = twoTimes.t2 - i2;
406 i2 = i1 + (i2 - i1) % c_TimeAxisSize;
407 for (int k = i1 + 1; k < i2; k++) m[k % c_TimeAxisSize] = 1;
408
409 double relSigma = twoTimes.sigma / meanTimeDifference;
410 double sig2 = relSigma * relSigma;
411
412 for (int jj = i1; jj < i2 + 1; jj++) {
413 int j = jj % c_TimeAxisSize;
414 for (int kk = i1; kk < i2 + 1; kk++) {
415 int k = kk % c_TimeAxisSize;
416 A(j, k) += m[j] * m[k] / sig2;
417 }
418 }
419 for (int k = 0; k < c_TimeAxisSize; k++) b[k] += m[k] / sig2;
420 }
421
422 // save as histograms
423
424 string forWhat = "scrod " + to_string(scrodID) + " channel " + to_string(chan);
425 if (m_saveMatrix) {
426 saveAsHistogram(A, "matA_ch" + to_string(chan), "Matrix for " + forWhat);
427 saveAsHistogram(b, "vecB_ch" + to_string(chan), "Right side for " + forWhat);
428 }
429
430 // invert matrix A and solve the equation: x = A^{-1}b
431
432 double det = 0;
433 A.Invert(&det);
434 if (det == 0) {
435 B2INFO("... channel " << chan << " failed");
436 return false;
437 }
438
439 vector<double> x(c_TimeAxisSize, 0.0);
440 for (int k = 0; k < c_TimeAxisSize; k++) {
441 for (int j = 0; j < c_TimeAxisSize; j++) {
442 x[k] += A(k, j) * b[j];
443 }
444 }
445
446 // calculate chi^2
447
448 double chi2 = 0;
449 int ndf = -c_TimeAxisSize;
450
451 for (const auto& twoTimes : ntuple) {
452 if (!twoTimes.good) continue;
453
454 vector<double> m(c_TimeAxisSize, 0.0);
455 int i1 = int(twoTimes.t1);
456 m[i1 % c_TimeAxisSize] = 1.0 - (twoTimes.t1 - i1);
457 int i2 = int(twoTimes.t2);
458 m[i2 % c_TimeAxisSize] = twoTimes.t2 - i2;
459 i2 = i1 + (i2 - i1) % c_TimeAxisSize;
460 for (int k = i1 + 1; k < i2; k++) m[k % c_TimeAxisSize] = 1;
461 double s = -1.0;
462 for (int k = 0; k < c_TimeAxisSize; k++) s += m[k] * x[k];
463 double relSigma = twoTimes.sigma / meanTimeDifference;
464 double sig2 = relSigma * relSigma;
465 chi2 += s * s / sig2;
466 ndf++;
467 }
468 Hchi2.SetBinContent(chan + 1, chi2 / ndf);
469 Hndf.SetBinContent(chan + 1, ndf);
470
471 // constrain sum of x to 2*syncTimeBase and calculate sample times
472
473 double sum = 0;
474 for (auto xi : x) sum += xi;
475 if (sum == 0) {
476 B2ERROR("sum == 0");
477 return false;
478 }
479 double DeltaT = 2 * m_syncTimeBase / sum;
480 for (auto& xi : x) xi *= DeltaT;
481 HDeltaT.SetBinContent(chan + 1, DeltaT);
482
483 vector<double> err;
484 for (int k = 0; k < c_TimeAxisSize; k++) err.push_back(sqrt(A(k, k)) * DeltaT);
485
486 vector<double> sampleTimes;
487 sampleTimes.push_back(0);
488 for (auto xi : x) sampleTimes.push_back(xi + sampleTimes.back());
489
490 // save results as histograms
491
492 if (m_saveMatrix) saveAsHistogram(A, "invA_ch" + to_string(chan), "Inverted matrix for " + forWhat);
493 saveAsHistogram(x, err, "dt_ch" + to_string(chan), "Sample time bins for " + forWhat,
494 "sample number", "#Delta t [ns]");
495 saveAsHistogram(sampleTimes, "sampleTimes_ch" + to_string(chan),
496 "Time base corrections for " + forWhat, "sample number", "t [ns]");
497
498 // calibrated cal pulse time difference
499
500 string name = "timeDiffcal_ch" + to_string(chan);
501 string title = "Calibrated cal pulse time difference vs. sample for " + forWhat;
502 TH2F Hcor(name.c_str(), title.c_str(), c_TimeAxisSize, 0, c_TimeAxisSize,
503 100, DeltaT - 0.5, DeltaT + 0.5);
504 Hcor.SetXTitle("sample number");
505 Hcor.SetYTitle("time difference [ns]");
506 Hcor.SetStats(kTRUE);
507
508 TOPSampleTimes timeBase;
509 timeBase.setTimeAxis(sampleTimes, sampleTimes.back() / 2);
510
511 for (const auto& twoTimes : ntuple) {
512 if (!twoTimes.good) continue;
513 double dt = timeBase.getDeltaTime(0, twoTimes.t2, twoTimes.t1);
514 int sample = int(twoTimes.t1) % c_TimeAxisSize;
515 Hcor.Fill(sample, dt);
516 }
517 Hcor.Write();
518
519 B2INFO("... channel " << chan << " OK (chi^2/ndf = " << chi2 / ndf
520 << ", ndf = " << ndf << ")");
521
522 return true;
523 }
524
525
526 bool TOPTimeBaseCalibratorModule::iterativeTBC(const std::vector<TwoTimes>& ntuple,
527 unsigned scrodID, unsigned chan,
528 double meanTimeDifference,
529 TH1F& Hchi2, TH1F& Hndf,
530 TH1F& HDeltaT)
531 {
532 std::vector<double> xval(c_TimeAxisSize + 1, 0.0);
533 double wx = 2 * m_syncTimeBase / static_cast<double>(c_TimeAxisSize);
534 for (int i = 0; i < c_TimeAxisSize + 1; i++) xval[i] = i * wx;
535
536 B2INFO("TimeBaseCalibration starts for channel#" << chan);
537
538 double pre_chi2 = 10000000.0;
539 unsigned num_small_dev = 0;
540
541 for (unsigned j = 0; j < m_numIterations; j++) {
542
543 Iteration(ntuple, xval);
544 double this_chi2 = Chisq(ntuple, xval);
545 if (this_chi2 < 0)continue;
546 double deltaChi2 = pre_chi2 - this_chi2;
547 if (deltaChi2 < -m_dchi2_min) break;
548 if (fabs(deltaChi2) < m_deltamin) num_small_dev++;
549 if (num_small_dev > m_conv_iter) break;
550 pre_chi2 = this_chi2;
551 }
552
553 // calculate chi^2
554
555 double chi2 = Chisq(ntuple, xval);
556 Hchi2.SetBinContent(chan + 1, chi2);
557 Hndf.SetBinContent(chan + 1, m_good);
558
559 // constrain sum of x to 2*syncTimeBase and calculate sample times, not necessary here
560
561 double sum = 0;
562 for (auto xi : xval) sum += xi;
563 if (sum == 0) {
564 B2ERROR("sum == 0");
565 return false;
566 }
567
568 double DeltaT = meanTimeDifference * (2 * m_syncTimeBase / static_cast<double>(c_TimeAxisSize));
569 HDeltaT.SetBinContent(chan + 1, DeltaT);
570
571 std::vector<double> timeInterval;
572 for (int i = 0; i < c_TimeAxisSize; i++)timeInterval.push_back(xval[i + 1] - xval[i]);
573
574
575 std::vector<double> sampleTimes;
576 for (auto xi : xval) sampleTimes.push_back(xi);
577
578 // save results as histograms
579 std::string forWhat = "scrod " + to_string(scrodID) + " channel " + to_string(chan);
580 saveAsHistogram(timeInterval, "dt_ch" + to_string(chan), "Sample time bins for " + forWhat,
581 "sample number", "#Delta t [ns]");
582 saveAsHistogram(sampleTimes, "sampleTimes_ch" + to_string(chan),
583 "Time base corrections for " + forWhat, "sample number", "t [ns]");
584
585 // calibrated cal pulse time difference
586 std::string name = "timeDiffcal_ch" + to_string(chan);
587 std::string title = "Calibrated cal pulse time difference vs. sample for " + forWhat;
588 TH2F Hcor(name.c_str(), title.c_str(), c_TimeAxisSize, 0, c_TimeAxisSize,
589 100, DeltaT - 0.5, DeltaT + 0.5);
590 Hcor.SetXTitle("sample number");
591 Hcor.SetYTitle("time difference [ns]");
592 Hcor.SetStats(kTRUE);
593
594 TOPSampleTimes timeBase;
595 timeBase.setTimeAxis(sampleTimes, sampleTimes.back() / 2);
596
597 for (const auto& twoTimes : ntuple) {
598 if (!twoTimes.good) continue;
599 double dt = timeBase.getDeltaTime(0, twoTimes.t2, twoTimes.t1);
600 int sample = int(twoTimes.t1) % c_TimeAxisSize;
601 Hcor.Fill(sample, dt);
602 }
603 Hcor.Write();
604
605 B2INFO("... channel " << chan << " OK (chi^2/ndf = " << chi2
606 << ", ndf = " << m_good << ")");
607
608 return true;
609 }
610
611 void TOPTimeBaseCalibratorModule::Iteration(const std::vector<TwoTimes>& ntuple, std::vector<double>& xval)
612 {
613 for (int i = 0; i < c_TimeAxisSize; i++) {
614 double wdth = xval[i + 1] - xval[i];
615 if (wdth < m_min_binwidth) {
616 xval[i] = xval[i] - 0.5 * fabs(wdth) - 0.5 * m_min_binwidth;
617 xval[i + 1] = xval[i + 1] + 0.5 * fabs(wdth) + 0.5 * m_min_binwidth;
618 }
619 if (wdth > m_max_binwidth) {
620 xval[i] = xval[i] - 0.5 * fabs(wdth) - 0.5 * m_max_binwidth;
621 xval[i + 1] = xval[i + 1] + 0.5 * fabs(wdth) + 0.5 * m_max_binwidth;
622 }
623 }
624
625 if (xval[0] != 0)
626 for (int i = 0; i < c_TimeAxisSize; i++) xval[i] = xval[i] - xval[0];
627
628 std::vector<double> xxval(c_TimeAxisSize + 1, 0.0);
629 for (int i = 0; i < c_TimeAxisSize + 1; i++) xxval[i] = xval[i];
630
631 double chi2_0 = Chisq(ntuple, xxval);
632 if (chi2_0 < 0) B2ERROR("iTBC chisq_0<0! xval has problem.");
633
634 std::vector<double> dr_chi2(c_TimeAxisSize + 1, 0.0);
635 TH1D hdrsamp_try("hdrsamp_try", "dchi2/dx distribution", 100, -0.01, 0.01);
636
637 for (int smp = 1; smp < c_TimeAxisSize; smp++) {
638 xxval[smp] = xval[smp] + m_dev_step;
639 double chi2_ch = Chisq(ntuple, xxval);
640 if (chi2_ch < 0)continue;
641 dr_chi2[smp] = (chi2_ch - chi2_0) / m_dev_step;
642 hdrsamp_try.Fill(dr_chi2[smp]);
643 xxval[smp] = xval[smp];
644 }
645
646 for (int smp = 1; smp < c_TimeAxisSize; smp++) {
647 double vx_it_step = dr_chi2[smp] * m_xstep;
648 xval[smp] = xval[smp] - vx_it_step;
649 }
650
651 //save rms of dchi2/dxval.
652 m_dchi2dxv = hdrsamp_try.GetRMS();
653 //change m_xstep
655 }
656
657 double TOPTimeBaseCalibratorModule::Chisq(const std::vector<TwoTimes>& ntuple, const std::vector<double>& xxval)
658 {
659 double sum1 = 0.0;
660 double sum2 = 0.0; //sum od dt and dt**2
661
662 m_good = 0;
663
664 for (const auto& twoTimes : ntuple) {
665 if (!twoTimes.good) continue;
666
667 int i1 = int(twoTimes.t1);
668 double fr = twoTimes.t1 - i1;
669 double samp0 = i1 % 256;
670 double ctdc1 = xxval[samp0] + fr * (xxval[samp0 + 1] - xxval[samp0]);
671 int i2 = int(twoTimes.t2);
672 fr = twoTimes.t2 - i2;
673 double samp1 = i2 % 256;
674 double ctdc2 = xxval[samp1] + fr * (xxval[samp1 + 1] - xxval[samp1]);
675 double cdt = 0.0;
676 if (samp1 > samp0) cdt = ctdc2 - ctdc1;
677 else cdt = ctdc2 - ctdc1 + m_syncTimeBase * 2;
678
679 if (cdt < m_dt_max && cdt > m_dt_min) {
680 sum1 += cdt;
681 sum2 += cdt * cdt;
682 m_good++;
683 }
684 }
685
686 double chi2 = -1.0;
687
688 if (m_good > 10) {
689 double mean = sum1 / m_good;
690 chi2 = (sum2 - m_good * mean * mean) / m_good / m_sigm2_exp;
691 }
692 return chi2;
693 }
694
695
696
697 void TOPTimeBaseCalibratorModule::saveAsHistogram(const std::vector<double>& vec,
698 const std::string& name,
699 const std::string& title,
700 const std::string& xTitle,
701 const std::string& yTitle)
702 {
703 if (vec.empty()) return;
704
705 TH1F h(name.c_str(), title.c_str(), vec.size(), 0, vec.size());
706 h.SetXTitle(xTitle.c_str());
707 h.SetYTitle(yTitle.c_str());
708 if (name.find("Fit") != string::npos) h.SetLineColor(2);
709
710 for (unsigned i = 0; i < vec.size(); i++) h.SetBinContent(i + 1, vec[i]);
711
712 h.Write();
713 }
714
715
716 void TOPTimeBaseCalibratorModule::saveAsHistogram(const vector<double>& vec,
717 const vector<double>& err,
718 const string& name,
719 const string& title,
720 const string& xTitle,
721 const string& yTitle)
722 {
723 if (vec.empty()) return;
724
725 TH1F h(name.c_str(), title.c_str(), vec.size(), 0, vec.size());
726 h.SetXTitle(xTitle.c_str());
727 h.SetYTitle(yTitle.c_str());
728
729 for (unsigned i = 0; i < vec.size(); i++) h.SetBinContent(i + 1, vec[i]);
730 for (unsigned i = 0; i < err.size(); i++) h.SetBinError(i + 1, err[i]);
731
732 h.Write();
733 }
734
735
737 const string& name,
738 const string& title)
739 {
740 int n = M.GetNrows();
741 TH2F h(name.c_str(), title.c_str(), n, 0, n, n, 0, n);
742 h.SetXTitle("columns");
743 h.SetYTitle("rows");
744
745 for (int j = 0; j < n; j++) {
746 for (int k = 0; k < n; k++) {
747 h.SetBinContent(j + 1, n - k, M(j, k));
748 }
749 }
750
751 h.Write();
752 }
753
754
755
757} // end Belle2 namespace
758
void setDescription(const std::string &description)
Sets the description of the module.
Definition Module.cc:214
Module()
Constructor.
Definition Module.cc:30
const TOPNominalTDC & getNominalTDC() const
Returns nominal time-to-digit conversion parameters.
double getSampleWidth() const
Returns time difference between two samples.
Class to store unpacked raw data (hits in feature-extraction format) It provides also calculation of ...
Definition TOPRawDigit.h:24
double getFWHM() const
Returns signal full width half maximum.
Calibration constants of a single ASIC channel: time axis (sample times)
double getDeltaTime(int window, double sample2, double sample1) const
Returns time difference between sample2 and sample1.
double m_deltamin
minimum chisq change in an iteration.
double m_xstep
unit for an interaction of delta(X_s)
TH2F m_calPulseSecond
pulse height versus width of the second calibration pulse
double m_sigm2_exp
(sigma_0(dT))**2 for nomarlization of chisq = sum{dT^2/sigma^2}
std::vector< TwoTimes > m_ntuples[c_NumChannels]
channel wise data
double m_min_binwidth
minimum time interval of one sample
TH2F m_goodHits
pulse height versus width of all good hits
double m_syncTimeBase
synchronization time (two ASIC windows)
double m_new_xstep
m_xstep = m_new_xstep if m_dchi2dxv < m_change_step
double m_dt_min
minimum Delta T of raw calpulse
bool m_saveMatrix
if true, save also matrix and its inverse in a root file
double m_dt_max
maximum Delta T of raw calpulse
double m_dev_step
a step size to calculate the value of d(chisq)/dxval
unsigned m_conv_iter
Number of iteration with chisq changes less than deltamin.
double m_dchi2dxv
rms of 255 dchi2/dxval values
double m_max_binwidth
maximum time interval of one sample
unsigned m_minHits
minimal required hits per channel
double m_minTimeDiff
lower bound for time difference [samples]
double m_change_xstep
update m_xstep if m_dchi2dxv < m_change_step
StoreArray< TOPDigit > m_digits
collection of digits
unsigned m_numIterations
Number of Iterations of iTBC.
std::string m_directoryName
directory name for the output root files
int m_good
good events used for chisq cal.
TH2F m_calPulseFirst
pulse height versus width of the first calibration pulse
double m_dchi2_min
quit if chisq increase in iterations is larger than this value.
bool m_useFallingEdge
if true, use falling edge instead of rising
double m_maxTimeDiff
upper bound for time difference [samples]
@ c_NumBoardstacks
number of boardstacks per module
@ c_TimeAxisSize
number of samples to calibrate
const TOPGeometry * getGeometry() const
Returns pointer to geometry object using basf2 units.
static TOPGeometryPar * Instance()
Static method to obtain the pointer to its instance.
const FrontEndMapper & getFrontEndMapper() const
Returns front-end mapper (mapping of SCROD's to positions within TOP modules)
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition Module.h:559
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition Module.h:649
double sqrt(double a)
sqrt for double
Definition beamHelpers.h:28
bool matrixInversion(const std::vector< TwoTimes > &ntuple, unsigned scrodID, unsigned scrodChannel, double meanTimeDifference, TH1F &Hchi2, TH1F &Hndf, TH1F &HDeltaT)
Method by matrix inversion.
bool calibrateChannel(std::vector< TwoTimes > &ntuple, unsigned scrodID, unsigned scrodChannel, TH1F &Hchi2, TH1F &Hndf, TH1F &HDeltaT)
calibrate single channel
virtual void initialize() override
Initialize the Module.
virtual void event() override
Event processor.
virtual void terminate() override
Termination action.
void Iteration(const std::vector< TwoTimes > &ntuple, std::vector< double > &xval)
Iteration function called by iterativeTBC()
double Chisq(const std::vector< TwoTimes > &ntuple, const std::vector< double > &xxval)
Return the chisqure of one set of TBC constants (xval) in iTBC calculation.
bool iterativeTBC(const std::vector< TwoTimes > &ntuple, unsigned scrodID, unsigned scrodChannel, double meanTimeDifference, TH1F &Hchi2, TH1F &Hndf, TH1F &HDeltaT)
Method by iteration of chi2 minimization.
static void saveAsHistogram(const std::vector< double > &vec, const std::string &name, const std::string &title, const std::string &xTitle="", const std::string &yTitle="")
Save vector to histogram and write it out.
void setTimeAxis(double syncTimeBase)
Sets equidistant time axis (uncalibrated).
Abstract base class for different kinds of events.
STL namespace.
Structure to hold some of the calpulse data.
Structure to hold calpulse raw times expressed in samples since sample 0 of window 0.