Belle II Software  release-08-01-10
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 
15 using namespace std;
16 using namespace Belle2;
17 using namespace ECL;
18 
19 REG_MODULE(ECLLOM);
20 
21 
22 ECLLOMModule::ECLLOMModule()
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++) {
36  m_CoincidenceCounterMatrix[i][j] = 0;
37  m_SumCoincidenceCounterMatrix[i][j] = 0;
38  }
39  }
40 }
41 
42 ECLLOMModule::~ECLLOMModule()
43 {
44 
45 }
46 
47 void ECLLOMModule::initialize()
48 {
49  if (!(m_MCParticles.isRequired() && m_TrgEclWaveforms.isRequired())) {
50  //Fatal is not neccessary 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 
95 void ECLLOMModule::beginRun()
96 {
97 }
98 
99 void ECLLOMModule::event()
100 {
101  clear_lom_data();
102  get_MCparticles();
103  get_waveforms();
104  calculate_discr_output();
105  calculate_amplitudes();
106  // LOM logic
107  for (int iSample = 300; iSample < 500; iSample++) { //300-500 window where BhaBha is expected
108  m_isBhabhaPatternFE = calculate_FE_quality(iSample);
109  m_isBhabhaPatternBE = calculate_BE_quality(iSample);
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  //coinsidence at first tick
117  m_isBhabha = true;
118  m_BhNum++;
119  }
120  }
121  }
122  m_testtree->Fill();
123  m_evtNum++;
124 }
125 
126 void ECLLOMModule::endRun()
127 {
128 }
129 
130 void ECLLOMModule::terminate()
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 
146 void ECLLOMModule::get_MCparticles()
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 
171 void ECLLOMModule::get_waveforms()
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 
198 void ECLLOMModule::calculate_discr_output()
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 
255 bool ECLLOMModule::calculate_BE_quality(int iSample)
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 
271 bool ECLLOMModule::calculate_FE_quality(int iSample)
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 
286 void ECLLOMModule::calculate_coincidence(int iSample)
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 
308 void ECLLOMModule::clear_lom_data()
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 
344 void ECLLOMModule::calculate_amplitudes()
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++) {
357  if (m_FESum_Amplitude[i] > m_FESum_MaxAmp) {
358  m_FESum_MaxAmp = m_FESum_Amplitude[i];
359  m_FESum_MaxId = i;
360  }
361  if (m_BESum_Amplitude[i] > m_BESum_MaxAmp) {
362  m_BESum_MaxAmp = m_BESum_Amplitude[i];
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 }
Digitize result.
#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.