Belle II Software  release-05-02-19
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 47 of file KLMHitRateCounter.h.

Member Function Documentation

◆ accumulate()

void accumulate ( unsigned  timeStamp)
overridevirtual

Accumulate hits.

Parameters
[in]timeStampTime stamp.

Implements HitRateBase.

Definition at line 37 of file KLMHitRateCounter.cc.

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 }

◆ initialize()

void initialize ( TTree *  tree)
overridevirtual

Class initializer.

Parameters
[in,out]treeData tree.

Implements HitRateBase.

Definition at line 16 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 71 of file KLMHitRateCounter.cc.


The documentation for this class was generated from the following files:
Belle2::Background::KLMHitRateCounter::m_digits
StoreArray< KLMDigit > m_digits
KLM digits.
Definition: KLMHitRateCounter.h:119
Belle2::KLMElementNumbers::getTotalModuleNumber
static constexpr int getTotalModuleNumber()
Get total number of modules.
Definition: KLMElementNumbers.h:251
Belle2::Background::KLMHitRateCounter::m_ElementNumbers
const KLMElementNumbers * m_ElementNumbers
KLM element numbers.
Definition: KLMHitRateCounter.h:122
LogVar
Class to store variables with their name which were sent to the logging service.
Definition: LogVariableStream.h:24
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