Belle II Software  release-06-00-14
TOPHitRateCounter.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 include
10 #include <background/modules/BeamBkgHitRateMonitor/TOPHitRateCounter.h>
11 #include <top/geometry/TOPGeometryPar.h>
12 
13 // framework aux
14 #include <framework/gearbox/Unit.h>
15 #include <framework/logging/Logger.h>
16 
17 using namespace std;
18 
19 namespace Belle2 {
24  namespace Background {
25 
26  void TOPHitRateCounter::initialize(TTree* tree)
27  {
28  // register collection(s) as optional, your detector might be excluded in DAQ
29  m_digits.isOptional();
30 
31  // set branch address
32  tree->Branch("top", &m_rates, "slotRates[16]/F:averageRate/F:numEvents/I:valid/O");
33 
34  // create control histograms
35  m_hits = new TH1F("top_hits", "time distribution of hits; digit.time [ns]",
36  1000, -500, 500);
37  m_hitsInWindow = new TH1F("top_hitsInWindow",
38  "time distribution of hits; digit.time [ns]",
39  100, m_timeOffset - m_timeWindow / 2,
40  m_timeOffset + m_timeWindow / 2);
41 
42  // check parameters
43  if (m_timeWindow <= 0) B2FATAL("invalid time window for TOP: " << m_timeWindow);
44 
45  // set fractions of active channels
46  setActiveFractions();
47  }
48 
49  void TOPHitRateCounter::clear()
50  {
51  m_buffer.clear();
52  }
53 
54  void TOPHitRateCounter::accumulate(unsigned timeStamp)
55  {
56  // check if data are available
57  if (not m_digits.isValid()) return;
58 
59  // get buffer element
60  auto& rates = m_buffer[timeStamp];
61 
62  // increment event counter
63  rates.numEvents++;
64 
65  // accumulate hits (weighted by efficiency correction)
66  const auto* topgp = TOP::TOPGeometryPar::Instance();
67  for (const auto& digit : m_digits) {
68  if (digit.getHitQuality() != TOPDigit::c_Good) continue;
69  m_hits->Fill(digit.getTime());
70 
71  if (fabs(digit.getTime() - m_timeOffset) > m_timeWindow / 2) continue;
72  m_hitsInWindow->Fill(digit.getTime());
73 
74  auto effi = topgp->getRelativePixelEfficiency(digit.getModuleID(),
75  digit.getPixelID());
76  float wt = std::min(1.0 / effi, 10.0);
77  rates.slotRates[digit.getModuleID() - 1] += wt;
78  rates.averageRate += wt;
79  }
80 
81  // set flag to true to indicate the rates are valid
82  rates.valid = true;
83 
84  }
85 
86  void TOPHitRateCounter::normalize(unsigned timeStamp)
87  {
88  // copy buffer element
89  m_rates = m_buffer[timeStamp];
90 
91  if (not m_rates.valid) return;
92 
93  // normalize
94  m_rates.normalize();
95 
96  // convert rates to [MHz/PMT]
97  for (auto& slotRate : m_rates.slotRates) {
98  slotRate /= m_timeWindow / Unit::us * 32; // 32 PMT's per slot
99  }
100  m_rates.averageRate /= m_timeWindow / Unit::us * 32 * 16; // 16 slots
101 
102  // check if DB object has changed and if so set fractions of active channels
103  if (m_channelMask.hasChanged()) setActiveFractions();
104 
105  // correct rates for masked-out channels
106  for (int m = 0; m < 16; m++) {
107  double fraction = m_activeFractions[m];
108  if (fraction > 0) {
109  m_rates.slotRates[m] /= fraction;
110  } else {
111  m_rates.slotRates[m] = 0;
112  }
113  }
114  m_rates.averageRate /= m_activeTotal; // we don't expect full TOP to be masked-out
115 
116  }
117 
118  void TOPHitRateCounter::setActiveFractions()
119  {
120  if (not m_channelMask.isValid()) {
121  for (auto& fraction : m_activeFractions) fraction = 1;
122  m_activeTotal = 1;
123  B2WARNING("TOPHitRateCounter: no valid channel mask - active fractions set to 1");
124  return;
125  }
126 
127  for (int m = 0; m < 16; m++) {
128  m_activeFractions[m] = m_channelMask->getActiveFraction(m + 1);
129  }
130  m_activeTotal = m_channelMask->getActiveFraction();
131 
132  }
133 
134 
135  } // Background namespace
137 } // Belle2 namespace
138 
Abstract base class for different kinds of events.