Belle II Software  release-08-01-10
KLMHitRateCounter Class Reference

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

#include <KLMHitRateCounter.h>

Inheritance diagram for KLMHitRateCounter:
Collaboration diagram for KLMHitRateCounter:

Classes

struct  TreeStruct
 Tree data structure. More...
 

Public Member Functions

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

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.

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.

◆ initialize()

void initialize ( TTree *  tree)
overridevirtual

Class initializer.

Parameters
[in,out]treeData tree.

Implements HitRateBase.

Definition at line 14 of file KLMHitRateCounter.cc.

◆ 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.


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