10 #include <ecl/modules/eclLOM/eclLOMModule.h>
13 #include <Math/Boost.h>
22 ECLLOMModule::ECLLOMModule()
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"));
34 for (
int i = 0; i < 16; i++) {
35 for (
int j = 0; j < 16; j++) {
36 m_CoincidenceCounterMatrix[i][j] = 0;
37 m_SumCoincidenceCounterMatrix[i][j] = 0;
42 ECLLOMModule::~ECLLOMModule()
47 void ECLLOMModule::initialize()
49 if (!(m_MCParticles.isRequired() && m_TrgEclWaveforms.isRequired())) {
52 B2FATAL(
"Not all collections found, exiting processing");
55 m_testfile =
new TFile(m_testFileName.c_str(),
"RECREATE");
56 m_testtree =
new TTree(
"lom_tree",
"");
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");
68 m_testtree->Branch(
"Bhabha", &m_isBhabha);
69 m_testtree->Branch(
"BhNum", &m_BhNum);
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");
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");
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");
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);
95 void ECLLOMModule::beginRun()
99 void ECLLOMModule::event()
104 calculate_discr_output();
105 calculate_amplitudes();
107 for (
int iSample = 300; iSample < 500; iSample++) {
108 m_isBhabhaPatternFE = calculate_FE_quality(iSample);
109 m_isBhabhaPatternBE = calculate_BE_quality(iSample);
110 calculate_coincidence(iSample);
112 for (
int iFESector = 0; iFESector < 16; iFESector++) {
114 int iBESector = (iFESector + 8) % 16;
115 if (m_SumCoincidenceMatrix[iFESector][iBESector] == 1 && m_isBhabhaPatternFE && m_isBhabhaPatternBE) {
126 void ECLLOMModule::endRun()
130 void ECLLOMModule::terminate()
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]);
139 std::cout << std::endl;
146 void ECLLOMModule::get_MCparticles()
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();
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();
171 void ECLLOMModule::get_waveforms()
174 int n_trg_wf = m_TrgEclWaveforms.getEntries();
177 for (
int i = 0; i < n_trg_wf; i++) {
180 int tc_theta_id = TCWaveform->getThetaID();
181 int tc_phi_id = TCWaveform->getPhiID();
183 TCWaveform->fillWaveform(m_wf);
185 int iSectorIndex = (tc_phi_id - 1) / 2;
186 if (tc_theta_id == 1 && !m_includeInnerFE)
continue;
188 for (
int iSample = 0; iSample < 64; iSample++) {
189 if (tc_theta_id <= 3) {
190 m_FE_Waveform_100ns[iSectorIndex][iSample] += m_wf[iSample];
192 if (tc_theta_id == 16 || tc_theta_id == 17) m_BE_Waveform_100ns[iSectorIndex][iSample] += m_wf[iSample];
198 void ECLLOMModule::calculate_discr_output()
200 for (
int iSector = 0; iSector < 16; iSector++) {
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;
207 for (
int iSector = 0; iSector < 16; iSector++) {
208 for (
int iSample = 0; iSample < 63; iSample++) {
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];
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];
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];
222 int TimeOfDiscr = int(m_discrTime / 10);
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];
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;
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;
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;
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;
255 bool ECLLOMModule::calculate_BE_quality(
int iSample)
260 for (
int iBESector = 0; iBESector < 16; iBESector++) {
261 if (m_BEQual_Discr[iBESector][iSample]) {
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);
271 bool ECLLOMModule::calculate_FE_quality(
int iSample)
275 for (
int iFESector = 0; iFESector < 16; iFESector++) {
276 if (m_FEQual_Discr[iFESector][iSample]) {
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);
286 void ECLLOMModule::calculate_coincidence(
int iSample)
288 for (
int iFESector = 0; iFESector < 16; iFESector++) {
289 for (
int iBESector = 0; iBESector < 16; iBESector++) {
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]++;
295 m_CoincidenceMatrix[iFESector][iBESector] = 0;
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]++;
302 m_SumCoincidenceMatrix[iFESector][iBESector] = 0;
308 void ECLLOMModule::clear_lom_data()
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;
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;
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;
344 void ECLLOMModule::calculate_amplitudes()
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];
356 for (
int i = 0; i < 16; i++) {
357 if (m_FESum_Amplitude[i] > m_FESum_MaxAmp) {
358 m_FESum_MaxAmp = m_FESum_Amplitude[i];
361 if (m_BESum_Amplitude[i] > m_BESum_MaxAmp) {
362 m_BESum_MaxAmp = m_BESum_Amplitude[i];
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]);
368 if (m_FE_Amplitude[i] > m_thresholdFE) m_h1FEHits->Fill(i);
369 if (m_BE_Amplitude[i] > m_thresholdBE) m_h1BEHits->Fill(i);
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Abstract base class for different kinds of events.