Belle II Software  release-08-01-10
TOPLaserCalibratorCollectorModule.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/TOPLaserCalibratorCollector/TOPLaserCalibratorCollectorModule.h>
10 
11 //ROOT
12 #include <TTree.h>
13 #include <TMath.h>
14 
15 //TOP
16 #include <top/dataobjects/TOPDigit.h>
17 #include <top/dataobjects/TOPRawDigit.h>
18 #include <top/dataobjects/TOPSimHit.h>
19 
20 #include <algorithm>
21 #include <vector>
22 
23 using namespace Belle2;
24 
25 REG_MODULE(TOPLaserCalibratorCollector);
26 
28 {
29  // Set module properties
30  setDescription("Collector module for the TOP ChannelT0 calibration and the quality monitoring using laser and pulser data");
32 
33  addParam("useReferencePulse", m_useReferencePulse, "Use the pulser as reference", bool(false));
34  addParam("refChannel", m_refChannel, "Channel to be used as reference");
35  addParam("refSlot", m_refSlot, "Slot to be used as reference");
36  addParam("pulserDeltaT", m_pulserDeltaT, "Approximate time difference between the two calpulses, in ns", float(21.8));
37  addParam("pulserDeltaTTolerance", m_pulserDeltaTTolerance, "Window around the nominal deltaT used to select a double pulse, in ns ",
38  float(2.));
39  addParam("storeMCTruth", m_storeMCTruth, " Store the TOPSimHits information instead of the TOPDigit one. \
40  If this option is used, useReferencePulse will be automatically set to false.");
41 
42 }
43 
44 
45 
47 {
48 
49  // Create the output trees
50 
51  auto hitTree = new TTree("hitTree", "hitTree");
52  hitTree->Branch<short>("channel", &m_channel);
53  hitTree->Branch<short>("slot", &m_slot);
54  hitTree->Branch<float>("hitTime", &m_hitTime);
55  hitTree->Branch<float>("dVdt", &m_dVdt);
56  hitTree->Branch<float>("refTime", &m_refTime);
57  hitTree->Branch<float>("amplitude", &m_amplitude);
58  hitTree->Branch<float>("width", &m_width);
59  hitTree->Branch<short>("sample", &m_sample);
60  hitTree->Branch<short>("window", &m_window);
61  hitTree->Branch<int>("event", &m_event);
62  hitTree->Branch<bool>("refTimeValid", &m_refTimeValid);
63 
64  registerObject<TTree>("hitTree", hitTree);
65 
66  m_TOPDigitArray.isRequired();
67 
68  if (m_storeMCTruth)
69  m_useReferencePulse = false;
70 
71  m_event = 0;
72 }
73 
74 
78 {
79  float refTimes[16] = {0.}; // Reference time for each slot
80  std::vector<bool> refTimesValid(16, true);
81 
82  // first loop over TOPDigits to find all the pairs of digits that satisfy the double-pulse conditions
83  if (m_useReferencePulse) {
84 
85  std::vector<float> calPulseTimes[16];
86  for (const auto& digit : m_TOPDigitArray) {
87  if (digit.getHitQuality() != TOPDigit::c_CalPulse or digit.getChannel() != m_refChannel)
88  continue; // remove photons and everything not on the ref channel
89  calPulseTimes[digit.getModuleID() - 1].push_back(digit.getTime());
90  }
91 
92  for (int i = 0; i < 16; i++) {
93  refTimesValid[i] = false;
94  auto& calTimes = calPulseTimes[i];
95  if (calTimes.size() < 2) continue;
96  std::sort(calTimes.begin(), calTimes.end());
97  for (unsigned k = 0; k < calTimes.size() - 1; k++) {
98  auto t1 = calTimes[k];
99  auto t2 = calTimes[k + 1];
100  if (fabs(fabs(t2 - t1) - m_pulserDeltaT) < m_pulserDeltaTTolerance) {
101  refTimes[i] = t1;
102  refTimesValid[i] = true;
103  break;
104  }
105  }
106  }
107 
108  }
109 
110  TTree* hitTree = getObjectPtr<TTree>("hitTree");
111 
112  // then fill the tree with good digits (photons and cal pulses)
113  for (const auto& digit : m_TOPDigitArray) {
114  if (digit.getHitQuality() == TOPDigit::c_Junk) continue; // remove the bad hits
115  m_channel = digit.getChannel();
116  m_slot = digit.getModuleID(); // this is 1-based
117  m_dVdt = 0.5 * TMath::Sqrt(-2.*TMath::Log(0.5)) * digit.getPulseHeight() / digit.getPulseWidth();
118 
119  if (m_refSlot > 0) {
120  m_refTime = refTimes[m_refSlot - 1];
121  m_refTimeValid = refTimesValid[m_refSlot - 1];
122  } else {
123  m_refTime = refTimes[m_slot - 1];
124  m_refTimeValid = refTimesValid[m_slot - 1];
125  }
126 
127  m_amplitude = digit.getPulseHeight() ;
128  m_width = digit.getPulseWidth();
129 
130  if (m_storeMCTruth) {
131  const auto* simHit = digit.getRelated<TOPSimHit>();
132  if (not simHit) continue; // no tree entry if MC truth doesn't exist
133  m_hitTime = simHit->getTime();
134  } else {
135  m_hitTime = digit.getTime() - m_refTime;
136  }
137 
138  m_window = -1;
139  const auto* rawDigit = digit.getRelated<TOPRawDigit>();
140  if (rawDigit) m_window = rawDigit->getASICWindow(); // window from which the feature is extracted
141  m_sample = digit.getModulo256Sample(); // sample number refered in TBC
142  hitTree->Fill();
143  }
144  m_event++;
145 }
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:208
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:80
TOPLaserCalibratorCollectorModule()
Constructor: Sets the description, the properties and the parameters of the module.
StoreArray< TOPDigit > m_TOPDigitArray
Required input array of TOPDigits.
bool m_storeMCTruth
Store the TOPSimHits information instead of the TOPDigit one.
float m_hitTime
Hit time with respect to the reference pulse (ns)
float m_pulserDeltaTTolerance
Window around the nominal deltaT used to select a double pulse, in ns.
void collect() override
Main mathod, called for each event.
bool m_useReferencePulse
Use the electronic pulser as reference.
bool m_refTimeValid
true when the time of the reference pulse is valid
float m_pulserDeltaT
Approximate time difference between the two calpulses, in ns.
Class to store unpacked raw data (hits in feature-extraction format) It provides also calculation of ...
Definition: TOPRawDigit.h:24
Class to store simulated hits of Cherenkov photons on PMT's input for digitization module (TOPDigitiz...
Definition: TOPSimHit.h:27
REG_MODULE(arichBtest)
Register the Module.
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:560
Abstract base class for different kinds of events.