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