Belle II Software development
ECLHitRateCounter.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9// Own header.
10#include <background/modules/BeamBkgHitRateMonitor/ECLHitRateCounter.h>
11
12#include <ecl/dbobjects/ECLCrystalCalib.h>
13#include <framework/database/DBObjPtr.h>
14#include <framework/geometry/B2Vector3.h>
15
16#include <numeric>
17
18using namespace std;
19
20
21namespace Belle2 {
26 namespace Background {
27
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 }
46
48 {
49 m_buffer.clear();
50 }
51
52 void ECLHitRateCounter::accumulate(unsigned timeStamp)
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 }
101
102 void ECLHitRateCounter::normalize(unsigned timeStamp)
103 {
104 // copy buffer element
105 m_rates = m_buffer[timeStamp];
106
107 if (not m_rates.validDspRate) return;
108
109 // normalize
111
112 //get average ret per crystal in segment
113 for (int i = 0; i < 16; i++) {
115 }
116 }
117
118
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 m_segmentMap.insert(std::pair<int, int>(cid, 2));
137 m_crystalsInSegment[2] += 1;
138 } else {
139 m_segmentMap.insert(std::pair<int, int>(cid, 3));
140 m_crystalsInSegment[3] += 1;
141 }
142 } else if (cid < 7777) {
143 if (z > 0) {
144 if (phi > 0.7853 && phi < 2.356) {
145 m_segmentMap.insert(std::pair<int, int>(cid, 4));
146 m_crystalsInSegment[4] += 1;
147 } else if (phi >= 2.356 || phi <= -2.356) {
148 m_segmentMap.insert(std::pair<int, int>(cid, 5));
149 m_crystalsInSegment[5] += 1;
150 } else if (phi > -2.356 && phi < -0.7853) {
151 m_segmentMap.insert(std::pair<int, int>(cid, 6));
152 m_crystalsInSegment[6] += 1;
153 } else {
154 m_segmentMap.insert(std::pair<int, int>(cid, 7));
155 m_crystalsInSegment[7] += 1;
156 }
157 } else {
158 if (phi > 0.7853 && phi < 2.356) {
159 m_segmentMap.insert(std::pair<int, int>(cid, 8));
160 m_crystalsInSegment[8] += 1;
161 } else if (phi >= 2.356 || phi <= -2.356) {
162 m_segmentMap.insert(std::pair<int, int>(cid, 9));
163 m_crystalsInSegment[9] += 1;
164 } else if (phi > -2.356 && phi < -0.7853) {
165 m_segmentMap.insert(std::pair<int, int>(cid, 10));
166 m_crystalsInSegment[10] += 1;
167 } else {
168 m_segmentMap.insert(std::pair<int, int>(cid, 11));
169 m_crystalsInSegment[11] += 1;
170 }
171 }
172 } else {
173 if (phi > 0.7853 && phi < 2.356) {
174 m_segmentMap.insert(std::pair<int, int>(cid, 12));
175 m_crystalsInSegment[12] += 1;
176 } else if (phi >= 2.356 || phi <= -2.356) {
177 m_segmentMap.insert(std::pair<int, int>(cid, 13));
178 m_crystalsInSegment[13] += 1;
179 } else if (phi > -2.356 && phi < -0.7853) {
180 m_segmentMap.insert(std::pair<int, int>(cid, 14));
181 m_crystalsInSegment[14] += 1;
182 } else {
183 m_segmentMap.insert(std::pair<int, int>(cid, 15));
184 m_crystalsInSegment[15] += 1;
185 }
186 }
187 }
188 }
189
190 } // Background namespace
192} // Belle2 namespace
DataType Phi() const
The azimuth angle.
Definition: B2Vector3.h:151
DataType Z() const
access variable Z (= .at(2) without boundary check)
Definition: B2Vector3.h:435
std::map< int, int > m_segmentMap
map with keys containing ECL CellID and values containing segment number
int m_crystalsInSegment[16]
array cotaining the number of crystals in given segment
std::map< unsigned, TreeStruct > m_buffer
average rates in time stamps
StoreArray< ECLDsp > m_dsps
collection of ECL waveforms
virtual void initialize(TTree *tree) override
Class initializer: set branch addresses and other staf.
Belle2::ECL::ECLGeometryPar * m_geometry
pointer to ECLGeometryPar
std::vector< float > m_waveformNoise
vector used to store ECL electronic noise constants foe each crystal
virtual void accumulate(unsigned timeStamp) override
Accumulate hits.
int findECLSegment(int cellid)
Find the correcsponding ECL segment based on the cellID.
void segmentECL()
Performs ECL segmentation; Done once per run; Populates a map which connects each ECL crystal with a ...
virtual void normalize(unsigned timeStamp) override
Normalize accumulated hits (e.g.
virtual void clear() override
Clear time-stamp buffer to prepare for 'accumulate'.
std::vector< float > m_ADCtoEnergy
vector used to store ECL calibration constants for each crystal
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
static ECLGeometryPar * Instance()
Static method to get a reference to the ECLGeometryPar instance.
ROOT::Math::XYZVector GetCrystalPos(int cid)
The Position of crystal.
void Mapping(int cid)
Mapping theta, phi Id.
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
Abstract base class for different kinds of events.
STL namespace.
bool validDspRate
status for rates calculated from waveforms, true if waveforms for all crystals are recorded
float averageDspBkgRate[16]
average background rate per crystal in given segment, calculated using ECL waveforms [hits/second]
void normalize()
normalize accumulated rates based on ECL waveforms