Belle II Software development
KLMHitRateCounter.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/KLMHitRateCounter.h>
11
12using namespace Belle2::Background;
13
15{
16 // register collection(s) as optional, your detector might be excluded in DAQ
17 m_digits.isOptional();
20
21 // set branch address
22 std::string branches;
23 branches =
24 "moduleRates[" +
26 "]/F:averageRate/F:numEvents/I:valid/O";
27 tree->Branch("klm", &m_rates, branches.c_str());
28}
29
31{
32 m_buffer.clear();
33}
34
35void KLMHitRateCounter::accumulate(unsigned timeStamp)
36{
37 // check if data are available
38 if (not m_digits.isValid())
39 return;
40
41 // get buffer element
42 auto& rates = m_buffer[timeStamp];
43
44 // increment event counter
45 rates.numEvents++;
46
47 // accumulate hits
48 for (const KLMDigit& klmDigit : m_digits) {
49 if (!klmDigit.inRPC() && !klmDigit.isGood())
50 continue;
51 uint16_t moduleNumber =
53 klmDigit.getSubdetector(), klmDigit.getSection(), klmDigit.getSector(),
54 klmDigit.getLayer());
55 int module = m_ModuleArrayIndex->getIndex(moduleNumber);
56 if (module >= 0 and module < KLMElementNumbers::getTotalModuleNumber()) {
57 rates.moduleRates[module] += 1;
58 } else {
59 B2ERROR("KLMHitRateCounter: module number out of range"
60 << LogVar("module", module));
61 }
62 rates.averageRate += 1;
63 }
64
65 // set flag to true to indicate the rates are valid
66 rates.valid = true;
67}
68
69void KLMHitRateCounter::normalize(unsigned timeStamp)
70{
71 // copy buffer element
72 m_rates = m_buffer[timeStamp];
73
74 if (not m_rates.valid) return;
75
76 // normalize
78
79 /* Normalize the hit rate per 1 strip. */
80 for (int i = 0; i < KLMElementNumbers::getTotalModuleNumber(); ++i) {
81 uint16_t module = m_ModuleArrayIndex->getNumber(i);
82 int activeStrips = m_ChannelStatus->getActiveStripsInModule(module);
83 if (activeStrips == 0)
84 m_rates.moduleRates[i] = 0;
85 else
86 m_rates.moduleRates[i] /= activeStrips;
87 }
88}
std::map< unsigned, TreeStruct > m_buffer
Buffer.
StoreArray< KLMDigit > m_digits
KLM digits.
DBObjPtr< KLMChannelStatus > m_ChannelStatus
KLM channel status.
virtual void initialize(TTree *tree) override
Class initializer.
const KLMElementNumbers * m_ElementNumbers
KLM element numbers.
virtual void accumulate(unsigned timeStamp) override
Accumulate hits.
const KLMModuleArrayIndex * m_ModuleArrayIndex
KLM module array index.
virtual void normalize(unsigned timeStamp) override
Normalize accumulated hits (i.e.
virtual void clear() override
Clear time-stamp buffer to prepare for 'accumulate'.
KLM digit (class representing a digitized hit in RPCs or scintillators).
Definition: KLMDigit.h:29
uint16_t getNumber(uint16_t index) const
Get element number.
uint16_t getIndex(uint16_t number) const
Get element index.
static const KLMElementNumbers & Instance()
Instantiation.
static constexpr int getTotalModuleNumber()
Get total number of modules.
KLMModuleNumber moduleNumber(int subdetector, int section, int sector, int layer) const
Get module number.
static const KLMModuleArrayIndex & Instance()
Instantiation.
Class to store variables with their name which were sent to the logging service.
float moduleRates[KLMElementNumbers::getTotalModuleNumber()]
Hit rates in each module.
void normalize()
Normalize accumulated hits to single event.