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 seting 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 71 of file CDCHitRateCounter.cc.

72 {
73 // check if data are available
74 if (not m_digits.isValid()) return;
75
76 // get buffer element
77 auto& rates = m_buffer[timeStamp];
78
79 // increment event counter
80 rates.numEvents++;
81
83 std::map<const CDCHit*, bool> CDCHitToBackgroundFlag;
85 StoreWrappedObjPtr<std::vector<CDCWireHit>> storeVector("CDCWireHitVector");
86 if (not storeVector) {
87 B2FATAL("CDCWireHitVector is unaccessible in DataStore."
88 "Need TFCDC_WireHitParameter module before.");
89 }
90 const std::vector<CDCWireHit>& cdcWireHitVector = *storeVector;
91 for (const auto& cdcWireHit : cdcWireHitVector) {
92 const CDCHit* cdcHit = cdcWireHit.getHit();
93 CDCHitToBackgroundFlag[cdcHit] = cdcWireHit->hasBackgroundFlag();
94 }
95 }
96
97 CDC::CDCGeometryPar& geometryPar = CDC::CDCGeometryPar::Instance();
98
106 for (CDCHit& hit : m_digits) {
107 const WireID wireID(hit.getID());
108 if (m_enableBadWireTreatment && geometryPar.isBadWire(wireID))
109 continue;
110
111 const unsigned short iSuperLayer = hit.getISuperLayer();
112 const unsigned short iLayer = hit.getICLayer();
113 const unsigned short iWireInLayer = hit.getIWire();
114 const unsigned short iPhiBin = getIPhiBin(iSuperLayer, iWireInLayer);
115 const unsigned short adc = hit.getADCCount();
116 const short tdc = hit.getTDCCount();
117
119 if (CDCHitToBackgroundFlag[&hit]) {
121 unsigned short newStatus = (hit.getStatus() | 0x100);
122 hit.setStatus(newStatus);
123 }
124 continue;
125 }
126 } else {
127 if (iSuperLayer == 0 && adc < 15)
128 continue;
129 if (iSuperLayer != 0 && adc < 18)
130 continue;
131 }
132
133 if (not isInTimeWindow(iSuperLayer, tdc))
134 continue;
135
136
137 rates.averageRate += 1;
138 rates.superLayerHitRate[iSuperLayer] += 1;
139 rates.layerHitRate[iLayer] += 1;
140 rates.layerPhiHitRate[iLayer][iPhiBin] += 1;
141 }
142
143 // set flag to true to indicate the rates are valid
144 rates.valid = true;
145 }
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 66 of file CDCHitRateCounter.cc.

67 {
68 m_buffer.clear();
69 }

◆ 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 230 of file CDCHitRateCounter.cc.

231 {
232 static const unsigned short nlayer_in_SL[f_nSuperLayer] = { 8, 6, 6,
233 6, 6, 6,
234 6, 6, 6
235 };
236 static const unsigned short nwire_in_layer[f_nSuperLayer] = { 160, 160, 192,
237 224, 256, 288,
238 320, 352, 384
239 };
240
241 CDC::CDCGeometryPar& geometryPar = CDC::CDCGeometryPar::Instance();
242
244 int contLayerID = 0;
245
246 for (int iSL = 0 ; iSL < f_nSuperLayer ; ++iSL) {
248 for (int iL = 0 ; iL < nlayer_in_SL[iSL] ; ++iL) {
249 m_nActiveWireInLayer[contLayerID] = 0;
250 for (int iPhi = 0 ; iPhi < f_nPhiDivision ; ++iPhi)
251 m_nActiveWireInLayerPhi[contLayerID][iPhi] = 0;
252
253 for (int i = 0 ; i < nwire_in_layer[iSL] ; ++i) {
254 WireID wireID(iSL, iL, i);
255 if (not geometryPar.isBadWire(wireID)) {
256 ++m_nActiveWireInLayerPhi[contLayerID][getIPhiBin(iSL, i)];
257 ++m_nActiveWireInLayer[contLayerID];
258 }
259 }
261 ++contLayerID;
262 }
264 }
265
266
267 //B2INFO("CDC, # of Active wires / # of total wires");
268 std::cout << "CDC, # of Active wires / # of total wires" << std::endl;
269 int contLayerID_2 = 0;
270 for (int iSL = 0 ; iSL < f_nSuperLayer ; ++iSL) {
271 for (int iL = 0 ; iL < nlayer_in_SL[iSL] ; ++iL) {
272 //B2INFO("Layer "<< contLayerID_2 << ": "
273 std::cout << "Layer " << contLayerID_2 << ": "
274 << m_nActiveWireInLayer[contLayerID_2] << " / "
275 //<< nwire_in_layer[iSL] );
276 << nwire_in_layer[iSL] << std::endl;
277 ++contLayerID_2;
278 }
279 }
280
281 }
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 202 of file CDCHitRateCounter.cc.

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

◆ getIPhiBin()

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

get the bin ID of the division.

Definition at line 284 of file CDCHitRateCounter.cc.

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

◆ initialize()

void initialize ( TTree *  tree)
overridevirtual

Class initializer: set branch addresses and other staf.

Parameters
treea valid TTree pointer

Implements HitRateBase.

Definition at line 36 of file CDCHitRateCounter.cc.

37 {
38 // register collection(s) as optional, your detector might be excluded in DAQ
39 m_digits.isOptional();
40
41 // set branch address
42 //tree->Branch("cdc", &m_rates, "averageRate/F:numEvents/I:valid/O");
43 stringstream leafList;
44 leafList
45 << "averageRate/F:"
46 << "superLayerHitRate[" << f_nSuperLayer << "]/F:"
47 << "layerHitRate[" << f_nLayer << "]/F:"
48 << "layerPhiHitRate[" << f_nLayer << "][" << f_nPhiDivision << "]/F:"
49 << "timeWindowForSmallCell/I:"
50 << "timeWindowForNormalCell/I:"
51 << "nActiveWireInTotal/I:"
52 << "nActiveWireInSuperLayer[" << f_nSuperLayer << "]/I:"
53 << "nActiveWireInLayer[" << f_nLayer << "]/I:"
54 << "nActiveWireInLayerPhi[" << f_nLayer << "][" << f_nPhiDivision << "]/I:"
55 << "numEvents/I:"
56 << "valid/O";
57 tree->Branch("cdc", &m_rates, leafList.str().c_str());
58
61 } else {
63 }
64 }
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 147 of file CDCHitRateCounter.cc.

148 {
149 // copy buffer element
150 m_rates = m_buffer[timeStamp];
151
152 if (not m_rates.valid) return;
153
156
157 // normalize : hit rate in the time stamp in kHz
159
162 if (m_nActiveWireInTotal == 0) {
164 } else {
166 }
168 for (int iSL = 0 ; iSL < f_nSuperLayer ; ++iSL)
169 if (m_nActiveWireInSuperLayer[iSL] == 0) {
170 m_rates.superLayerHitRate[iSL] = 0;
171 } else {
173 }
175 for (int iLayer = 0 ; iLayer < f_nLayer ; ++iLayer)
176 if (m_nActiveWireInLayer[iLayer] == 0) {
177 m_rates.layerHitRate[iLayer] = 0;
178 } else {
179 m_rates.layerHitRate[iLayer] /= m_nActiveWireInLayer[iLayer];
180 }
182 for (int iLayer = 0 ; iLayer < f_nLayer ; ++iLayer)
183 for (int iPhi = 0 ; iPhi < f_nPhiDivision ; ++iPhi)
184 if (m_nActiveWireInLayerPhi[iLayer][iPhi] == 0) {
185 m_rates.layerPhiHitRate[iLayer][iPhi] = 0;
186 } else {
187 m_rates.layerPhiHitRate[iLayer][iPhi] /= m_nActiveWireInLayerPhi[iLayer][iPhi];
188 }
189
192 for (int i = 0 ; i < f_nSuperLayer ; ++i)
194 for (int i = 0 ; i < f_nLayer ; ++i)
196 for (int i = 0 ; i < f_nLayer ; ++i)
197 for (int j = 0 ; j < f_nPhiDivision ; ++j)
199 }
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: