Belle II Software  release-05-02-19
TOPTimeBaseCalibratorModule.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2016 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Marko Staric, Hiromichi Kichimi and Xiaolong Wang *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 // Own include
12 #include <top/modules/TOPTimeBaseCalibrator/TOPTimeBaseCalibratorModule.h>
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 
34 using namespace std;
35 
36 namespace Belle2 {
42  using namespace TOP;
43 
44  //-----------------------------------------------------------------
45  // Register module
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 interation 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 
99  }
100 
101  TOPTimeBaseCalibratorModule::~TOPTimeBaseCalibratorModule()
102  {
103  }
104 
105  void TOPTimeBaseCalibratorModule::initialize()
106  {
107 
108  // input
109  m_digits.isRequired();
110 
111  // checks
112  const auto* geo = TOPGeometryPar::Instance()->getGeometry();
113  if (!geo->isModuleIDValid(m_moduleID))
114  B2ERROR("Invalid module ID: " << m_moduleID);
115 
116  // check for existance and mkdir if not
117  if (m_directoryName.empty()) m_directoryName = "./";
118  if (m_directoryName != "./") gSystem->mkdir(m_directoryName.c_str(), kTRUE);
119 
120  // synchronization time corresponding to two ASIC windows
121  m_syncTimeBase = geo->getNominalTDC().getSyncTimeBase();
122 
123  // control histograms
124  m_goodHits = TH2F("goodHits", "pulse height vs. pulse width of all good hits",
125  100, 0, 10, 100, 0, 2000);
126  m_goodHits.SetXTitle("pulse width (FWHM) [ns]");
127  m_goodHits.SetYTitle("pulse height [ADC counts]");
128  m_calPulseFirst = TH2F("calPulseFirst",
129  "pulse height vs. pulse width of the first calibration pulse",
130  100, 0, 10, 100, 0, 2000);
131  m_calPulseFirst.SetXTitle("pulse width (FWHM) [ns]");
132  m_calPulseFirst.SetYTitle("pulse height [ADC counts]");
133  m_calPulseSecond = TH2F("calPulseSecond",
134  "pulse height vs. pulse width of the second calibration pulse",
135  100, 0, 10, 100, 0, 2000);
136  m_calPulseSecond.SetXTitle("pulse width (FWHM) [ns]");
137  m_calPulseSecond.SetYTitle("pulse height [ADC counts]");
138 
139  }
140 
141 
142  void TOPTimeBaseCalibratorModule::beginRun()
143  {
144  }
145 
146  void TOPTimeBaseCalibratorModule::event()
147  {
148  const auto* geo = TOPGeometryPar::Instance()->getGeometry();
149  double sampleWidth = geo->getNominalTDC().getSampleWidth();
150 
151  vector<Hit> hits[c_NumChannels];
152 
153  for (const auto& digit : m_digits) {
154  if (digit.getModuleID() != m_moduleID) continue;
155  if (digit.getHitQuality() != TOPDigit::c_Junk) {
156  m_goodHits.Fill(digit.getPulseWidth(), digit.getPulseHeight());
157  }
158  if (digit.getHitQuality() != TOPDigit::c_CalPulse) continue;
159  double rawTime = digit.getRawTime();
160  double errScaleFactor = 1;
161  if (m_useFallingEdge) {
162  const auto* rawDigit = digit.getRelated<TOPRawDigit>();
163  if (!rawDigit) {
164  B2ERROR("No relation to TOPRawDigit - can't determine falling edge time error");
165  continue;
166  }
167  // rawTime may include corrections due to window number discontinuity,
168  // therefore one must add the width and not just use getCFDFallingTime()
169  rawTime += rawDigit->getFWHM();
170  errScaleFactor = rawDigit->getCFDFallingTimeError(1.0) / rawDigit->getCFDLeadingTimeError(1.0);
171  }
172  double t = rawTime + digit.getFirstWindow() * c_WindowSize;
173  if (t < 0) {
174  B2ERROR("Got negative sample number - digit ignored");
175  continue;
176  }
177  double et = digit.getTimeError() / sampleWidth * errScaleFactor;
178  if (et <= 0) {
179  B2ERROR("Time error is not given - digit ignored");
180  continue;
181  }
182  unsigned channel = digit.getChannel();
183  if (channel < c_NumChannels) {
184  hits[channel].push_back(Hit(t, et, digit.getPulseHeight(), digit.getPulseWidth()));
185  }
186  }
187 
188  for (unsigned channel = 0; channel < c_NumChannels; channel++) {
189  const auto& channelHits = hits[channel];
190  if (channelHits.size() == 2) {
191  double t0 = channelHits[0].time;
192  double t1 = channelHits[1].time;
193  auto diff = fabs(t0 - t1); // since not sorted yet
194  if (diff < m_minTimeDiff) continue;
195  if (diff > m_maxTimeDiff) continue;
196  double sig0 = channelHits[0].timeErr;
197  double sig1 = channelHits[1].timeErr;
198  double sigma = sqrt(sig0 * sig0 + sig1 * sig1);
199  m_ntuples[channel].push_back(TwoTimes(t0, t1, sigma));
200  if (t0 < t1) { // check, since not sorted yet
201  m_calPulseFirst.Fill(channelHits[0].pulseWidth, channelHits[0].pulseHeight);
202  m_calPulseSecond.Fill(channelHits[1].pulseWidth, channelHits[1].pulseHeight);
203  } else {
204  m_calPulseSecond.Fill(channelHits[0].pulseWidth, channelHits[0].pulseHeight);
205  m_calPulseFirst.Fill(channelHits[1].pulseWidth, channelHits[1].pulseHeight);
206  }
207  } else if (channelHits.size() > 2) {
208  B2WARNING("More than two cal pulses per channel found - ignored");
209  }
210  }
211 
212  }
213 
214 
215  void TOPTimeBaseCalibratorModule::endRun()
216  {
217  }
218 
219  void TOPTimeBaseCalibratorModule::terminate()
220  {
221 
222  const auto& feMapper = TOPGeometryPar::Instance()->getFrontEndMapper();
223 
224  for (unsigned bs = 0; bs < c_NumBoardstacks; bs ++) {
225 
226  // determine scrod ID
227 
228  const auto* feMap = feMapper.getMap(m_moduleID, bs);
229  if (!feMap) {
230  B2ERROR("No front-end map available for slot " << m_moduleID << ", BS" << bs);
231  continue;
232  }
233  unsigned scrodID = feMap->getScrodID();
234 
235  // open output root file
236 
237  string fileName = m_directoryName + "/tbcSlot";
238  if (m_moduleID < 10) {
239  fileName += "0" + to_string(m_moduleID);
240  } else {
241  fileName += to_string(m_moduleID);
242  }
243  fileName += "_" + to_string(bs) + "-scrod" +
244  to_string(scrodID) + ".root";
245  TFile* fout = TFile::Open(fileName.c_str(), "recreate");
246  if (!fout) {
247  B2ERROR("Can't open the output file " << fileName);
248  continue;
249  }
250 
251  // write control histograms
252 
253  m_goodHits.Write();
254  m_calPulseFirst.Write();
255  m_calPulseSecond.Write();
256 
257  B2INFO("Fitting time base corrections for SCROD " << scrodID
258  << ", output file: " << fileName);
259 
260  // book histograms
261 
262  TH1F Hchan("channels", "channels with cal pulse",
263  c_NumScrodChannels, 0, c_NumScrodChannels);
264  Hchan.SetXTitle("channel number");
265  Hchan.SetYTitle("double cal pulse counts");
266  TH1F Hsuccess("success", "successfuly fitted channels",
267  c_NumScrodChannels, 0, c_NumScrodChannels);
268  Hsuccess.SetXTitle("channel number");
269  TH1F Hchi2("chi2", "normalized chi2",
270  c_NumScrodChannels, 0, c_NumScrodChannels);
271  Hchi2.SetXTitle("channel number");
272  Hchi2.SetYTitle("chi2/ndf");
273  TH1F Hndf("ndf", "degrees of freedom",
274  c_NumScrodChannels, 0, c_NumScrodChannels);
275  Hndf.SetXTitle("channel number");
276  Hndf.SetYTitle("ndf");
277  TH1F HDeltaT("DeltaT", "Fitted double pulse delay time",
278  c_NumScrodChannels, 0, c_NumScrodChannels);
279  HDeltaT.SetXTitle("channel number");
280  HDeltaT.SetYTitle("#Delta T [ns]");
281 
282  // loop over channels of a single SCROD
283 
284  for (unsigned chan = 0; chan < c_NumScrodChannels; chan++) {
285  unsigned channel = chan + bs * c_NumScrodChannels;
286  auto& ntuple = m_ntuples[channel];
287  Hchan.SetBinContent(chan + 1, ntuple.size());
288  if (ntuple.size() > m_minHits) {
289  bool ok = calibrateChannel(ntuple, scrodID, chan, Hchi2, Hndf, HDeltaT);
290  if (ok) Hsuccess.Fill(chan);
291  } else {
292  B2INFO("... channel " << chan << " statistics too low ("
293  << ntuple.size() << " double cal pulses)");
294  }
295  }
296 
297  // write histograms and close the file
298 
299  Hchan.Write();
300  Hsuccess.Write();
301  Hchi2.Write();
302  Hndf.Write();
303  HDeltaT.Write();
304  fout->Close();
305 
306  }
307  }
308 
309 
310  bool TOPTimeBaseCalibratorModule::calibrateChannel(vector<TwoTimes>& ntuple,
311  unsigned scrodID, unsigned chan,
312  TH1F& Hchi2, TH1F& Hndf,
313  TH1F& HDeltaT)
314  {
315 
316  // determine outlayer removal cuts
317 
318  string name = "timeDiff_ch" + to_string(chan);
319  string forWhat = "scrod " + to_string(scrodID) + " channel " + to_string(chan);
320  string title = "Cal pulse time difference vs. sample for " + forWhat;
321  TProfile prof(name.c_str(), title.c_str(), c_TimeAxisSize, 0, c_TimeAxisSize,
322  m_minTimeDiff, m_maxTimeDiff, "S");
323  prof.SetXTitle("sample number");
324  prof.SetYTitle("time difference [samples]");
325 
326  name = "sampleOccup_ch" + to_string(chan);
327  title = "Occupancy for " + forWhat;
328  TH1F hist(name.c_str(), title.c_str(), c_TimeAxisSize, 0, c_TimeAxisSize);
329  hist.SetXTitle("sample number");
330  hist.SetYTitle("entries per sample");
331 
332  vector<double> profMeans(c_TimeAxisSize, (m_maxTimeDiff + m_minTimeDiff) / 2);
333  double sigma = (m_maxTimeDiff - m_minTimeDiff) / 6;
334 
335  for (int iter = 0; iter < 5; iter++) {
336  prof.Reset();
337  double x = 0;
338  double x2 = 0;
339  int n = 0;
340  for (auto& twoTimes : ntuple) {
341  int sample = int(twoTimes.t1) % c_TimeAxisSize;
342  double diff = twoTimes.t2 - twoTimes.t1;
343  double Ddiff = diff - profMeans[sample];
344  if (fabs(Ddiff) < 3 * sigma) {
345  prof.Fill(sample, diff);
346  hist.Fill(sample);
347  x += Ddiff;
348  x2 += Ddiff * Ddiff;
349  n++;
350  twoTimes.good = true;
351  } else {
352  twoTimes.good = false;
353  }
354  }
355  for (int i = 0; i < c_TimeAxisSize; i++) {
356  profMeans[i] = prof.GetBinContent(i + 1);
357  }
358  if (n == 0) return false;
359  x2 /= n;
360  x /= n;
361  sigma = sqrt(x2 - x * x);
362  }
363  prof.Write();
364  hist.Write();
365 
366  // calculate average time difference
367 
368  double meanTimeDifference = 0;
369  for (auto& x : profMeans) meanTimeDifference += x;
370  meanTimeDifference /= profMeans.size();
371  meanTimeDifference -= int(meanTimeDifference / c_TimeAxisSize) * c_TimeAxisSize;
372 
373  // perform the fit
374 
375  switch (m_method) {
376  case 0:
377  B2INFO("... channel " << chan << ": "
378  << ntuple.size() << " double cal pulses)");
379  return false;
380  case 1:
381  return matrixInversion(ntuple, scrodID, chan, meanTimeDifference,
382  Hchi2, Hndf, HDeltaT);
383  case 2:
384  return iterativeTBC(ntuple, scrodID, chan, meanTimeDifference,
385  Hchi2, Hndf, HDeltaT);
386  case 3:
387  B2ERROR("Singuler value decomposition not implemented yet");
388  return false;
389  default:
390  B2ERROR("Unknown method " << m_method);
391  return false;
392  }
393 
394  }
395 
396 
397  bool TOPTimeBaseCalibratorModule::matrixInversion(const vector<TwoTimes>& ntuple,
398  unsigned scrodID, unsigned chan,
399  double meanTimeDifference,
400  TH1F& Hchi2, TH1F& Hndf,
401  TH1F& HDeltaT)
402  {
403 
404  // Ax = b: construct matrix A and right side vector b
405 
406  TMatrixDSym A(c_TimeAxisSize);
407  vector<double> b(c_TimeAxisSize, 0.0);
408 
409  for (const auto& twoTimes : ntuple) {
410  if (!twoTimes.good) continue;
411 
412  vector<double> m(c_TimeAxisSize, 0.0);
413  int i1 = int(twoTimes.t1);
414  m[i1 % c_TimeAxisSize] = 1.0 - (twoTimes.t1 - i1);
415  int i2 = int(twoTimes.t2);
416  m[i2 % c_TimeAxisSize] = twoTimes.t2 - i2;
417  i2 = i1 + (i2 - i1) % c_TimeAxisSize;
418  for (int k = i1 + 1; k < i2; k++) m[k % c_TimeAxisSize] = 1;
419 
420  double relSigma = twoTimes.sigma / meanTimeDifference;
421  double sig2 = relSigma * relSigma;
422 
423  for (int jj = i1; jj < i2 + 1; jj++) {
424  int j = jj % c_TimeAxisSize;
425  for (int kk = i1; kk < i2 + 1; kk++) {
426  int k = kk % c_TimeAxisSize;
427  A(j, k) += m[j] * m[k] / sig2;
428  }
429  }
430  for (int k = 0; k < c_TimeAxisSize; k++) b[k] += m[k] / sig2;
431  }
432 
433  // save as histograms
434 
435  string forWhat = "scrod " + to_string(scrodID) + " channel " + to_string(chan);
436  saveAsHistogram(A, "matA_ch" + to_string(chan), "Matrix for " + forWhat);
437  saveAsHistogram(b, "vecB_ch" + to_string(chan), "Right side for " + forWhat);
438 
439  // invert matrix A and solve the equation: x = A^{-1}b
440 
441  double det = 0;
442  A.Invert(&det);
443  if (det == 0) {
444  B2INFO("... channel " << chan << " failed");
445  return false;
446  }
447 
448  vector<double> x(c_TimeAxisSize, 0.0);
449  for (int k = 0; k < c_TimeAxisSize; k++) {
450  for (int j = 0; j < c_TimeAxisSize; j++) {
451  x[k] += A(k, j) * b[j];
452  }
453  }
454 
455  // calculate chi^2
456 
457  double chi2 = 0;
458  int ndf = -c_TimeAxisSize;
459 
460  for (const auto& twoTimes : ntuple) {
461  if (!twoTimes.good) continue;
462 
463  vector<double> m(c_TimeAxisSize, 0.0);
464  int i1 = int(twoTimes.t1);
465  m[i1 % c_TimeAxisSize] = 1.0 - (twoTimes.t1 - i1);
466  int i2 = int(twoTimes.t2);
467  m[i2 % c_TimeAxisSize] = twoTimes.t2 - i2;
468  i2 = i1 + (i2 - i1) % c_TimeAxisSize;
469  for (int k = i1 + 1; k < i2; k++) m[k % c_TimeAxisSize] = 1;
470  double s = -1.0;
471  for (int k = 0; k < c_TimeAxisSize; k++) s += m[k] * x[k];
472  double relSigma = twoTimes.sigma / meanTimeDifference;
473  double sig2 = relSigma * relSigma;
474  chi2 += s * s / sig2;
475  ndf++;
476  }
477  Hchi2.SetBinContent(chan + 1, chi2 / ndf);
478  Hndf.SetBinContent(chan + 1, ndf);
479 
480  // constrain sum of x to 2*syncTimeBase and calculate sample times
481 
482  double sum = 0;
483  for (auto xi : x) sum += xi;
484  if (sum == 0) {
485  B2ERROR("sum == 0");
486  return false;
487  }
488  double DeltaT = 2 * m_syncTimeBase / sum;
489  for (auto& xi : x) xi *= DeltaT;
490  HDeltaT.SetBinContent(chan + 1, DeltaT);
491 
492  vector<double> err;
493  for (int k = 0; k < c_TimeAxisSize; k++) err.push_back(sqrt(A(k, k)) * DeltaT);
494 
495  vector<double> sampleTimes;
496  sampleTimes.push_back(0);
497  for (auto xi : x) sampleTimes.push_back(xi + sampleTimes.back());
498 
499  // save results as histograms
500 
501  saveAsHistogram(A, "invA_ch" + to_string(chan), "Inverted matrix for " + forWhat);
502  saveAsHistogram(x, err, "dt_ch" + to_string(chan), "Sample time bins for " + forWhat,
503  "sample number", "#Delta t [ns]");
504  saveAsHistogram(sampleTimes, "sampleTimes_ch" + to_string(chan),
505  "Time base corrections for " + forWhat, "sample number", "t [ns]");
506 
507  // calibrated cal pulse time difference
508 
509  string name = "timeDiffcal_ch" + to_string(chan);
510  string title = "Calibrated cal pulse time difference vs. sample for " + forWhat;
511  TH2F Hcor(name.c_str(), title.c_str(), c_TimeAxisSize, 0, c_TimeAxisSize,
512  100, DeltaT - 0.5, DeltaT + 0.5);
513  Hcor.SetXTitle("sample number");
514  Hcor.SetYTitle("time difference [ns]");
515  Hcor.SetStats(kTRUE);
516 
517  TOPSampleTimes timeBase;
518  timeBase.setTimeAxis(sampleTimes, sampleTimes.back() / 2);
519 
520  for (const auto& twoTimes : ntuple) {
521  if (!twoTimes.good) continue;
522  double dt = timeBase.getDeltaTime(0, twoTimes.t2, twoTimes.t1);
523  int sample = int(twoTimes.t1) % c_TimeAxisSize;
524  Hcor.Fill(sample, dt);
525  }
526  Hcor.Write();
527 
528  B2INFO("... channel " << chan << " OK (chi^2/ndf = " << chi2 / ndf
529  << ", ndf = " << ndf << ")");
530 
531  return true;
532  }
533 
534 
535  bool TOPTimeBaseCalibratorModule::iterativeTBC(const std::vector<TwoTimes>& ntuple,
536  unsigned scrodID, unsigned chan,
537  double meanTimeDifference,
538  TH1F& Hchi2, TH1F& Hndf,
539  TH1F& HDeltaT)
540  {
541  std::vector<double> xval(c_TimeAxisSize + 1, 0.0);
542  double wx = 2 * m_syncTimeBase / c_TimeAxisSize;
543  for (int i = 0; i < c_TimeAxisSize + 1; i++) xval[i] = i * wx;
544 
545  B2INFO("TimeBaseCalibration starts for channel#" << chan);
546 
547  double pre_chi2 = 10000000.0;
548  unsigned num_small_dev = 0;
549 
550  for (unsigned j = 0; j < m_numIterations; j++) {
551 
552  Iteration(ntuple, xval);
553  double this_chi2 = Chisq(ntuple, xval);
554  if (this_chi2 < 0)continue;
555  double deltaChi2 = pre_chi2 - this_chi2;
556  if (deltaChi2 < -m_dchi2_min) break;
557  if (fabs(deltaChi2) < m_deltamin) num_small_dev++;
558  if (num_small_dev > m_conv_iter) break;
559  pre_chi2 = this_chi2;
560  }
561 
562  // calculate chi^2
563 
564  double chi2 = Chisq(ntuple, xval);
565  Hchi2.SetBinContent(chan + 1, chi2);
566  Hndf.SetBinContent(chan + 1, m_good);
567 
568  // constrain sum of x to 2*syncTimeBase and calculate sample times, not necessary here
569 
570  double sum = 0;
571  for (auto xi : xval) sum += xi;
572  if (sum == 0) {
573  B2ERROR("sum == 0");
574  return false;
575  }
576 
577  double DeltaT = meanTimeDifference * (2 * m_syncTimeBase / c_TimeAxisSize);
578  HDeltaT.SetBinContent(chan + 1, DeltaT);
579 
580  std::vector<double> timeInterval;
581  for (int i = 0; i < c_TimeAxisSize; i++)timeInterval.push_back(xval[i + 1] - xval[i]);
582 
583 
584  std::vector<double> sampleTimes;
585  for (auto xi : xval) sampleTimes.push_back(xi);
586 
587  // save results as histograms
588  std::string forWhat = "scrod " + to_string(scrodID) + " channel " + to_string(chan);
589  saveAsHistogram(timeInterval, "dt_ch" + to_string(chan), "Sample time bins for " + forWhat,
590  "sample number", "#Delta t [ns]");
591  saveAsHistogram(sampleTimes, "sampleTimes_ch" + to_string(chan),
592  "Time base corrections for " + forWhat, "sample number", "t [ns]");
593 
594  // calibrated cal pulse time difference
595  std::string name = "timeDiffcal_ch" + to_string(chan);
596  std::string title = "Calibrated cal pulse time difference vs. sample for " + forWhat;
597  TH2F Hcor(name.c_str(), title.c_str(), c_TimeAxisSize, 0, c_TimeAxisSize,
598  100, DeltaT - 0.5, DeltaT + 0.5);
599  Hcor.SetXTitle("sample number");
600  Hcor.SetYTitle("time difference [ns]");
601  Hcor.SetStats(kTRUE);
602 
603  TOPSampleTimes timeBase;
604  timeBase.setTimeAxis(sampleTimes, sampleTimes.back() / 2);
605 
606  for (const auto& twoTimes : ntuple) {
607  if (!twoTimes.good) continue;
608  double dt = timeBase.getDeltaTime(0, twoTimes.t2, twoTimes.t1);
609  int sample = int(twoTimes.t1) % c_TimeAxisSize;
610  Hcor.Fill(sample, dt);
611  }
612  Hcor.Write();
613 
614  B2INFO("... channel " << chan << " OK (chi^2/ndf = " << chi2
615  << ", ndf = " << m_good << ")");
616 
617  return true;
618  }
619 
620  void TOPTimeBaseCalibratorModule::Iteration(const std::vector<TwoTimes>& ntuple, std::vector<double>& xval)
621  {
622  for (int i = 0; i < c_TimeAxisSize; i++) {
623  double wdth = xval[i + 1] - xval[i];
624  if (wdth < m_min_binwidth) {
625  xval[i] = xval[i] - 0.5 * fabs(wdth) - 0.5 * m_min_binwidth;
626  xval[i + 1] = xval[i + 1] + 0.5 * fabs(wdth) + 0.5 * m_min_binwidth;
627  }
628  if (wdth > m_max_binwidth) {
629  xval[i] = xval[i] - 0.5 * fabs(wdth) - 0.5 * m_max_binwidth;
630  xval[i + 1] = xval[i + 1] + 0.5 * fabs(wdth) + 0.5 * m_max_binwidth;
631  }
632  }
633 
634  if (xval[0] != 0)
635  for (int i = 0; i < c_TimeAxisSize; i++) xval[i] = xval[i] - xval[0];
636 
637  std::vector<double> xxval(c_TimeAxisSize + 1, 0.0);
638  for (int i = 0; i < c_TimeAxisSize + 1; i++) xxval[i] = xval[i];
639 
640  double chi2_0 = Chisq(ntuple, xxval);
641  if (chi2_0 < 0) B2ERROR("iTBC chisq_0<0! xval has problem.");
642 
643  std::vector<double> dr_chi2(c_TimeAxisSize + 1, 0.0);
644  TH1D hdrsamp_try("hdrsamp_try", "dchi2/dx distribution", 100, -0.01, 0.01);
645 
646  for (int smp = 1; smp < c_TimeAxisSize; smp++) {
647  xxval[smp] = xval[smp] + m_dev_step;
648  double chi2_ch = Chisq(ntuple, xxval);
649  if (chi2_ch < 0)continue;
650  dr_chi2[smp] = (chi2_ch - chi2_0) / m_dev_step;
651  hdrsamp_try.Fill(dr_chi2[smp]);
652  xxval[smp] = xval[smp];
653  }
654 
655  for (int smp = 1; smp < c_TimeAxisSize; smp++) {
656  double vx_it_step = dr_chi2[smp] * m_xstep;
657  xval[smp] = xval[smp] - vx_it_step;
658  }
659 
660  //save rms of dchi2/dxval.
661  m_dchi2dxv = hdrsamp_try.GetRMS();
662  //change m_xstep
663  if (fabs(m_dchi2dxv) < m_change_xstep) m_xstep = m_new_xstep;
664  }
665 
666  double TOPTimeBaseCalibratorModule::Chisq(const std::vector<TwoTimes>& ntuple, const std::vector<double>& xxval)
667  {
668  double sum1 = 0.0;
669  double sum2 = 0.0; //sum od dt and dt**2
670 
671  m_good = 0;
672 
673  for (const auto& twoTimes : ntuple) {
674  if (!twoTimes.good) continue;
675 
676  std::vector<double> m(c_TimeAxisSize, 0.0);
677 
678  int i1 = int(twoTimes.t1);
679  double fr = twoTimes.t1 - i1;
680  double samp0 = i1 % 256;
681  double ctdc1 = xxval[samp0] + fr * (xxval[samp0 + 1] - xxval[samp0]);
682  int i2 = int(twoTimes.t2);
683  fr = twoTimes.t2 - i2;
684  double samp1 = i2 % 256;
685  double ctdc2 = xxval[samp1] + fr * (xxval[samp1 + 1] - xxval[samp1]);
686  double cdt = 0.0;
687  if (samp1 > samp0) cdt = ctdc2 - ctdc1;
688  else cdt = ctdc2 - ctdc1 + m_syncTimeBase * 2;
689 
690  if (cdt < m_dt_max && cdt > m_dt_min) {
691  sum1 += cdt;
692  sum2 += cdt * cdt;
693  m_good++;
694  }
695  }
696 
697  double chi2 = -1.0;
698 
699  if (m_good > 10) {
700  double mean = sum1 / m_good;
701  chi2 = (sum2 - m_good * mean * mean) / m_good / m_sigm2_exp;
702  }
703  return chi2;
704  }
705 
706 
707 
708  void TOPTimeBaseCalibratorModule::saveAsHistogram(const std::vector<double>& vec,
709  const std::string& name,
710  const std::string& title,
711  const std::string& xTitle,
712  const std::string& yTitle) const
713  {
714  if (vec.empty()) return;
715 
716  TH1F h(name.c_str(), title.c_str(), vec.size(), 0, vec.size());
717  h.SetXTitle(xTitle.c_str());
718  h.SetYTitle(yTitle.c_str());
719  if (name.find("Fit") != string::npos) h.SetLineColor(2);
720 
721  for (unsigned i = 0; i < vec.size(); i++) h.SetBinContent(i + 1, vec[i]);
722 
723  h.Write();
724  }
725 
726 
727  void TOPTimeBaseCalibratorModule::saveAsHistogram(const vector<double>& vec,
728  const vector<double>& err,
729  const string& name,
730  const string& title,
731  const string& xTitle,
732  const string& yTitle) const
733  {
734  if (vec.empty()) return;
735 
736  TH1F h(name.c_str(), title.c_str(), vec.size(), 0, vec.size());
737  h.SetXTitle(xTitle.c_str());
738  h.SetYTitle(yTitle.c_str());
739 
740  for (unsigned i = 0; i < vec.size(); i++) h.SetBinContent(i + 1, vec[i]);
741  for (unsigned i = 0; i < err.size(); i++) h.SetBinError(i + 1, err[i]);
742 
743  h.Write();
744  }
745 
746 
747  void TOPTimeBaseCalibratorModule::saveAsHistogram(const TMatrixDSym& M,
748  const string& name,
749  const string& title) const
750  {
751  int n = M.GetNrows();
752  TH2F h(name.c_str(), title.c_str(), n, 0, n, n, 0, n);
753  h.SetXTitle("columns");
754  h.SetYTitle("rows");
755 
756  for (int j = 0; j < n; j++) {
757  for (int k = 0; k < n; k++) {
758  h.SetBinContent(j + 1, n - k, M(j, k));
759  }
760  }
761 
762  h.Write();
763  }
764 
765 
766 
768 } // end Belle2 namespace
769 
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::TOPTimeBaseCalibratorModule
Time base calibrator.
Definition: TOPTimeBaseCalibratorModule.h:87
Belle2::Module
Base class for Modules.
Definition: Module.h:74
Belle2::TOPRawDigit
Class to store unpacked raw data (hits in feature-extraction format) It provides also calculation of ...
Definition: TOPRawDigit.h:34
Belle2::TOPSampleTimes
Calibration constants of a singe ASIC channel: time axis (sample times)
Definition: TOPSampleTimes.h:34
Belle2::TwoTimes
Structure to hold calpulse raw times expressed in samples since sample 0 of window 0.
Definition: TOPTimeBaseCalibratorModule.h:61
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::TOPSampleTimes::setTimeAxis
void setTimeAxis(double syncTimeBase)
Sets equidistant time axis (uncalibrated).
Definition: TOPSampleTimes.cc:25
Belle2::TOPSampleTimes::getDeltaTime
double getDeltaTime(int window, double sample2, double sample1) const
Returns time difference between sample2 and sample1.
Definition: TOPSampleTimes.h:167
Belle2::Hit
Structure to hold some of the calpulse data.
Definition: TOPTimeBaseCalibratorModule.h:40