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