Belle II Software development
ECLHitRateCounter Class Reference

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

#include <ECLHitRateCounter.h>

Inheritance diagram for ECLHitRateCounter:
HitRateBase

Classes

struct  TreeStruct
 tree structure More...
 

Public Member Functions

 ECLHitRateCounter ()
 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 segmentECL ()
 Performs ECL segmentation; Done once per run; Populates a map which connects each ECL crystal with a segment number (0-15); Segments 0-3 are in the forward endcap, 4-7 are in the barrel with z<0, 8-11 are in the barrel with z>0, 12-15 are in the backward endcap; Segment 0 contains crystals with 45deg < phi < 135deg, Segment 1 contains crystals with 135deg < phi < 225deg, Segment 2 contains crystals with 225deg < phi < 315deg, Segment 3 contains crystals with phi < 45deg or phi > 315deg, With the same angular pattern continuing for barrel and BWD encap segments.
 
int findECLSegment (int cellid)
 Find the correcsponding ECL segment based on the cellID.
 
void findElectronicsNoise ()
 Find the electronics noise correction for each cellID Reads a file with a histogram containing electronics noise level of each crystal.
 

Private Attributes

TreeStruct m_rates
 tree variables
 
std::map< unsigned, TreeStructm_buffer
 average rates in time stamps
 
StoreArray< ECLDigitm_digits
 collection of digits
 
StoreArray< ECLDspm_dsps
 collection of ECL waveforms
 
std::vector< float > m_ADCtoEnergy
 vector used to store ECL calibration constants for each crystal
 
std::vector< float > m_waveformNoise
 vector used to store ECL electronic noise constants foe each crystal
 
Belle2::ECL::ECLGeometryParm_geometry {nullptr}
 pointer to ECLGeometryPar
 
std::map< int, int > m_segmentMap
 map with keys containing ECL CellID and values containing segment number
 
int m_crystalsInSegment [16] = {0}
 array containing the number of crystals in given segment
 

Detailed Description

Class for monitoring beam background hit rates of ECL.

Definition at line 28 of file ECLHitRateCounter.h.

Constructor & Destructor Documentation

◆ ECLHitRateCounter()

ECLHitRateCounter ( )
inline

Constructor.

Definition at line 59 of file ECLHitRateCounter.h.

60 {}

Member Function Documentation

◆ accumulate()

void accumulate ( unsigned timeStamp)
overridevirtual

Accumulate hits.

Parameters
timeStamptime stamp

Implements HitRateBase.

Definition at line 52 of file ECLHitRateCounter.cc.

53 {
54 // check if data are available
55 if (not m_dsps.isValid()) return;
56
57 //calculate rates using waveforms
58 //The background rate for a crystal is calculated as
59 //rate = rms_pedestal_squared / (average_photon_energy_squared * time_constant)
60 //where time_constant=2.53 us and average_photon_energy_squared = 1 MeV
61 if (m_dsps.getEntries() == 8736) {
62
63 // get buffer element
64 auto& rates = m_buffer[timeStamp];
65
66 // increment event counter
67 rates.numEvents++;
68
69 for (auto& aECLDsp : m_dsps) {
70
71 int nadc = aECLDsp.getNADCPoints();
72 int cellID = aECLDsp.getCellId();
73 int segmentNumber = findECLSegment(cellID);
74 int crysID = cellID - 1;
75 std::vector<int> dspAv = aECLDsp.getDspA();
76
77 //finding the pedestal value
78 double dspMean = (std::accumulate(dspAv.begin(), dspAv.begin() + nadc, 0.0)) / nadc;
79 double wpsum = 0;
80 for (int v = 0; v < nadc; v++) {
81 wpsum += pow(dspAv[v] - dspMean, 2);
82 }
83 double dspRMS = sqrt(wpsum / nadc);
84 double dspSigma = dspRMS * abs(m_ADCtoEnergy[crysID]);
85
86 //calculating the background rate per second
87 double dspBkgRate = ((pow(dspSigma, 2)) - (pow(m_waveformNoise[crysID], 2))) / (2.53 * 1e-12);
88 if (dspBkgRate < 0) {
89 dspBkgRate = 0;
90 }
91
92 //hit rate for segment in ECL, which is later normalized per 1Hz
93 rates.averageDspBkgRate[segmentNumber] += dspBkgRate;
94
95 }
96
97 // set flag to true to indicate the rates are valid
98 rates.validDspRate = true;
99 }
100 }
double sqrt(double a)
sqrt for double
Definition beamHelpers.h:28

◆ clear()

void clear ( )
overridevirtual

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

Implements HitRateBase.

Definition at line 47 of file ECLHitRateCounter.cc.

48 {
49 m_buffer.clear();
50 }

◆ findECLSegment()

int findECLSegment ( int cellid)
inlineprivate

Find the correcsponding ECL segment based on the cellID.

Parameters
cellidECL crystal CellID

Definition at line 108 of file ECLHitRateCounter.h.

109 {
110 return m_segmentMap.find(cellid)->second;
111 }

◆ initialize()

void initialize ( TTree * tree)
overridevirtual

Class initializer: set branch addresses and other staf.

Parameters
treea valid TTree pointer

Implements HitRateBase.

Definition at line 28 of file ECLHitRateCounter.cc.

29 {
30 segmentECL();
31 // register collection(s) as optional, your detector might be excluded in DAQ
32 m_dsps.isOptional();
33
34 // set branch address
35 tree->Branch("ecl", &m_rates, "averageDspBkgRate[16]/F:numEvents/I:validDspRate/O");
36
37 //ECL calibration
38 DBObjPtr<ECLCrystalCalib> ECLElectronicsCalib("ECLCrystalElectronics"), ECLECalib("ECLCrystalEnergy"),
39 ECLWaveformNoise("ECLWaveformRMS");
40 m_ADCtoEnergy.resize(8736);
41 if (ECLElectronicsCalib) for (int i = 0; i < 8736; i++) m_ADCtoEnergy[i] = ECLElectronicsCalib->getCalibVector()[i];
42 if (ECLECalib) for (int i = 0; i < 8736; i++) m_ADCtoEnergy[i] *= ECLECalib->getCalibVector()[i];
43 m_waveformNoise.resize(8736);
44 if (ECLWaveformNoise) for (int i = 0; i < 8736; i++) m_waveformNoise[i] = ECLWaveformNoise->getCalibVector()[i];
45 }

◆ normalize()

void normalize ( unsigned timeStamp)
overridevirtual

Normalize accumulated hits (e.g.

transform to rates)

Parameters
timeStamptime stamp

Implements HitRateBase.

Definition at line 102 of file ECLHitRateCounter.cc.

103 {
104 // copy buffer element
105 m_rates = m_buffer[timeStamp];
106
107 if (not m_rates.validDspRate) return;
108
109 // normalize
110 m_rates.normalize();
111
112 //get average ret per crystal in segment
113 for (int i = 0; i < 16; i++) {
114 m_rates.averageDspBkgRate[i] /= m_crystalsInSegment[i];
115 }
116 }

◆ segmentECL()

void segmentECL ( )
private

Performs ECL segmentation; Done once per run; Populates a map which connects each ECL crystal with a segment number (0-15); Segments 0-3 are in the forward endcap, 4-7 are in the barrel with z<0, 8-11 are in the barrel with z>0, 12-15 are in the backward endcap; Segment 0 contains crystals with 45deg < phi < 135deg, Segment 1 contains crystals with 135deg < phi < 225deg, Segment 2 contains crystals with 225deg < phi < 315deg, Segment 3 contains crystals with phi < 45deg or phi > 315deg, With the same angular pattern continuing for barrel and BWD encap segments.

Definition at line 119 of file ECLHitRateCounter.cc.

120 {
122 for (int cid = 1; cid < 8737; cid++) {
123 m_geometry->Mapping(cid);
124 const B2Vector3D position = m_geometry->GetCrystalPos(cid - 1);
125 const double phi = position.Phi();
126 const double z = position.Z();
127
128 if (cid < 1297) {
129 if (phi > 0.7853 && phi < 2.356) {
130 m_segmentMap.insert(std::pair<int, int>(cid, 0));
131 m_crystalsInSegment[0] += 1;
132 } else if (phi >= 2.356 || phi <= -2.356) {
133 m_segmentMap.insert(std::pair<int, int>(cid, 1));
134 m_crystalsInSegment[1] += 1;
135 // } else if (phi > -2.356 && phi < -0.7853) {
136 } else if (phi < -0.7853) {
137 m_segmentMap.insert(std::pair<int, int>(cid, 2));
138 m_crystalsInSegment[2] += 1;
139 } else {
140 m_segmentMap.insert(std::pair<int, int>(cid, 3));
141 m_crystalsInSegment[3] += 1;
142 }
143 } else if (cid < 7777) {
144 if (z > 0) {
145 if (phi > 0.7853 && phi < 2.356) {
146 m_segmentMap.insert(std::pair<int, int>(cid, 4));
147 m_crystalsInSegment[4] += 1;
148 } else if (phi >= 2.356 || phi <= -2.356) {
149 m_segmentMap.insert(std::pair<int, int>(cid, 5));
150 m_crystalsInSegment[5] += 1;
151 // } else if (phi > -2.356 && phi < -0.7853) {
152 } else if (phi < -0.7853) {
153 m_segmentMap.insert(std::pair<int, int>(cid, 6));
154 m_crystalsInSegment[6] += 1;
155 } else {
156 m_segmentMap.insert(std::pair<int, int>(cid, 7));
157 m_crystalsInSegment[7] += 1;
158 }
159 } else {
160 if (phi > 0.7853 && phi < 2.356) {
161 m_segmentMap.insert(std::pair<int, int>(cid, 8));
162 m_crystalsInSegment[8] += 1;
163 } else if (phi >= 2.356 || phi <= -2.356) {
164 m_segmentMap.insert(std::pair<int, int>(cid, 9));
165 m_crystalsInSegment[9] += 1;
166 // } else if (phi > -2.356 && phi < -0.7853) {
167 } else if (phi < -0.7853) {
168 m_segmentMap.insert(std::pair<int, int>(cid, 10));
169 m_crystalsInSegment[10] += 1;
170 } else {
171 m_segmentMap.insert(std::pair<int, int>(cid, 11));
172 m_crystalsInSegment[11] += 1;
173 }
174 }
175 } else {
176 if (phi > 0.7853 && phi < 2.356) {
177 m_segmentMap.insert(std::pair<int, int>(cid, 12));
178 m_crystalsInSegment[12] += 1;
179 } else if (phi >= 2.356 || phi <= -2.356) {
180 m_segmentMap.insert(std::pair<int, int>(cid, 13));
181 m_crystalsInSegment[13] += 1;
182 // } else if (phi > -2.356 && phi < -0.7853) {
183 } else if (phi < -0.7853) {
184 m_segmentMap.insert(std::pair<int, int>(cid, 14));
185 m_crystalsInSegment[14] += 1;
186 } else {
187 m_segmentMap.insert(std::pair<int, int>(cid, 15));
188 m_crystalsInSegment[15] += 1;
189 }
190 }
191 }
192 }
static ECLGeometryPar * Instance()
Static method to get a reference to the ECLGeometryPar instance.
B2Vector3< double > B2Vector3D
typedef for common usage with double
Definition B2Vector3.h:516

Member Data Documentation

◆ m_ADCtoEnergy

std::vector<float> m_ADCtoEnergy
private

vector used to store ECL calibration constants for each crystal

Definition at line 129 of file ECLHitRateCounter.h.

◆ m_buffer

std::map<unsigned, TreeStruct> m_buffer
private

average rates in time stamps

Definition at line 123 of file ECLHitRateCounter.h.

◆ m_crystalsInSegment

int m_crystalsInSegment[16] = {0}
private

array containing the number of crystals in given segment

Definition at line 135 of file ECLHitRateCounter.h.

135{0};

◆ m_digits

StoreArray<ECLDigit> m_digits
private

collection of digits

Definition at line 126 of file ECLHitRateCounter.h.

◆ m_dsps

StoreArray<ECLDsp> m_dsps
private

collection of ECL waveforms

Definition at line 127 of file ECLHitRateCounter.h.

◆ m_geometry

Belle2::ECL::ECLGeometryPar* m_geometry {nullptr}
private

pointer to ECLGeometryPar

Definition at line 133 of file ECLHitRateCounter.h.

133{nullptr};

◆ m_rates

TreeStruct m_rates
private

tree variables

Definition at line 120 of file ECLHitRateCounter.h.

◆ m_segmentMap

std::map<int, int> m_segmentMap
private

map with keys containing ECL CellID and values containing segment number

Definition at line 134 of file ECLHitRateCounter.h.

◆ m_waveformNoise

std::vector<float> m_waveformNoise
private

vector used to store ECL electronic noise constants foe each crystal

Definition at line 130 of file ECLHitRateCounter.h.


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