Belle II Software  release-08-01-10
TOPModuleT0DeltaTCollectorModule.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 #include <top/modules/collectors/TOPModuleT0DeltaTCollectorModule.h>
10 #include <string>
11 #include <TH1F.h>
12 #include <TH2F.h>
13 
14 using namespace std;
15 
16 namespace Belle2 {
22  //-----------------------------------------------------------------
24  //-----------------------------------------------------------------
25 
26  REG_MODULE(TOPModuleT0DeltaTCollector);
27 
28  //-----------------------------------------------------------------
29  // Implementation
30  //-----------------------------------------------------------------
31 
32  TOPModuleT0DeltaTCollectorModule::TOPModuleT0DeltaTCollectorModule()
33  {
34  // set module description and processing properties
35  setDescription("A collector for module T0 calibration with chi2 minimization of "
36  "time differences between slots (method DeltaT). Useful for initial "
37  "(rough) calibration, since the results are found slightly biased. "
38  "For the final (precise) calibration one has to use LL method.");
39  setPropertyFlags(c_ParallelProcessingCertified);
40 
41  // module parameters
42  addParam("numBins", m_numBins,
43  "number of bins of histograms of time difference", 400);
44  addParam("timeRange", m_timeRange,
45  "time range [ns] of histograms of time difference", 20.0);
46 
47  }
48 
49 
50  void TOPModuleT0DeltaTCollectorModule::prepare()
51  {
52 
53  m_timeZeros.isRequired();
54 
55  auto slotPairs = new TH2F("slots", "slot pairs: number of events",
56  16, 0.5, 16.5, 16, 0.5, 16.5);
57  slotPairs->SetXTitle("first slot number");
58  slotPairs->SetYTitle("second slot number");
59  registerObject<TH2F>("slots", slotPairs);
60 
61  double xmin = -m_timeRange / 2;
62  double xmax = m_timeRange / 2;
63  for (int slot1 = 1; slot1 <= 9; slot1++) {
64  for (int slot2 = slot1 + 7; slot2 <= slot1 + 9; slot2++) {
65  if (slot2 > 16) continue;
66  string name = "deltaT0_" + to_string(slot1) + "-" + to_string(slot2);
67  string title = "time difference: slot " + to_string(slot1) + " - slot "
68  + to_string(slot2);
69  auto h = new TH1F(name.c_str(), title.c_str(), m_numBins, xmin, xmax);
70  h->SetXTitle("time difference [ns]");
71  registerObject<TH1F>(name, h);
72  }
73  }
74 
75  }
76 
77 
78  void TOPModuleT0DeltaTCollectorModule::collect()
79  {
80 
81  if (m_timeZeros.getEntries() != 2) return;
82 
83  const auto* timeZero1 = m_timeZeros[0];
84  const auto* timeZero2 = m_timeZeros[1];
85  if (timeZero1->getModuleID() > timeZero2->getModuleID()) {
86  timeZero1 = m_timeZeros[1];
87  timeZero2 = m_timeZeros[0];
88  }
89 
90  int slot1 = timeZero1->getModuleID();
91  int slot2 = timeZero2->getModuleID();
92  if (slot1 > 9 or slot2 - slot1 < 7 or slot2 - slot1 > 9) return;
93  string name = "deltaT0_" + to_string(slot1) + "-" + to_string(slot2);
94  auto h = getObjectPtr<TH1F>(name);
95  if (not h) return;
96 
97  double t1 = timeZero1->getTime();
98  if (m_moduleT0->isCalibrated(slot1)) t1 += m_moduleT0->getT0(slot1);
99  double t2 = timeZero2->getTime();
100  if (m_moduleT0->isCalibrated(slot2)) t2 += m_moduleT0->getT0(slot2);
101  h->Fill(t1 - t2);
102 
103  auto slotPairs = getObjectPtr<TH2F>("slots");
104  slotPairs->Fill(slot1, slot2);
105  }
106 
108 }
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.