Belle II Software development
KLMHitRateCounter Class Reference

Class for monitoring beam background hit rates of EKLM. More...

#include <KLMHitRateCounter.h>

Inheritance diagram for KLMHitRateCounter:
HitRateBase

Classes

struct  TreeStruct
 Tree data structure. More...
 

Public Member Functions

 KLMHitRateCounter ()
 Constructor.
 
virtual void initialize (TTree *tree) override
 Class initializer.
 
virtual void clear () override
 Clear time-stamp buffer to prepare for 'accumulate'.
 
virtual void accumulate (unsigned timeStamp) override
 Accumulate hits.
 
virtual void normalize (unsigned timeStamp) override
 Normalize accumulated hits (i.e.
 

Private Attributes

TreeStruct m_rates
 Tree data.
 
std::map< unsigned, TreeStructm_buffer
 Buffer.
 
StoreArray< KLMDigitm_digits
 KLM digits.
 
const KLMElementNumbersm_ElementNumbers = nullptr
 KLM element numbers.
 
const KLMModuleArrayIndexm_ModuleArrayIndex = nullptr
 KLM module array index.
 
DBObjPtr< KLMChannelStatusm_ChannelStatus
 KLM channel status.
 

Detailed Description

Class for monitoring beam background hit rates of EKLM.

Definition at line 37 of file KLMHitRateCounter.h.

Constructor & Destructor Documentation

◆ KLMHitRateCounter()

KLMHitRateCounter ( )
inline

Constructor.

Definition at line 75 of file KLMHitRateCounter.h.

75{};

Member Function Documentation

◆ accumulate()

void accumulate ( unsigned  timeStamp)
overridevirtual

Accumulate hits.

Parameters
[in]timeStampTime stamp.

Implements HitRateBase.

Definition at line 35 of file KLMHitRateCounter.cc.

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}
std::map< unsigned, TreeStruct > m_buffer
Buffer.
StoreArray< KLMDigit > m_digits
KLM digits.
const KLMElementNumbers * m_ElementNumbers
KLM element numbers.
const KLMModuleArrayIndex * m_ModuleArrayIndex
KLM module array index.
uint16_t getIndex(uint16_t number) const
Get element index.
static constexpr int getTotalModuleNumber()
Get total number of modules.
KLMModuleNumber moduleNumber(int subdetector, int section, int sector, int layer) const
Get module number.
Class to store variables with their name which were sent to the logging service.

◆ clear()

void clear ( )
overridevirtual

Clear time-stamp buffer to prepare for 'accumulate'.

Implements HitRateBase.

Definition at line 30 of file KLMHitRateCounter.cc.

31{
32 m_buffer.clear();
33}

◆ initialize()

void initialize ( TTree *  tree)
overridevirtual

Class initializer.

Parameters
[in,out]treeData tree.

Implements HitRateBase.

Definition at line 14 of file KLMHitRateCounter.cc.

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}
static const KLMElementNumbers & Instance()
Instantiation.
static const KLMModuleArrayIndex & Instance()
Instantiation.

◆ normalize()

void normalize ( unsigned  timeStamp)
overridevirtual

Normalize accumulated hits (i.e.

transform to rates).

Parameters
[in]timeStampTime stamp.

Implements HitRateBase.

Definition at line 69 of file KLMHitRateCounter.cc.

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}
DBObjPtr< KLMChannelStatus > m_ChannelStatus
KLM channel status.
uint16_t getNumber(uint16_t index) const
Get element number.
float moduleRates[KLMElementNumbers::getTotalModuleNumber()]
Hit rates in each module.
void normalize()
Normalize accumulated hits to single event.

Member Data Documentation

◆ m_buffer

std::map<unsigned, TreeStruct> m_buffer
private

Buffer.

Definition at line 106 of file KLMHitRateCounter.h.

◆ m_ChannelStatus

DBObjPtr<KLMChannelStatus> m_ChannelStatus
private

KLM channel status.

Definition at line 118 of file KLMHitRateCounter.h.

◆ m_digits

StoreArray<KLMDigit> m_digits
private

KLM digits.

Definition at line 109 of file KLMHitRateCounter.h.

◆ m_ElementNumbers

const KLMElementNumbers* m_ElementNumbers = nullptr
private

KLM element numbers.

Definition at line 112 of file KLMHitRateCounter.h.

◆ m_ModuleArrayIndex

const KLMModuleArrayIndex* m_ModuleArrayIndex = nullptr
private

KLM module array index.

Definition at line 115 of file KLMHitRateCounter.h.

◆ m_rates

TreeStruct m_rates
private

Tree data.

Definition at line 103 of file KLMHitRateCounter.h.


The documentation for this class was generated from the following files: