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