Belle II Software development
CDCHitRateCounter Class Reference

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

#include <CDCHitRateCounter.h>

Inheritance diagram for CDCHitRateCounter:
HitRateBase

Classes

struct  TreeStruct
 tree structure More...
 

Public Member Functions

 CDCHitRateCounter (const int timeWindowLowerEdge_smallCell, const int timeWindowUpperEdge_smallCell, const int timeWindowLowerEdge_normalCell, const int timeWindowUpperEdge_normalCell, const bool enableBadWireTreatment, const bool enableBackgroundHitFilter, const bool enableMarkBackgroundHit)
 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

bool isInTimeWindow (const int SL, const short tdc)
 return true if the hit is in the given time window
 
void countActiveWires_countAll ()
 set m_nActiveWireInTotal, m_nActiveWireInLayer[] and m_nActiveWireInSuperLayer[].
 
void countActiveWires ()
 set m_nActiveWireInTotal, m_nActiveWireInLayer[] and m_nActiveWireInSuperLayer[].
 
unsigned short getIPhiBin (unsigned short iSL, unsigned short iWireInLayer)
 get the bin ID of the division.
 

Private Attributes

TreeStruct m_rates
 tree variables
 
std::map< unsigned, TreeStructm_buffer
 average rates in time stamps
 
StoreArray< CDCHitm_digits
 collection of digits
 
const int m_timeWindowLowerEdge_smallCell
 lower edge of the time window for small cells [tdc count = 0.982536 ns]
 
const int m_timeWindowUpperEdge_smallCell
 upper edge of the time window for small cells [tdc count = 0.982536 ns]
 
const int m_timeWindowLowerEdge_normalCell
 lower edge of the time window for normal cells [tdc count = 0.982536 ns]
 
const int m_timeWindowUpperEdge_normalCell
 upper edge of the time window for normal cells [tdc count = 0.982536 ns]
 
const bool m_enableBadWireTreatment
 flag to enable the bad wire treatment.
 
const bool m_enableBackgroundHitFilter
 flag to enable the CDC background hit (crosstakl, noise) filter.
 
const bool m_enableMarkBackgroundHit
 flag to enable to mark background flag on CDCHit (set 0x100 bit for CDCHit::m_status).
 
int m_nActiveWireInTotal
 the number of wires used in this hit-rate calculation in the whole CDC
 
int m_nActiveWireInSuperLayer [f_nSuperLayer] = {0}
 the number of wires used in this hit-rate calculation in each suler layer
 
int m_nActiveWireInLayer [f_nLayer] = {0}
 the number of wires used in this hit-rate calculation in each layer
 
int m_nActiveWireInLayerPhi [f_nLayer][f_nPhiDivision] = {{0}}
 the number of wires used in this hit-rate calculation in each phi bin in each layer
 

Static Private Attributes

static const int f_nSuperLayer = 9
 the number of super layers
 
static const int f_nLayer = 56
 the number of layers
 
static const int f_nPhiDivision = 8
 the number of division in phi
 

Detailed Description

Class for monitoring beam background hit rates of CDC.

Definition at line 28 of file CDCHitRateCounter.h.

Constructor & Destructor Documentation

◆ CDCHitRateCounter()

CDCHitRateCounter ( const int  timeWindowLowerEdge_smallCell,
const int  timeWindowUpperEdge_smallCell,
const int  timeWindowLowerEdge_normalCell,
const int  timeWindowUpperEdge_normalCell,
const bool  enableBadWireTreatment,
const bool  enableBackgroundHitFilter,
const bool  enableMarkBackgroundHit 
)
inline

Constructor.

Parameters
timeWindowLowerEdge_smallCelllower edge of the timing window for the layers with small cells (SL0)
timeWindowUpperEdge_smallCellupper edge of the timing window for the layers with small cells (SL0)
timeWindowLowerEdge_normalCelllower edge of the timing window for the layers with normal cells (SL1-8)
timeWindowUpperEdge_normalCellupper edge of the timing window for the layers with normal cells (SL1-8)
enableBadWireTreatmentflag to enable the bad wire treatment. default: true
enableBackgroundHitFilterflag to enable the CDC background hit (crosstakl, noise) filter. default: true
enableMarkBackgroundHitflag to enable to mark background flag on CDCHit (set 0x100 bit for CDCHit::m_status). default: false

Definition at line 108 of file CDCHitRateCounter.h.

114 :
115 m_timeWindowLowerEdge_smallCell(timeWindowLowerEdge_smallCell),
116 m_timeWindowUpperEdge_smallCell(timeWindowUpperEdge_smallCell),
117 m_timeWindowLowerEdge_normalCell(timeWindowLowerEdge_normalCell),
118 m_timeWindowUpperEdge_normalCell(timeWindowUpperEdge_normalCell),
119 m_enableBadWireTreatment(enableBadWireTreatment),
120 m_enableBackgroundHitFilter(enableBackgroundHitFilter),
121 m_enableMarkBackgroundHit(enableMarkBackgroundHit)
122 {
125 B2FATAL("invalid setting of CDC time window");
126 }
127 }
const bool m_enableBackgroundHitFilter
flag to enable the CDC background hit (crosstakl, noise) filter.
const int m_timeWindowLowerEdge_normalCell
lower edge of the time window for normal cells [tdc count = 0.982536 ns]
const int m_timeWindowUpperEdge_smallCell
upper edge of the time window for small cells [tdc count = 0.982536 ns]
const bool m_enableBadWireTreatment
flag to enable the bad wire treatment.
const int m_timeWindowUpperEdge_normalCell
upper edge of the time window for normal cells [tdc count = 0.982536 ns]
const bool m_enableMarkBackgroundHit
flag to enable to mark background flag on CDCHit (set 0x100 bit for CDCHit::m_status).
const int m_timeWindowLowerEdge_smallCell
lower edge of the time window for small cells [tdc count = 0.982536 ns]

Member Function Documentation

◆ accumulate()

void accumulate ( unsigned  timeStamp)
overridevirtual

Accumulate hits.

Parameters
timeStamptime stamp

Implements HitRateBase.

Definition at line 69 of file CDCHitRateCounter.cc.

70 {
71 // check if data are available
72 if (not m_digits.isValid()) return;
73
74 // get buffer element
75 auto& rates = m_buffer[timeStamp];
76
77 // increment event counter
78 rates.numEvents++;
79
81 std::map<const CDCHit*, bool> CDCHitToBackgroundFlag;
83 StoreWrappedObjPtr<std::vector<CDCWireHit>> storeVector("CDCWireHitVector");
84 if (not storeVector) {
85 B2FATAL("CDCWireHitVector is unaccessible in DataStore."
86 "Need TFCDC_WireHitParameter module before.");
87 }
88 const std::vector<CDCWireHit>& cdcWireHitVector = *storeVector;
89 for (const auto& cdcWireHit : cdcWireHitVector) {
90 const CDCHit* cdcHit = cdcWireHit.getHit();
91 CDCHitToBackgroundFlag[cdcHit] = cdcWireHit->hasBackgroundFlag();
92 }
93 }
94
95 CDC::CDCGeometryPar& geometryPar = CDC::CDCGeometryPar::Instance();
96
104 for (CDCHit& hit : m_digits) {
105 const WireID wireID(hit.getID());
106 if (m_enableBadWireTreatment && geometryPar.isBadWire(wireID))
107 continue;
108
109 const unsigned short iSuperLayer = hit.getISuperLayer();
110 const unsigned short iLayer = hit.getICLayer();
111 const unsigned short iWireInLayer = hit.getIWire();
112 const unsigned short iPhiBin = getIPhiBin(iSuperLayer, iWireInLayer);
113 const unsigned short adc = hit.getADCCount();
114 const short tdc = hit.getTDCCount();
115
117 if (CDCHitToBackgroundFlag[&hit]) {
119 unsigned short newStatus = (hit.getStatus() | 0x100);
120 hit.setStatus(newStatus);
121 }
122 continue;
123 }
124 } else {
125 if (iSuperLayer == 0 && adc < 15)
126 continue;
127 if (iSuperLayer != 0 && adc < 18)
128 continue;
129 }
130
131 if (not isInTimeWindow(iSuperLayer, tdc))
132 continue;
133
134
135 rates.averageRate += 1;
136 rates.superLayerHitRate[iSuperLayer] += 1;
137 rates.layerHitRate[iLayer] += 1;
138 rates.layerPhiHitRate[iLayer][iPhiBin] += 1;
139 }
140
141 // set flag to true to indicate the rates are valid
142 rates.valid = true;
143 }
std::map< unsigned, TreeStruct > m_buffer
average rates in time stamps
bool isInTimeWindow(const int SL, const short tdc)
return true if the hit is in the given time window
StoreArray< CDCHit > m_digits
collection of digits
unsigned short getIPhiBin(unsigned short iSL, unsigned short iWireInLayer)
get the bin ID of the division.
static CDCGeometryPar & Instance(const CDCGeometry *=nullptr)
Static method to get a reference to the CDCGeometryPar instance.

◆ clear()

void clear ( )
overridevirtual

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

Implements HitRateBase.

Definition at line 64 of file CDCHitRateCounter.cc.

65 {
66 m_buffer.clear();
67 }

◆ countActiveWires()

void countActiveWires ( )
private

set m_nActiveWireInTotal, m_nActiveWireInLayer[] and m_nActiveWireInSuperLayer[].

called in initialize function. count the number of wires excluding dead wires (bad channels).

end i loop

end iL loop

end iSL loop

Definition at line 228 of file CDCHitRateCounter.cc.

229 {
230 static const unsigned short nlayer_in_SL[f_nSuperLayer] = { 8, 6, 6,
231 6, 6, 6,
232 6, 6, 6
233 };
234 static const unsigned short nwire_in_layer[f_nSuperLayer] = { 160, 160, 192,
235 224, 256, 288,
236 320, 352, 384
237 };
238
239 CDC::CDCGeometryPar& geometryPar = CDC::CDCGeometryPar::Instance();
240
242 int contLayerID = 0;
243
244 for (int iSL = 0 ; iSL < f_nSuperLayer ; ++iSL) {
246 for (int iL = 0 ; iL < nlayer_in_SL[iSL] ; ++iL) {
247 m_nActiveWireInLayer[contLayerID] = 0;
248 for (int iPhi = 0 ; iPhi < f_nPhiDivision ; ++iPhi)
249 m_nActiveWireInLayerPhi[contLayerID][iPhi] = 0;
250
251 for (int i = 0 ; i < nwire_in_layer[iSL] ; ++i) {
252 WireID wireID(iSL, iL, i);
253 if (not geometryPar.isBadWire(wireID)) {
254 ++m_nActiveWireInLayerPhi[contLayerID][getIPhiBin(iSL, i)];
255 ++m_nActiveWireInLayer[contLayerID];
256 }
257 }
259 ++contLayerID;
260 }
262 }
263
264
265 //B2INFO("CDC, # of Active wires / # of total wires");
266 std::cout << "CDC, # of Active wires / # of total wires" << std::endl;
267 int contLayerID_2 = 0;
268 for (int iSL = 0 ; iSL < f_nSuperLayer ; ++iSL) {
269 for (int iL = 0 ; iL < nlayer_in_SL[iSL] ; ++iL) {
270 //B2INFO("Layer "<< contLayerID_2 << ": "
271 std::cout << "Layer " << contLayerID_2 << ": "
272 << m_nActiveWireInLayer[contLayerID_2] << " / "
273 //<< nwire_in_layer[iSL] );
274 << nwire_in_layer[iSL] << std::endl;
275 ++contLayerID_2;
276 }
277 }
278
279 }
int m_nActiveWireInLayer[f_nLayer]
the number of wires used in this hit-rate calculation in each layer
int m_nActiveWireInTotal
the number of wires used in this hit-rate calculation in the whole CDC
int m_nActiveWireInLayerPhi[f_nLayer][f_nPhiDivision]
the number of wires used in this hit-rate calculation in each phi bin in each layer
static const int f_nSuperLayer
the number of super layers
static const int f_nPhiDivision
the number of division in phi
int m_nActiveWireInSuperLayer[f_nSuperLayer]
the number of wires used in this hit-rate calculation in each suler layer

◆ countActiveWires_countAll()

void countActiveWires_countAll ( )
private

set m_nActiveWireInTotal, m_nActiveWireInLayer[] and m_nActiveWireInSuperLayer[].

called in initialize function. count the number of all the wires including dead wires (bad channels).

Definition at line 200 of file CDCHitRateCounter.cc.

201 {
202 static const unsigned short nlayer_in_SL[f_nSuperLayer] = { 8, 6, 6,
203 6, 6, 6,
204 6, 6, 6
205 };
206 static const unsigned short nwire_in_layer[f_nSuperLayer] = { 160, 160, 192,
207 224, 256, 288,
208 320, 352, 384
209 };
210
212 int contLayerID = 0;
213 for (int iSL = 0 ; iSL < f_nSuperLayer ; ++iSL) {
215 for (int i = 0 ; i < nlayer_in_SL[iSL] ; ++i) {
216 m_nActiveWireInTotal += nwire_in_layer[iSL];
217 m_nActiveWireInSuperLayer[iSL] += nwire_in_layer[iSL];
218 m_nActiveWireInLayer[contLayerID] = nwire_in_layer[iSL];
219 for (int j = 0 ; j < f_nPhiDivision ; ++j)
220 m_nActiveWireInLayerPhi[contLayerID][j] = nwire_in_layer[iSL] / f_nPhiDivision;
221 ++contLayerID;
222 }
223 }
224 }

◆ getIPhiBin()

unsigned short getIPhiBin ( unsigned short  iSL,
unsigned short  iWireInLayer 
)
private

get the bin ID of the division.

Definition at line 282 of file CDCHitRateCounter.cc.

283 {
284 static const unsigned short nwire_in_layer[f_nSuperLayer] = { 160, 160, 192,
285 224, 256, 288,
286 320, 352, 384
287 };
288 return iWireInLayer / (nwire_in_layer[iSL] / f_nPhiDivision);
289 }

◆ initialize()

void initialize ( TTree *  tree)
overridevirtual

Class initializer: set branch addresses and other staf.

Parameters
treea valid TTree pointer

Implements HitRateBase.

Definition at line 34 of file CDCHitRateCounter.cc.

35 {
36 // register collection(s) as optional, your detector might be excluded in DAQ
37 m_digits.isOptional();
38
39 // set branch address
40 //tree->Branch("cdc", &m_rates, "averageRate/F:numEvents/I:valid/O");
41 stringstream leafList;
42 leafList
43 << "averageRate/F:"
44 << "superLayerHitRate[" << f_nSuperLayer << "]/F:"
45 << "layerHitRate[" << f_nLayer << "]/F:"
46 << "layerPhiHitRate[" << f_nLayer << "][" << f_nPhiDivision << "]/F:"
47 << "timeWindowForSmallCell/I:"
48 << "timeWindowForNormalCell/I:"
49 << "nActiveWireInTotal/I:"
50 << "nActiveWireInSuperLayer[" << f_nSuperLayer << "]/I:"
51 << "nActiveWireInLayer[" << f_nLayer << "]/I:"
52 << "nActiveWireInLayerPhi[" << f_nLayer << "][" << f_nPhiDivision << "]/I:"
53 << "numEvents/I:"
54 << "valid/O";
55 tree->Branch("cdc", &m_rates, leafList.str().c_str());
56
59 } else {
61 }
62 }
static const int f_nLayer
the number of layers
void countActiveWires_countAll()
set m_nActiveWireInTotal, m_nActiveWireInLayer[] and m_nActiveWireInSuperLayer[].
void countActiveWires()
set m_nActiveWireInTotal, m_nActiveWireInLayer[] and m_nActiveWireInSuperLayer[].

◆ isInTimeWindow()

bool isInTimeWindow ( const int  SL,
const short  tdc 
)
inlineprivate

return true if the hit is in the given time window

Parameters
SLsuper layer ID which the wire of the hit belongs to
tdcTDC value of the hit

Definition at line 175 of file CDCHitRateCounter.h.

176 {
177 if (SL == 0) {
179 } else {
181 }
182 }

◆ normalize()

void normalize ( unsigned  timeStamp)
overridevirtual

Normalize accumulated hits (e.g.

transform to rates)

Parameters
timeStamptime stamp

copy nActiveWire

Implements HitRateBase.

Definition at line 145 of file CDCHitRateCounter.cc.

146 {
147 // copy buffer element
148 m_rates = m_buffer[timeStamp];
149
150 if (not m_rates.valid) return;
151
154
155 // normalize : hit rate in the time stamp in kHz
157
160 if (m_nActiveWireInTotal == 0) {
162 } else {
164 }
166 for (int iSL = 0 ; iSL < f_nSuperLayer ; ++iSL)
167 if (m_nActiveWireInSuperLayer[iSL] == 0) {
168 m_rates.superLayerHitRate[iSL] = 0;
169 } else {
171 }
173 for (int iLayer = 0 ; iLayer < f_nLayer ; ++iLayer)
174 if (m_nActiveWireInLayer[iLayer] == 0) {
175 m_rates.layerHitRate[iLayer] = 0;
176 } else {
177 m_rates.layerHitRate[iLayer] /= m_nActiveWireInLayer[iLayer];
178 }
180 for (int iLayer = 0 ; iLayer < f_nLayer ; ++iLayer)
181 for (int iPhi = 0 ; iPhi < f_nPhiDivision ; ++iPhi)
182 if (m_nActiveWireInLayerPhi[iLayer][iPhi] == 0) {
183 m_rates.layerPhiHitRate[iLayer][iPhi] = 0;
184 } else {
185 m_rates.layerPhiHitRate[iLayer][iPhi] /= m_nActiveWireInLayerPhi[iLayer][iPhi];
186 }
187
190 for (int i = 0 ; i < f_nSuperLayer ; ++i)
192 for (int i = 0 ; i < f_nLayer ; ++i)
194 for (int i = 0 ; i < f_nLayer ; ++i)
195 for (int j = 0 ; j < f_nPhiDivision ; ++j)
197 }
float superLayerHitRate[f_nSuperLayer]
SuperLayer average hit rate in kHz.
float layerPhiHitRate[f_nLayer][f_nPhiDivision]
Layer (in the phi bin) average hit rate in kHz.
int nActiveWireInLayerPhi[f_nLayer][f_nPhiDivision]
number of wires used in this analysis in each phi bin in each layer
float layerHitRate[f_nLayer]
Layer average hit rate in kHz.
int nActiveWireInTotal
number of wires used in this analysis in the whole CDC
int nActiveWireInSuperLayer[f_nSuperLayer]
number of wires used in this analysis in each super layer
int nActiveWireInLayer[f_nLayer]
number of wires used in this analysis in each layer
float averageRate
total detector average hit rate in KHz
void normalize()
normalize accumulated hits to hit rate in kHz
int timeWindowForSmallCell
time window for the small cells in 2*508.887 MHz clock ( 1 clock = 0.982536 ns)
int timeWindowForNormalCell
time window for the normal cells in 2*508.887 MHz clock ( 1 clock = 0.982536 ns)

Member Data Documentation

◆ f_nLayer

const int f_nLayer = 56
staticprivate

the number of layers

Definition at line 32 of file CDCHitRateCounter.h.

◆ f_nPhiDivision

const int f_nPhiDivision = 8
staticprivate

the number of division in phi

Definition at line 33 of file CDCHitRateCounter.h.

◆ f_nSuperLayer

const int f_nSuperLayer = 9
staticprivate

the number of super layers

Definition at line 31 of file CDCHitRateCounter.h.

◆ m_buffer

std::map<unsigned, TreeStruct> m_buffer
private

average rates in time stamps

Definition at line 159 of file CDCHitRateCounter.h.

◆ m_digits

StoreArray<CDCHit> m_digits
private

collection of digits

Definition at line 162 of file CDCHitRateCounter.h.

◆ m_enableBackgroundHitFilter

const bool m_enableBackgroundHitFilter
private

flag to enable the CDC background hit (crosstakl, noise) filter.

default: true

Definition at line 188 of file CDCHitRateCounter.h.

◆ m_enableBadWireTreatment

const bool m_enableBadWireTreatment
private

flag to enable the bad wire treatment.

default: true

Definition at line 187 of file CDCHitRateCounter.h.

◆ m_enableMarkBackgroundHit

const bool m_enableMarkBackgroundHit
private

flag to enable to mark background flag on CDCHit (set 0x100 bit for CDCHit::m_status).

default: false

Definition at line 190 of file CDCHitRateCounter.h.

◆ m_nActiveWireInLayer

int m_nActiveWireInLayer[f_nLayer] = {0}
private

the number of wires used in this hit-rate calculation in each layer

Definition at line 195 of file CDCHitRateCounter.h.

◆ m_nActiveWireInLayerPhi

int m_nActiveWireInLayerPhi[f_nLayer][f_nPhiDivision] = {{0}}
private

the number of wires used in this hit-rate calculation in each phi bin in each layer

Definition at line 196 of file CDCHitRateCounter.h.

◆ m_nActiveWireInSuperLayer

int m_nActiveWireInSuperLayer[f_nSuperLayer] = {0}
private

the number of wires used in this hit-rate calculation in each suler layer

Definition at line 194 of file CDCHitRateCounter.h.

◆ m_nActiveWireInTotal

int m_nActiveWireInTotal
private
Initial value:
=
0

the number of wires used in this hit-rate calculation in the whole CDC

Definition at line 192 of file CDCHitRateCounter.h.

◆ m_rates

TreeStruct m_rates
private

tree variables

Definition at line 156 of file CDCHitRateCounter.h.

◆ m_timeWindowLowerEdge_normalCell

const int m_timeWindowLowerEdge_normalCell
private

lower edge of the time window for normal cells [tdc count = 0.982536 ns]

Definition at line 167 of file CDCHitRateCounter.h.

◆ m_timeWindowLowerEdge_smallCell

const int m_timeWindowLowerEdge_smallCell
private

lower edge of the time window for small cells [tdc count = 0.982536 ns]

Definition at line 165 of file CDCHitRateCounter.h.

◆ m_timeWindowUpperEdge_normalCell

const int m_timeWindowUpperEdge_normalCell
private

upper edge of the time window for normal cells [tdc count = 0.982536 ns]

Definition at line 168 of file CDCHitRateCounter.h.

◆ m_timeWindowUpperEdge_smallCell

const int m_timeWindowUpperEdge_smallCell
private

upper edge of the time window for small cells [tdc count = 0.982536 ns]

Definition at line 166 of file CDCHitRateCounter.h.


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