Belle II Software development
eclLOMModule.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 <ecl/modules/eclLOM/eclLOMModule.h>
11
12/* ROOT headers. */
13#include <Math/Boost.h>
14
15using namespace std;
16using namespace Belle2;
17using namespace ECL;
18
19REG_MODULE(ECLLOM);
20
21
23{
24 setDescription("module to emulate Luminosity Online Monitor");
25 addParam("thresholdFE", m_thresholdFE, "Threshold for Forward Endcap [GeV]", 3.0);
26 addParam("thresholdBE", m_thresholdBE, "Threshold for Backward Endcap [GeV]", 1.0);
27 addParam("thresholdBkg", m_thresholdBkg, "Threshold when sector considered as lighted [GeV]", 0.5);
28 addParam("discrTime", m_discrTime, "Duration of '1' (positive) signal from discriminators [ns]", 1000.0);
29 addParam("includeInnerFE", m_includeInnerFE, "Flag to include/exclude Inner part of the Forward Endcap", false);
30 addParam("saveSignal", m_saveSignal, "Flag to store or not signals' waveforms", false);
31 addParam("testFileName", m_testFileName, "output file", std::string("lomtest.root"));
32
33 m_evtNum = 0;
34 for (int i = 0; i < 16; i++) {
35 for (int j = 0; j < 16; j++) {
38 }
39 }
40}
41
43{
44
45}
46
48{
49 if (!(m_MCParticles.isRequired() && m_TrgEclWaveforms.isRequired())) {
50 //Fatal is not necessary here as the storeArrays should just look
51 //empty if not registered but let's make sure everything is present
52 B2FATAL("Not all collections found, exiting processing");
53 }
54
55 m_testfile = new TFile(m_testFileName.c_str(), "RECREATE");
56 m_testtree = new TTree("lom_tree", "");
57
58 m_testtree->Branch("ev", &m_evtNum);
59 m_testtree->Branch("BE_amp[16]", m_BE_Amplitude, "BE_amp[16]/D");
60 m_testtree->Branch("FE_amp[16]", m_FE_Amplitude, "FE_amp[16]/D");
61 m_testtree->Branch("FESum_MaxAmp", &m_FESum_MaxAmp);
62 m_testtree->Branch("BESum_MaxAmp", &m_BESum_MaxAmp);
63 m_testtree->Branch("FESum_MaxId", &m_FESum_MaxId);
64 m_testtree->Branch("BESum_MaxId", &m_BESum_MaxId);
65 m_testtree->Branch("BE_Pedal[16]", m_BE_Pedal, "BE_Pedal[16]/D");
66 m_testtree->Branch("FE_Pedal[16]", m_FE_Pedal, "FE_Pedal[16]/D");
67
68 m_testtree->Branch("Bhabha", &m_isBhabha);
69 m_testtree->Branch("BhNum", &m_BhNum);
70
71 m_testtree->Branch("mc_en[2]", m_mcen, "mc_en[2]/D");
72 m_testtree->Branch("mc_th[2]", m_mcth, "mc_th[2]/D");
73 m_testtree->Branch("mc_ph[2]", m_mcph, "mc_ph[2]/D");
74
75 m_testtree->Branch("com_en[2]", m_com_en, "com_en[2]/D");
76 m_testtree->Branch("com_th[2]", m_com_th, "com_th[2]/D");
77 m_testtree->Branch("com_ph[2]", m_com_ph, "com_ph[2]/D");
78
79 if (m_saveSignal) {
80 m_testtree->Branch("BE_wf[16][64]", m_BE_Waveform_100ns, "BE_wf[16][64]/D");
81 m_testtree->Branch("FE_wf[16][64]", m_FE_Waveform_100ns, "FE_wf[16][64]/D");
82 }
83
84 //additional histograms. Represent data over the whole dataset:
85 m_h2Coin = new TH2D("Coins", "Coincidence Matrix", 16, 0, 16, 16, 0, 16);
86 m_h2SumCoin = new TH2D("SumCoins", "Sum Coincidence Matrix", 16, 0, 16, 16, 0, 16);
87 m_h2FEAmp = new TH2D("FE_AmpId", "", 16, 0, 16, 100, 0, 8);
88 m_h2BEAmp = new TH2D("BE_AmpId", "", 16, 0, 16, 100, 0, 8);
89 m_h1BEHits = new TH1D("BE_Fired", "", 16, 0, 16);
90 m_h1FEHits = new TH1D("FE_Fired", "", 16, 0, 16);
91
92 m_NSamples = 631;
93}
94
96{
97}
98
100{
106 // LOM logic
107 for (int iSample = 300; iSample < 500; iSample++) { //300-500 window where BhaBha is expected
110 calculate_coincidence(iSample);
111 // generate bhabha signal
112 for (int iFESector = 0; iFESector < 16; iFESector++) {
113 // check opposite running sum:
114 int iBESector = (iFESector + 8) % 16;
115 if (m_SumCoincidenceMatrix[iFESector][iBESector] == 1 && m_isBhabhaPatternFE && m_isBhabhaPatternBE) {
116 //coincidence at first tick
117 m_isBhabha = true;
118 m_BhNum++;
119 }
120 }
121 }
122 m_testtree->Fill();
123 m_evtNum++;
124}
125
127{
128}
129
131{
132
133 for (int i = 0; i < 16; i++) {
134 for (int j = 0; j < 16; j++) {
135 std::cout << m_CoincidenceCounterMatrix[i][j] << " ";
136 m_h2Coin->SetBinContent(i + 1, j + 1, m_CoincidenceCounterMatrix[i][j]);
137 m_h2SumCoin->SetBinContent(i + 1, j + 1, m_SumCoincidenceCounterMatrix[i][j]);
138 }
139 std::cout << std::endl;
140 }
141 m_testfile->Write();
142 m_testfile->Close();
143}
144
145
147{
148 int nm_MCParticles = m_MCParticles.getEntries();
149 if (nm_MCParticles >= 4) {
150 for (int ind = 2; ind < 4; ind++) {
151 m_mcen[ind - 2] = m_MCParticles[ind]->getEnergy();
152 m_mcth[ind - 2] = m_MCParticles[ind]->getMomentum().Theta();
153 m_mcph[ind - 2] = m_MCParticles[ind]->getMomentum().Phi();
154 }
155
156 ROOT::Math::PxPyPzEVector SummP(m_MCParticles[0]->get4Vector() + m_MCParticles[1]->get4Vector());
157 ROOT::Math::XYZVector Boost_backV = SummP.BoostToCM();
158 ROOT::Math::PxPyPzEVector ComP[2];
159 ComP[0] = m_MCParticles[2]->get4Vector();
160 ComP[1] = m_MCParticles[3]->get4Vector();
161 ComP[0] = ROOT::Math::Boost(Boost_backV) * ComP[0];
162 ComP[1] = ROOT::Math::Boost(Boost_backV) * ComP[1];
163 for (int ind = 0; ind < 2; ind++) {
164 m_com_en[ind] = ComP[ind].E();
165 m_com_th[ind] = ComP[ind].Theta();
166 m_com_ph[ind] = ComP[ind].Phi();
167 }
168 }
169}
170
172{
173 //int n_trg_digi = TrgEclDigiArray.getEntries();
174 int n_trg_wf = m_TrgEclWaveforms.getEntries();
175 // calculate signals of endcap sectors for LOM input
176 // as sum of corresponding TC signals
177 for (int i = 0; i < n_trg_wf; i++) {
178 const TRGECLWaveform* TCWaveform = m_TrgEclWaveforms[i];
179 //int m_tcid = TCWaveform->getTCID();
180 int tc_theta_id = TCWaveform->getThetaID(); //FE:1,2,3 BE:16,17 Checked for rel 4 02 08
181 int tc_phi_id = TCWaveform->getPhiID(); // 1 - 32
182 double m_wf[64];
183 TCWaveform->fillWaveform(m_wf);
184
185 int iSectorIndex = (tc_phi_id - 1) / 2; // from 0 to 15
186 if (tc_theta_id == 1 && !m_includeInnerFE) continue;
187
188 for (int iSample = 0; iSample < 64; iSample++) {
189 if (tc_theta_id <= 3) { //Forward Endcap
190 m_FE_Waveform_100ns[iSectorIndex][iSample] += m_wf[iSample];
191 } else { // Backward Endcap
192 if (tc_theta_id == 16 || tc_theta_id == 17) m_BE_Waveform_100ns[iSectorIndex][iSample] += m_wf[iSample];
193 }
194 }
195 }
196}
197
199{
200 for (int iSector = 0; iSector < 16; iSector++) { // Calculating pedestals
201 for (int iSample = 15; iSample < 36; iSample++) {
202 m_BE_Pedal[iSector] += m_BE_Waveform_100ns[iSector][iSample] / 20;
203 m_FE_Pedal[iSector] += m_FE_Waveform_100ns[iSector][iSample] / 20;
204 }
205 }
206 double dAdT; // convert 100 ns signal to 10 ns
207 for (int iSector = 0; iSector < 16; iSector++) { // Linear interpolation from 100ns to 10ns
208 for (int iSample = 0; iSample < 63; iSample++) {
209 // forward
210 dAdT = (m_FE_Waveform_100ns[iSector][iSample + 1] - m_FE_Waveform_100ns[iSector][iSample]) / 10.0;
211 m_FE_Waveform_100ns[iSector][iSample] -= m_FE_Pedal[iSector]; //remove pedestals
212 for (int j = 0; j < 10; j++) m_FE_Waveform_10ns[iSector][iSample * 10 + j] = m_FE_Waveform_100ns[iSector][iSample] + j * dAdT;
213 m_FE_Waveform_10ns[iSector][630] = m_FE_Waveform_100ns[iSector][63] - m_FE_Pedal[iSector];
214 //backward
215 dAdT = (m_BE_Waveform_100ns[iSector][iSample + 1] - m_BE_Waveform_100ns[iSector][iSample]) / 10.0;
216 m_BE_Waveform_100ns[iSector][iSample] -= m_BE_Pedal[iSector];
217 for (int j = 0; j < 10; j++) m_BE_Waveform_10ns[iSector][iSample * 10 + j] = m_BE_Waveform_100ns[iSector][iSample] + j * dAdT;
218 m_BE_Waveform_10ns[iSector][630] = m_BE_Waveform_100ns[iSector][63] - m_BE_Pedal[iSector];
219 }
220 }
221 // calculate running sums for 10ns signal
222 int TimeOfDiscr = int(m_discrTime / 10); //discriminator duration in samples
223 for (int iSector = 0; iSector < 16; iSector++) {
224 for (int iSample = 1; iSample < m_NSamples; iSample++) {
225 int iNextSector = (iSector + 1) % 16;
226 m_BESum_Waveform_10ns[iSector][iSample] = m_BE_Waveform_10ns[iSector][iSample] + m_BE_Waveform_10ns[iNextSector][iSample];
227 m_FESum_Waveform_10ns[iSector][iSample] = m_FE_Waveform_10ns[iSector][iSample] + m_FE_Waveform_10ns[iNextSector][iSample];
228
229 //filling Discriminators' signals
230 if (m_FESum_Waveform_10ns[iSector][iSample] > m_thresholdFE && m_FESum_Discr[iSector][iSample] == 0) {
231 for (int j = iSample; j < iSample + TimeOfDiscr; j++) {
232 if (j < m_NSamples) m_FESum_Discr[iSector][j] = 1;
233 }
234 }
235 if (m_BESum_Waveform_10ns[iSector][iSample] > m_thresholdBE && m_BESum_Discr[iSector][iSample] == 0) {
236 for (int j = iSample; j < iSample + TimeOfDiscr; j++) {
237 if (j < m_NSamples) m_BESum_Discr[iSector][j] = 1;
238 }
239 }
240 if (m_FE_Waveform_10ns[iSector][iSample] > m_thresholdBkg && m_FEQual_Discr[iSector][iSample] == 0) {
241 for (int j = iSample; j < iSample + TimeOfDiscr; j++) {
242 if (j < m_NSamples) m_FEQual_Discr[iSector][j] = 1;
243 }
244 }
245 if (m_BE_Waveform_10ns[iSector][iSample] > m_thresholdBkg && m_BEQual_Discr[iSector][iSample] == 0) {
246 for (int j = iSample; j < iSample + TimeOfDiscr; j++) {
247 if (j < m_NSamples) m_BEQual_Discr[iSector][j] = 1;
248 }
249 }
250 }
251 }
252}
253
254
256{
257 int nhit = 0;
258 int First = 0;
259 // calculate quality signal for backward endcap
260 for (int iBESector = 0; iBESector < 16; iBESector++) {
261 if (m_BEQual_Discr[iBESector][iSample]) {
262 nhit++;
263 if (nhit == 1) First = iBESector;
264 if (nhit == 2 && !((iBESector + 1) % 16 == First || (First + 1) % 16 == iBESector)) return (false);
265 if (nhit >= 3) return (false);
266 }
267 }
268 return (true);
269}
270
272{
273 int nhit = 0;
274 int First = 0;
275 for (int iFESector = 0; iFESector < 16; iFESector++) {
276 if (m_FEQual_Discr[iFESector][iSample]) {
277 nhit++;
278 if (nhit == 1) First = iFESector;
279 if (nhit == 2 && !((iFESector + 1) % 16 == First || (First + 1) % 16 == iFESector)) return (false);
280 if (nhit >= 3) return (false);
281 }
282 }
283 return (true);
284}
285
287{
288 for (int iFESector = 0; iFESector < 16; iFESector++) {
289 for (int iBESector = 0; iBESector < 16; iBESector++) {
290
291 if (m_FE_Waveform_10ns[iFESector][iSample] > m_thresholdFE && m_BE_Waveform_10ns[iBESector][iSample] > m_thresholdBE) {
292 if (m_CoincidenceMatrix[iFESector][iBESector] == 0) m_CoincidenceCounterMatrix[iFESector][iBESector]++;
293 m_CoincidenceMatrix[iFESector][iBESector]++;
294 } else {
295 m_CoincidenceMatrix[iFESector][iBESector] = 0;
296 }
297
298 if (m_FESum_Discr[iFESector][iSample] && m_BESum_Discr[iBESector][iSample]) {
299 if (m_SumCoincidenceMatrix[iFESector][iBESector] == 0) m_SumCoincidenceCounterMatrix[iFESector][iBESector]++;
300 m_SumCoincidenceMatrix[iFESector][iBESector]++;
301 }
302 m_SumCoincidenceMatrix[iFESector][iBESector] = 0;
303 }
304 }
305}
306
307
309{
310 for (int isector = 0; isector < 16; isector++) {
311 for (int iSample = 0; iSample < 64; iSample++) {
312 m_BE_Waveform_100ns[isector][iSample] = 0;
313 m_FE_Waveform_100ns[isector][iSample] = 0;
314 }
315 for (int iSample = 0; iSample < m_NSamples; iSample++) {
316 m_BE_Waveform_10ns[isector][iSample] = 0;
317 m_FE_Waveform_10ns[isector][iSample] = 0;
318 m_BESum_Waveform_10ns[isector][iSample] = 0;
319 m_FESum_Waveform_10ns[isector][iSample] = 0;
320 m_FESum_Discr[isector][iSample] = 0;
321 m_BESum_Discr[isector][iSample] = 0;
322 m_FEQual_Discr[isector][iSample] = 0;
323 m_BEQual_Discr[isector][iSample] = 0;
324 }
325 m_BE_Pedal[isector] = 0;
326 m_FE_Pedal[isector] = 0;
327 m_BE_Amplitude[isector] = 0;
328 m_FE_Amplitude[isector] = 0;
329 m_BESum_Amplitude[isector] = 0;
330 m_FESum_Amplitude[isector] = 0;
331 for (int jsector = 0; jsector < 16; jsector++) {
332 m_CoincidenceMatrix[isector][jsector] = 0;
333 m_SumCoincidenceMatrix[isector][jsector] = 0;
334 }
335 }
336 m_isBhabha = 0;
337 m_BhNum = 0;
338 m_FESum_MaxAmp = 0;
339 m_BESum_MaxAmp = 0;
340 m_FESum_MaxId = -1;
341 m_BESum_MaxId = -1;
342}
343
345{
346 for (int iSample = 0; iSample < m_NSamples; iSample++) {
347 for (int isector = 0; isector < 16; isector++) {
348 if (m_FE_Waveform_10ns[isector][iSample] > m_FE_Amplitude[isector]) m_FE_Amplitude[isector] = m_FE_Waveform_10ns[isector][iSample];
349 if (m_FESum_Waveform_10ns[isector][iSample] > m_FESum_Amplitude[isector]) m_FESum_Amplitude[isector] =
350 m_FESum_Waveform_10ns[isector][iSample];
351 if (m_BE_Waveform_10ns[isector][iSample] > m_BE_Amplitude[isector]) m_BE_Amplitude[isector] = m_BE_Waveform_10ns[isector][iSample];
352 if (m_BESum_Waveform_10ns[isector][iSample] > m_BESum_Amplitude[isector]) m_BESum_Amplitude[isector] =
353 m_BESum_Waveform_10ns[isector][iSample];
354 }
355 }
356 for (int i = 0; i < 16; i++) {
359 m_FESum_MaxId = i;
360 }
363 m_BESum_MaxId = i;
364 }
365 if (m_FE_Amplitude[i] > 0.5) m_h2FEAmp->Fill(i, m_FE_Amplitude[i]);
366 if (m_BE_Amplitude[i] > 0.5) m_h2BEAmp->Fill(i, m_BE_Amplitude[i]);
367
368 if (m_FE_Amplitude[i] > m_thresholdFE) m_h1FEHits->Fill(i);
369 if (m_BE_Amplitude[i] > m_thresholdBE) m_h1BEHits->Fill(i);
370 }
371}
bool m_isBhabhaPatternBE
Quality signal for Backward endcap.
Definition: eclLOMModule.h:160
double m_FESum_MaxAmp
Maximum running sum amplitude in an event for Forward endcap.
Definition: eclLOMModule.h:134
int m_BESum_MaxId
Id of a sector with maximum amplitude in Backward endcap.
Definition: eclLOMModule.h:137
int m_BhNum
Number of Bha-bha signals in an event.
Definition: eclLOMModule.h:127
TH2D * m_h2BEAmp
Store sectors amplitudes for Backward endcap over all events.
Definition: eclLOMModule.h:143
TH2D * m_h2FEAmp
Store sectors amplitudes for Forward endcap over all events.
Definition: eclLOMModule.h:142
TH2D * m_h2SumCoin
Store number of coincedencies in running sums for i:j sectors (Forward:Backward) over all events.
Definition: eclLOMModule.h:141
double m_FE_Waveform_100ns[16][64]
Waveforms with 100ns sampling for Forward Endcap sectors.
Definition: eclLOMModule.h:110
double m_mcph[2]
Monte Carlo phi of the final state particles in main frame.
Definition: eclLOMModule.h:113
int m_NSamples
m_NSamples=631, number of samples for 10ns sampling.
Definition: eclLOMModule.h:150
int m_evtNum
Event number.
Definition: eclLOMModule.h:108
TFile * m_testfile
File to save output.
Definition: eclLOMModule.h:148
void calculate_discr_output()
Transforms waveforms into discriminators output.
double m_BESum_MaxAmp
Maximum running sum amplitude in an event for Backward endcap.
Definition: eclLOMModule.h:135
bool calculate_BE_quality(int iSample)
Return Quality (topology) flag at sample point, iSample, for Backward Endcap.
int m_CoincidenceCounterMatrix[16][16]
Stores number of concidences between waveforms exceeding threshold in i:j sectors (Forward:Backward).
Definition: eclLOMModule.h:163
int m_FESum_MaxId
Id of a sector with maximum amplitude in Forward endcap.
Definition: eclLOMModule.h:136
virtual void initialize() override
Initialize variables.
Definition: eclLOMModule.cc:47
double m_thresholdBE
Threshold [GeV] on signal for Backward Endcap .
Definition: eclLOMModule.h:101
TH1D * m_h1BEHits
Store number of events when Backward sector i has signal exceeding Bha-Bha threshold over all events.
Definition: eclLOMModule.h:144
virtual void event() override
event per event.
Definition: eclLOMModule.cc:99
void get_MCparticles()
Get MC particles parameters.
double m_thresholdFE
Threshold [GeV] on signal for Forward Endcap .
Definition: eclLOMModule.h:100
double m_BE_Waveform_10ns[16][631]
Waveforms with 10ns sampling for Backward Endcap sectors.
Definition: eclLOMModule.h:151
StoreArray< TRGECLWaveform > m_TrgEclWaveforms
Trigger waveforms.
Definition: eclLOMModule.h:122
bool m_FESum_Discr[16][631]
Discriminators values for running sums of Forward Endcap.
Definition: eclLOMModule.h:155
int m_SumCoincidenceMatrix[16][16]
Stores current coincidence duration [in samples] between running sums discriminators in i:j sectors (...
Definition: eclLOMModule.h:162
virtual void endRun() override
end run.
virtual ~ECLLOMModule()
Destructor.
Definition: eclLOMModule.cc:42
virtual void terminate() override
terminate.
void get_waveforms()
Get ECL waveforms comdined into sectors.
double m_discrTime
Discriminator's signal duration in ns.
Definition: eclLOMModule.h:103
double m_FESum_Amplitude[16]
Calculated amplitudes in running sums of Forward Endcap.
Definition: eclLOMModule.h:131
bool m_FEQual_Discr[16][631]
Discriminators values for Quality signal of Forward Endcap.
Definition: eclLOMModule.h:157
int m_CoincidenceMatrix[16][16]
Stores current coincidence duration [in samples] between waveforms exceeding threshold in i:j sectors...
Definition: eclLOMModule.h:161
void clear_lom_data()
Clear internal data.
double m_com_th[2]
Monte Carlo thetha of the final state particles in CMS frame.
Definition: eclLOMModule.h:115
TH1D * m_h1FEHits
Store number of events when Forward sector i has signal exceeding Bha-Bha threshold over all events.
Definition: eclLOMModule.h:145
bool m_BEQual_Discr[16][631]
Discriminators values for Quality signal of Backward Endcap.
Definition: eclLOMModule.h:158
double m_BE_Waveform_100ns[16][64]
Waveforms with 100ns sampling for Backward Endcap sectors.
Definition: eclLOMModule.h:109
double m_mcth[2]
Monte Carlo thetha of the final state particles in main frame.
Definition: eclLOMModule.h:112
virtual void beginRun() override
begin run.
Definition: eclLOMModule.cc:95
double m_FE_Pedal[16]
Calculated pedestal values for Forward Endcap.
Definition: eclLOMModule.h:133
TTree * m_testtree
Tree to store output.
Definition: eclLOMModule.h:149
double m_BESum_Amplitude[16]
Calculated amplitudes in running sums of Backward Endcap.
Definition: eclLOMModule.h:130
bool m_isBhabhaPatternFE
Quality signal for Forward endcap.
Definition: eclLOMModule.h:159
bool calculate_FE_quality(int iSample)
Return Quality (topology) flag at sample point, iSample, for Forward Endcap.
double m_FESum_Waveform_10ns[16][631]
Running sum's waveforms with 10ns sampling for Forward Endcap sectors.
Definition: eclLOMModule.h:154
ECLLOMModule()
Constructor.
Definition: eclLOMModule.cc:22
bool m_BESum_Discr[16][631]
Discriminators values for running sums of Backward Endcap.
Definition: eclLOMModule.h:156
bool m_isBhabha
Bha-bha signal for an event.
Definition: eclLOMModule.h:125
double m_BESum_Waveform_10ns[16][631]
Running sum's waveforms with 10ns sampling for Backward Endcap sectors.
Definition: eclLOMModule.h:153
TH2D * m_h2Coin
Store number of coincedencies for i:j sectors (Forward:Backward) over all events.
Definition: eclLOMModule.h:140
double m_mcen[2]
Monte Carlo energy of the final state particles in main frame.
Definition: eclLOMModule.h:111
double m_FE_Amplitude[16]
Calculated amplitudes in sectors of Forward Endcap.
Definition: eclLOMModule.h:129
void calculate_coincidence(int iSample)
Calculates Coincidence Matrix at sample point, iSample.
int m_SumCoincidenceCounterMatrix[16][16]
Stores number of concidences between running sums discriminators in i:j sectors (Forward:Backward).
Definition: eclLOMModule.h:164
double m_BE_Amplitude[16]
Calculated amplitudes in sectors of Backward Endcap.
Definition: eclLOMModule.h:128
double m_com_en[2]
Monte Carlo energy of the final state particles in CMS frame.
Definition: eclLOMModule.h:114
double m_FE_Waveform_10ns[16][631]
Waveforms with 10ns sampling for Forward Endcap sectors.
Definition: eclLOMModule.h:152
double m_BE_Pedal[16]
Calculated pedestal values for Backward Endcap.
Definition: eclLOMModule.h:132
bool m_includeInnerFE
Flag to include Inner part of the Forward Endcap.
Definition: eclLOMModule.h:104
bool m_saveSignal
Flag to save signal wavefroms into file.
Definition: eclLOMModule.h:105
StoreArray< MCParticle > m_MCParticles
MC particles.
Definition: eclLOMModule.h:119
void calculate_amplitudes()
Calculates amplitude [GeV] in an event for each sector.
double m_thresholdBkg
Threshold [GeV] on signal when sector considered as lighted.
Definition: eclLOMModule.h:102
std::string m_testFileName
Name of file to save output.
Definition: eclLOMModule.h:99
double m_com_ph[2]
Monte Carlo phi of the final state particles in CMS frame.
Definition: eclLOMModule.h:116
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:216
Digitize result.
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.
STL namespace.