Belle II Software development
ARICHHitRateCounter Class Reference

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

#include <ARICHHitRateCounter.h>

Inheritance diagram for ARICHHitRateCounter:
HitRateBase

Classes

struct  TreeStruct
 tree structure More...
 

Public Member Functions

 ARICHHitRateCounter ()
 Constructor.
 
virtual void initialize (TTree *tree) override
 Class initializer: set branch addresses and other staf.
 
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 (e.g.
 

Private Member Functions

void setActiveHapds ()
 Sets number of active hapds in each segment.
 

Private Attributes

TreeStruct m_rates
 tree variables
 
std::map< unsigned, TreeStructm_buffer
 average rates in time stamps
 
StoreArray< ARICHHitm_hits
 collection of digits
 
DBObjPtr< ARICHChannelMaskm_channelMask
 channel mask
 
DBObjPtr< ARICHModulesInfom_modulesInfo
 HAPD modules info.
 
double m_activeHapds [18] = {0}
 number of active HAPDS in each segment
 
double m_activeTotal = 0
 total number of active HAPDS
 
int m_segmentMap [420] = {0}
 mapping from module ID to segments
 

Detailed Description

Class for monitoring beam background hit rates of ARICH.

Definition at line 31 of file ARICHHitRateCounter.h.

Constructor & Destructor Documentation

◆ ARICHHitRateCounter()

ARICHHitRateCounter ( )
inline

Constructor.

Definition at line 60 of file ARICHHitRateCounter.h.

61 {}

Member Function Documentation

◆ accumulate()

void accumulate ( unsigned  timeStamp)
overridevirtual

Accumulate hits.

Parameters
timeStamptime stamp

Implements HitRateBase.

Definition at line 53 of file ARICHHitRateCounter.cc.

54 {
55 // check if data are available
56 if (not m_hits.isValid()) return;
57
58 // get buffer element
59 auto& rates = m_buffer[timeStamp];
60
61 // increment event counter
62 rates.numEvents++;
63
64 // count and weight hits accoring to channel efficiecny
65 for (const auto& hit : m_hits) {
66 if (hit.getModule() < 1 || hit.getModule() > 420) continue;
67 auto effi = m_modulesInfo->getChannelQE(hit.getModule(), hit.getChannel());
68 float wt = std::min(1.0 / effi, 100.);
69 rates.segmentRates[m_segmentMap[hit.getModule() - 1]] += wt;
70 rates.averageRate += wt;
71 }
72
73 // set flag to true to indicate the rates are valid
74 rates.valid = true;
75
76 }
std::map< unsigned, TreeStruct > m_buffer
average rates in time stamps
int m_segmentMap[420]
mapping from module ID to segments
StoreArray< ARICHHit > m_hits
collection of digits
DBObjPtr< ARICHModulesInfo > m_modulesInfo
HAPD modules info.

◆ clear()

void clear ( )
overridevirtual

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

Implements HitRateBase.

Definition at line 48 of file ARICHHitRateCounter.cc.

49 {
50 m_buffer.clear();
51 }

◆ initialize()

void initialize ( TTree *  tree)
overridevirtual

Class initializer: set branch addresses and other staf.

Parameters
treea valid TTree pointer

Implements HitRateBase.

Definition at line 24 of file ARICHHitRateCounter.cc.

25 {
26 // register collection(s) as optional, your detector might be excluded in DAQ
27 m_hits.isOptional();
28
29 // set branch address
30 tree->Branch("arich", &m_rates, "segmentRates[18]/F:averageRate/F:numEvents/I:valid/O");
31
32 // make map of modules to 18 segments
33 // (first hapd ring is 6 segments (by sector), for the rest, each segment merges 3 hapd rings (again by sector))
34 const int nModInRing[7] = {0, 42, 90, 144, 204, 270, 342};
35 int iRing = 0;
36 for (int i = 0; i < 420; i++) {
37 if (i == nModInRing[iRing + 1]) iRing++;
38 int segment = (i - nModInRing[iRing]) / (7 + iRing);
39 if (iRing > 0) segment += 6;
40 if (iRing > 3) segment += 6;
41 m_segmentMap[i] = segment;
42 }
43
44 // set fractions of active channels
46 }
void setActiveHapds()
Sets number of active hapds in each segment.

◆ normalize()

void normalize ( unsigned  timeStamp)
overridevirtual

Normalize accumulated hits (e.g.

transform to rates)

Parameters
timeStamptime stamp

Implements HitRateBase.

Definition at line 78 of file ARICHHitRateCounter.cc.

79 {
80 // copy buffer element
81 m_rates = m_buffer[timeStamp];
82
83 if (not m_rates.valid) return;
84
85 // normalize
87
88 // correct rates for masked-out channels
89 if (m_channelMask.hasChanged()) setActiveHapds();
90
91 for (int iSegment = 0; iSegment < 18; iSegment++) {
92 double nHapds = m_activeHapds[iSegment];
93 if (nHapds > 0) m_rates.segmentRates[iSegment] /= nHapds;
94 else m_rates.segmentRates[iSegment] = 0;
95 }
97 }
double m_activeTotal
total number of active HAPDS
DBObjPtr< ARICHChannelMask > m_channelMask
channel mask
double m_activeHapds[18]
number of active HAPDS in each segment
float segmentRates[18]
hit rates per HAPD [Hz] for 18 segments of arich
float averageRate
total detector average hit rate per HAPD [Hz]
void normalize()
normalize accumulated hits to single event

◆ setActiveHapds()

void setActiveHapds ( )
private

Sets number of active hapds in each segment.

Definition at line 99 of file ARICHHitRateCounter.cc.

100 {
101 for (auto& nactive : m_activeHapds) nactive = 0;
102
103 if (not m_channelMask.isValid()) {
104 for (unsigned imod = 1; imod < 421; imod++) m_activeHapds[m_segmentMap[imod - 1]] += 1.;
105 m_activeTotal = 420;
106 B2WARNING("ARICHHitRateCounter: no valid channel mask - all HAPDs set to active");
107 return;
108 }
109
110 int nactiveTotal = 0;
111 for (unsigned imod = 1; imod < 421; imod++) {
112 int nactive = 0;
113 for (unsigned ichn = 0; ichn < 144; ichn++) {
114 if (m_channelMask->isActive(imod, ichn)) nactive++;
115 }
116 nactiveTotal += nactive;
117 m_activeHapds[m_segmentMap[imod - 1]] += (float)nactive / 144.;
118 }
119 m_activeTotal = (float)nactiveTotal / 144.;
120 }

Member Data Documentation

◆ m_activeHapds

double m_activeHapds[18] = {0}
private

number of active HAPDS in each segment

Definition at line 107 of file ARICHHitRateCounter.h.

◆ m_activeTotal

double m_activeTotal = 0
private

total number of active HAPDS

Definition at line 108 of file ARICHHitRateCounter.h.

◆ m_buffer

std::map<unsigned, TreeStruct> m_buffer
private

average rates in time stamps

Definition at line 97 of file ARICHHitRateCounter.h.

◆ m_channelMask

DBObjPtr<ARICHChannelMask> m_channelMask
private

channel mask

Definition at line 103 of file ARICHHitRateCounter.h.

◆ m_hits

StoreArray<ARICHHit> m_hits
private

collection of digits

Definition at line 100 of file ARICHHitRateCounter.h.

◆ m_modulesInfo

DBObjPtr<ARICHModulesInfo> m_modulesInfo
private

HAPD modules info.

Definition at line 104 of file ARICHHitRateCounter.h.

◆ m_rates

TreeStruct m_rates
private

tree variables

Definition at line 94 of file ARICHHitRateCounter.h.

◆ m_segmentMap

int m_segmentMap[420] = {0}
private

mapping from module ID to segments

Definition at line 109 of file ARICHHitRateCounter.h.


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