Belle II Software  release-08-01-10
TRGECLDQMModule.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 #include <vector>
9 #include <sstream>
10 #include <string>
11 
12 #include <trg/ecl/modules/trgeclDQM/TRGECLDQMModule.h>
13 #include <trg/ecl/TrgEclMapping.h>
14 #include <trg/ecl/TrgEclCluster.h>
15 #include <trg/ecl/TrgEclDataBase.h>
16 
17 #include <framework/datastore/StoreArray.h>
18 
19 #include <TDirectory.h>
20 
21 #include <TH1D.h>
22 #include <TH2D.h>
23 
24 using namespace Belle2;
25 
26 REG_MODULE(TRGECLDQM);
27 
28 
30 {
31 
32 
33  setDescription("DQM for ECL Trigger system");
35 
36  TCId.clear();
37  TCHitWin.clear();
38  TCEnergy.clear();
39  TCTiming.clear();
40  RevoFAM.clear();
41  RevoTrg.clear();
42  FineTiming.clear();
43 
44  addParam("ECLNTC_grpclk", m_grpclknum, "The number of clocks to group N(TC)", 2);
45 }
46 
47 
49 {
50 
51 }
52 
53 
55 {
56  TDirectory* oldDir = gDirectory;
57  TDirectory* dirDQM = (TDirectory*)gDirectory->Get("TRG");
58  if (!dirDQM) {
59  dirDQM = oldDir->mkdir("TRG");
60  }
61  dirDQM->cd();
62 
63  h_TCId = new TH1D("h_TCId", "[TRGECL] Hit TC ID", 578, 0, 578);
64  h_TCthetaId = new TH1D("h_TCthetaId", "[TRGECL] Hit TC #theta ID", 19, 0, 19);
65  h_TCphiId_FWD = new TH1D("h_TCphiId_FWD", "[TRGECL] Hit TC #phi ID in FWD", 34, 0, 34);
66  h_TCphiId_BR = new TH1D("h_TCphiId_BR", "[TRGECL] Hit TC #phi ID in BR", 38, 0, 38);
67  h_TCphiId_BWD = new TH1D("h_TCphiId_BWD", "[TRGECL] Hit TC #phi ID in BWD", 34, 0, 34);
68  h_TotalEnergy = new TH1D("h_TotalEnergy", "[TRGECL] Total TC Energy (ADC)", 100, 0, 3000);
69  h_TCEnergy = new TH1D("h_TCEnergy", "[TRGECL] TC Energy (ADC)", 100, 0, 1500);
70  h_Narrow_TotalEnergy = new TH1D("h_Narrow_TotalEnergy", "[TRGECL] Total TC Energy (ADC)", 100, 0, 500);
71  h_Narrow_TCEnergy = new TH1D("h_Narrow_TCEnergy", "[TRGECL] TC Energy (ADC)", 100, 0, 100);
72  h_n_TChit_event = new TH1D("h_n_TChit_event", "[TRGECL] N(TC) ", 50, 0, 50);
73  h_n_TChit_clean = new TH1D("h_n_TChit_clean", "[TRGECL] N(TC) (Injection BG Clean)", 300, 0, 300);
74  h_n_TChit_injHER = new TH1D("h_n_TChit_injHER", "[TRGECL] N(TC) (HER Injection BG)", 300, 0, 300);
75  h_n_TChit_injLER = new TH1D("h_n_TChit_injLER", "[TRGECL] N(TC) (LER Injection BG)", 300, 0, 300);
76  h_nTChit_injtime = new TH2D("h_nTChit_injtime", "[TRGECL] N(TC) vs. Time since injection", 201, 0, 200, 100, 0, 50);
77  h_Cluster = new TH1D("h_Cluster", "[TRGECL] N(Cluster) ", 20, 0, 20);
78  h_TCTiming = new TH1D("h_TCTiming", "[TRGECL] TC Timing (ns)", 100, 3010, 3210);
79  h_TRGTiming = new TH1D("h_TRGTiming", "[TRGECL] TRG Timing (ns)", 100, 3010, 3210);
80  h_Cal_TCTiming = new TH1D("h_Cal_TCTiming", "[TRGECL] Cal TC Timing (ns)", 100, -400, 400);
81  h_Cal_TRGTiming = new TH1D("h_Cal_TRGTiming", "[TRGECL] TRG Timing (ns)", 100, -400, 400);
82  h_ECL_TriggerBit = new TH1D("h_ECL_TriggerBit", "[TRGECL] ECL Trigger Bit", 29, 0, 29);
83  h_Cluster_Energy_Sum = new TH1D("h_Cluster_Energy_Sum", "[TRGECL] Energy Sum of 2 Clusters (ADC)", 300, 0, 3000);
84 
85  h_nTChit_injtime->GetXaxis()->SetTitle("The number of TC hits");
86  h_nTChit_injtime->GetYaxis()->SetTitle("Time since injection [ms]");
87 
88  std::stringstream ntcclk_titlebase_ss;
89  ntcclk_titlebase_ss << "[TRGECL] N(TC_" << m_grpclknum << "_clk)";
90  std::string ntcclk_titlebase = ntcclk_titlebase_ss.str();
91  h_n_TChit_event_clkgrp = new TH1D("h_n_TChit_event_clkgrp", ntcclk_titlebase.c_str(), 50, 0, 50);
92  h_n_TChit_clean_clkgrp = new TH1D("h_n_TChit_clean_clkgrp", (ntcclk_titlebase + " (Injection BG Clean)").c_str(),
93  300, 0, 300);
94  h_n_TChit_injHER_clkgrp = new TH1D("h_n_TChit_injHER_clkgrp", (ntcclk_titlebase + " (HER Injection BG)").c_str(),
95  300, 0, 300);
96  h_n_TChit_injLER_clkgrp = new TH1D("h_n_TChit_injLER_clkgrp", (ntcclk_titlebase + " (LER Injection BG)").c_str(),
97  300, 0, 300);
98  h_nTChit_injtime_clkgrp = new TH2D("h_nTChit_injtime_clkgrp", (ntcclk_titlebase + " vs. Time since injection").c_str(), 201, 0, 200,
99  100, 0, 50);
100 
101  const std::vector<std::string> titleArr = {
102  "Forward", "Barrel", "Backward"
103  };
104  const std::vector<std::string> nameArr = {
105  "FWD", "BRL", "BWD"
106  };
107 
108  for (int i = 0; i < 3; i++) {
109  std::stringstream namebase_ss;
110  std::stringstream titlebase_ss;
111 
112  namebase_ss << "h_n_TChit_" << nameArr[i] << "_";
113  titlebase_ss << "[TRGECL] N_" << titleArr[i] << "(TC_" << m_grpclknum << "clk)";
114 
115  h_n_TChit_part_event_clkgrp[i] = new TH1D((namebase_ss.str() + "event_clkgrp").c_str(), titlebase_ss.str().c_str(), 50, 0, 50);
116  h_n_TChit_part_clean_clkgrp[i] = new TH1D((namebase_ss.str() + "clean_clkgrp").c_str(),
117  (titlebase_ss.str() + " (Injection BG Clean)").c_str(), 300, 0, 300);
118  h_n_TChit_part_injHER_clkgrp[i] = new TH1D((namebase_ss.str() + "injHER_clkgrp").c_str(),
119  (titlebase_ss.str() + " (HER Injection BG)").c_str(), 300, 0, 300);
120  h_n_TChit_part_injLER_clkgrp[i] = new TH1D((namebase_ss.str() + "injLER_clkgrp").c_str(),
121  (titlebase_ss.str() + " (LER Injection BG)").c_str(), 300, 0, 300);
122  h_nTChit_part_injtime_clkgrp[i] = new TH2D((namebase_ss.str() + "injtime_clkgrp").c_str(),
123  (titlebase_ss.str() + " vs. Time since injection").c_str(), 201, 0, 200, 100, 0, 50);
124  }
125 
126  const char* label[44] = {"Hit", "Timing Source(FWD)", "Timing Source(BR)", "Timing Source(BWD)", "physics Trigger", "2D Bhabha Veto", "3D Bhabha veto", "3D Bhabha Selection", "E Low", "E High", "E LOM", "Cluster Overflow", "Low multi bit 0", "Low multi bit 1", "Low multi bit 2", "Low multi bit 3", "Low multi bit 4", "Low multi bit 5", "Low multi bit 6", "Low multi bit 7", "Low multi bit 8", "Low multi bit 9", "Low multi bit 10", "Low multi bit 11", "Low multi bit 12", "Low multi bit 13", "mumu bit", "prescale bit", "ECL burst bit", "2D Bhabha bit 1", "2D Bhabha bit 2", "2D Bhabha bit 3", "2D Bhabha bit 4", "2D Bhabha bit 5", "2D Bhabha bit 6", "2D Bhabha bit 7", "2D Bhabha bit 8", "2D Bhabha bit 9", "2D Bhabha bit 10", "2D Bhabha bit 11", "2D Bhabha bit 12", "2D Bhabha bit 13", "2D Bhabha bit 14"};
127 
128  for (int j = 0; j < 29; j++) {
129  h_ECL_TriggerBit->GetXaxis()-> SetBinLabel(j + 1, label[j]);
130  }
131  h_ECL_TriggerBit->SetStats(0);
132 
133  oldDir->cd();
134 }
135 
136 
138 {
139  if (m_grpclknum <= 0 || m_grpclknum >= 9) {
140  B2WARNING("The number of grouping clocks for N(TC) is wrong (It should be 1 <= clk_num <= 8). The value will set to 2.");
141  m_grpclknum = 2;
142  }
143 
144  // calls back the defineHisto() function, but the HistoManager module has to be in the path
145  REG_HISTOGRAM
146 
147  trgeclHitArray.registerInDataStore();
148  trgeclEvtArray.registerInDataStore();
149  trgeclCluster.registerInDataStore();
150  trgeclSumArray.registerInDataStore();
151 }
152 
154 {
155 }
156 
158 {
159  // delete h_TCId;
160 }
161 
163 {
164  TCId.clear();
165  TCHitWin.clear();
166  TCEnergy.clear();
167  TCTiming.clear();
168  RevoFAM.clear();
169  RevoTrg.clear();
170  FineTiming.clear();
171 
172  // StoreArray<TRGECLUnpackerStore> trgeclHitArray;
173  /* cppcheck-suppress variableScope */
174  double HitTiming;
175  /* cppcheck-suppress variableScope */
176  double HitEnergy;
177  double HitRevoFam = 0;
178  double HitRevoTrg = 0;
179  double HitFineTiming = 0;
180  double HitRevoEvtTiming = 0;
181  double HitCalTiming = 0;
182  int CheckSum = 0;
183 
184  for (int iii = 0; iii < trgeclEvtArray.getEntries(); iii++) {
185  TRGECLUnpackerEvtStore* aTRGECLUnpackerEvtStore = trgeclEvtArray[iii];
186 
187  HitFineTiming = aTRGECLUnpackerEvtStore -> getEvtTime();
188  HitRevoTrg = aTRGECLUnpackerEvtStore -> getL1Revo();
189  HitRevoEvtTiming = aTRGECLUnpackerEvtStore -> getEvtRevo();
190  CheckSum = aTRGECLUnpackerEvtStore -> getEvtExist() ;
191 
192 
193  RevoTrg.push_back(HitRevoTrg);
194 
195 
196 
197  }
198  if (CheckSum == 0) {return;}
199 
200 
201 
202  for (int ii = 0; ii < trgeclHitArray.getEntries(); ii++) {
203  TRGECLUnpackerStore* aTRGECLUnpackerStore = trgeclHitArray[ii];
204  int TCID = (aTRGECLUnpackerStore->getTCId());
205  int hit_win = aTRGECLUnpackerStore -> getHitWin();
206  HitEnergy = aTRGECLUnpackerStore -> getTCEnergy();
207  HitTiming = aTRGECLUnpackerStore ->getTCTime();
208 
209  if (TCID < 1 || TCID > 576 || HitEnergy == 0) {continue;}
210  if (!(hit_win == 3 || hit_win == 4)) {continue;}
211  HitCalTiming = aTRGECLUnpackerStore ->getTCCALTime() ;
212  HitRevoFam = aTRGECLUnpackerStore-> getRevoFAM() ;
213 
214  TCId.push_back(TCID);
215  TCHitWin.push_back(hit_win);
216  TCEnergy.push_back(HitEnergy);
217  TCTiming.push_back(HitTiming);
218  RevoFAM.push_back(HitRevoFam);
219  FineTiming.push_back(HitCalTiming);
220  }
221  //
222  //
223  if (TCId.size() == 0) {return;}
224 
225  /* cppcheck-suppress variableScope */
226  int phy;
227  /* cppcheck-suppress variableScope */
228  int b1;
229  /* cppcheck-suppress variableScope */
230  int b2v;
231  /* cppcheck-suppress variableScope */
232  int b2s;
233  /* cppcheck-suppress variableScope */
234  int mu;
235  /* cppcheck-suppress variableScope */
236  int pre;
237  /* cppcheck-suppress variableScope */
238  int clover;
239  /* cppcheck-suppress variableScope */
240  int tsource;
241  /* cppcheck-suppress variableScope */
242  int b1type;
243  /* cppcheck-suppress variableScope */
244  int etot;
245  /* cppcheck-suppress variableScope */
246  int vlm;
247  /* cppcheck-suppress variableScope */
248  int eclburst;
249  // int s_hit_win= 0;
250  std::vector<int> trgbit ;
251  trgbit.resize(44, 0);
252  for (int iii = 0; iii < trgeclSumArray.getEntries(); iii++) {
253  TRGECLUnpackerSumStore* aTRGECLUnpackerSumStore = trgeclSumArray[iii];
254 
255  tsource = aTRGECLUnpackerSumStore ->getTimeType();
256  phy = aTRGECLUnpackerSumStore ->getPhysics();
257  b1 = aTRGECLUnpackerSumStore ->get2DBhabha();
258  b1type = aTRGECLUnpackerSumStore -> getBhabhaType();
259  b2v = aTRGECLUnpackerSumStore -> get3DBhabhaV();
260  b2s = aTRGECLUnpackerSumStore -> get3DBhabhaS() ;
261  etot = aTRGECLUnpackerSumStore -> getEtotType();
262  clover = aTRGECLUnpackerSumStore -> getICNOver();
263  vlm = aTRGECLUnpackerSumStore -> getLowMulti();
264  mu = aTRGECLUnpackerSumStore -> getMumu();
265  pre = aTRGECLUnpackerSumStore -> getPrescale();
266  eclburst = aTRGECLUnpackerSumStore -> getECLBST();
267 
268  //
269  trgbit[0] = 1;
270  trgbit[1] = tsource & 0x1;
271  trgbit[2] = (tsource >> 1) & 0x1;
272  trgbit[3] = (tsource >> 2) & 0x1;
273  trgbit[4] = phy;
274  trgbit[5] = b1;
275  trgbit[6] = b2v;
276  trgbit[7] = b2s;
277  trgbit[8] = etot & 0x1;
278  trgbit[9] = (etot >> 1) & 0x1;
279  trgbit[10] = (etot >> 2) & 0x1;
280  trgbit[11] = clover;
281 
282  for (int j = 0; j < 14; j++) {
283  trgbit[12 + j] = (vlm >> j) & 0x1;
284  }
285 
286  trgbit[26] = mu;
287  trgbit[27] = pre;
288  trgbit[28] = eclburst;
289 
290  trgbit[29] = b1type & 0x1;
291  trgbit[30] = (b1type >> 1) & 0x1;
292  trgbit[31] = (b1type >> 2) & 0x1;
293  trgbit[32] = (b1type >> 3) & 0x1;
294  trgbit[33] = (b1type >> 4) & 0x1;
295  trgbit[34] = (b1type >> 5) & 0x1;
296  trgbit[35] = (b1type >> 6) & 0x1;
297  trgbit[36] = (b1type >> 7) & 0x1;
298  trgbit[37] = (b1type >> 8) & 0x1;
299  trgbit[38] = (b1type >> 9) & 0x1;
300  trgbit[39] = (b1type >> 10) & 0x1;
301  trgbit[40] = (b1type >> 11) & 0x1;
302  trgbit[41] = (b1type >> 12) & 0x1;
303  trgbit[42] = (b1type >> 13) & 0x1;
304  trgbit[43] = (b1type >> 14) & 0x1;
305 
306 
307  }
308 
309  for (int j = 0; j < 29; j++) {
310  if (trgbit[j] == 0x1) {h_ECL_TriggerBit->Fill(j, 1);}
311  }
312 
313 
314  //----------------------
315  //Clustering
316  //----------------------
317  //
318 
319  TrgEclCluster _TCCluster ;
320  _TCCluster.setICN(TCId, TCEnergy, TCTiming);
321 
322  int c = _TCCluster.getNofCluster();
323  h_Cluster->Fill(c);
324  std::vector<double> ClusterTiming;
325  std::vector<double> ClusterEnergy;
326  std::vector<int> MaxTCId;
327  ClusterTiming.clear();
328  ClusterEnergy.clear();
329  MaxTCId.clear();
330 
331  for (int iii = 0; iii < trgeclCluster.getEntries(); iii++) {
332  TRGECLCluster* aTRGECLCluster = trgeclCluster[iii];
333  int maxTCId = aTRGECLCluster ->getMaxTCId();
334  double clusterenergy = aTRGECLCluster ->getEnergyDep();
335  double clustertiming = aTRGECLCluster -> getTimeAve();
336  ClusterTiming.push_back(clustertiming);
337  ClusterEnergy.push_back(clusterenergy);
338  MaxTCId.push_back(maxTCId);
339  }
340 
341 
342  std::vector<double> maxClusterEnergy;
343  std::vector<double> maxClusterTiming;
344  std::vector<int> maxCenterTCId;
345  maxClusterTiming.clear();
346  maxClusterEnergy.clear();
347  maxCenterTCId.clear();
348 
349  maxClusterEnergy.resize(2, 0.0);
350  maxClusterTiming.resize(2, 0.0);
351  maxCenterTCId.resize(2, 0.0);
352  const int cl_size = ClusterEnergy.size();
353  for (int icl = 0; icl < cl_size; icl++) {
354  if (maxClusterEnergy[0] < ClusterEnergy[icl]) {
355  maxClusterEnergy[0] = ClusterEnergy[icl];
356  maxClusterTiming[0] = ClusterTiming[icl];
357  maxCenterTCId[0] = MaxTCId[icl];
358  } else if (maxClusterEnergy[1] < ClusterEnergy[icl]) {
359  maxClusterEnergy[1] = ClusterEnergy[icl];
360  maxClusterTiming[1] = ClusterTiming[icl];
361  maxCenterTCId[1] = MaxTCId[icl];
362 
363  }
364 
365  }
366  TrgEclDataBase _database;
367 
368  std::vector<double> _3DBhabhaThreshold;
369  _3DBhabhaThreshold = {30, 45}; // /10 MeV
370 
371 
372  bool BtoBFlag = false;
373  bool BhabhaFlag = false;
374  int lut1 = _database.Get3DBhabhaLUT(maxCenterTCId[0]);
375  int lut2 = _database.Get3DBhabhaLUT(maxCenterTCId[1]);
376  int energy1 = 15 & lut1;
377  int energy2 = 15 & lut2;
378  lut1 >>= 4;
379  lut2 >>= 4;
380  int phi1 = 511 & lut1;
381  int phi2 = 511 & lut2;
382  lut1 >>= 9;
383  lut2 >>= 9;
384  int theta1 = lut1;
385  int theta2 = lut2;
386  int dphi = abs(phi1 - phi2);
387  if (dphi > 180) {dphi = 360 - dphi;}
388  int thetaSum = theta1 + theta2;
389  if (dphi > 160 && thetaSum > 165 && thetaSum < 190) {BtoBFlag = true;}
390 
391  if ((maxClusterEnergy[0] * 0.1) > _3DBhabhaThreshold[0] * energy1
392  && (maxClusterEnergy[1] * 0.1) > _3DBhabhaThreshold[0] * (energy2)
393  && ((maxClusterEnergy[0] * 0.1) > _3DBhabhaThreshold[1] * energy1
394  || (maxClusterEnergy[1] * 0.1) > _3DBhabhaThreshold[1] * (energy2))) {
395  if (BtoBFlag) {BhabhaFlag = true;}
396  }
397 
398 
399  if (BhabhaFlag) {
400  h_Cluster_Energy_Sum -> Fill((maxClusterEnergy[0] + maxClusterEnergy[1]) / 5.25);
401  }
402 
403 
404  const int NofTCHit = TCId.size();
405 
406  int nTCHitPerClk_total[8] = {0};
407  int nTCHitPerClk_part[3][8] = {0};
408  double totalEnergy = 0;
409  TrgEclMapping* a = new TrgEclMapping();
410  double max = 0;
411  double caltrgtiming = 0;
412  double diff = -1;
413  bool isHER;
414 
415  diff = m_trgTime->getTimeSinceLastInjectionInMicroSeconds() / 1000.;
416  isHER = m_trgTime->isHER();
417 
418  for (int ihit = 0; ihit < NofTCHit ; ihit ++) {
419  h_TCId -> Fill(TCId[ihit]);
420  h_TCthetaId -> Fill(a -> getTCThetaIdFromTCId(TCId[ihit]));
421  {
422  if (a->getTCThetaIdFromTCId(TCId[ihit]) < 4) {
423  h_TCphiId_FWD -> Fill(a->getTCPhiIdFromTCId(TCId[ihit]));
424  } else if (a->getTCThetaIdFromTCId(TCId[ihit]) > 3 && a->getTCThetaIdFromTCId(TCId[ihit]) < 16) {
425  h_TCphiId_BR -> Fill(a->getTCPhiIdFromTCId(TCId[ihit]));
426  } else {
427  h_TCphiId_BWD -> Fill(a->getTCPhiIdFromTCId(TCId[ihit]));
428 
429  }
430  }
431  h_TCEnergy -> Fill(TCEnergy[ihit]);
432  h_Narrow_TCEnergy -> Fill(TCEnergy[ihit]);
433  h_Cal_TCTiming -> Fill(FineTiming[ihit]);
434 
435  if (max < TCEnergy[ihit]) {
436  max = TCEnergy[ihit];
437  caltrgtiming = FineTiming[ihit];
438  }
439 
440  totalEnergy += TCEnergy[ihit];
441  double timing = 8 * HitRevoTrg - (128 * RevoFAM[ihit] + TCTiming[ihit]);
442  if (timing < 0) {timing = timing + 10240;}
443  h_TCTiming->Fill(timing);
444  nTCHitPerClk_total[TCHitWin[ihit]]++;
445 
446  if (TCId[ihit] <= 80) {
447  nTCHitPerClk_part[0][TCHitWin[ihit]]++;
448  } else if (80 < TCId[ihit] && TCId[ihit] <= 512) {
449  nTCHitPerClk_part[1][TCHitWin[ihit]]++;
450  } else {
451  nTCHitPerClk_part[2][TCHitWin[ihit]]++;
452  }
453  }
454 
455  const double revotime_in_us = 5.120 / m_hwclkdb->getAcceleratorRF();
456  int quotient;
457  double running_in_us, diff_in_us;
458 
459  diff_in_us = diff * 1000.;
460  quotient = diff_in_us / revotime_in_us;
461  running_in_us = diff_in_us - quotient * revotime_in_us;
462 
463  bool cond_clean, cond_injHER, cond_injLER;
464 
465  cond_clean = (6 < running_in_us && running_in_us < 8) && (50 < diff && diff < 70);
466 
467  cond_injHER = isHER && ((diff < 0.5) || ((diff < 20) && (2 < running_in_us && running_in_us < 3)));
468  cond_injLER = !isHER && ((diff < 0.5) || ((diff < 20) && (1 < running_in_us && running_in_us < 2)));
469 
470  h_n_TChit_event -> Fill(NofTCHit);
471  h_nTChit_injtime->Fill(NofTCHit, diff);
472 
473  if (cond_clean) {
474  h_n_TChit_clean->Fill(NofTCHit);
475  } else if (cond_injHER) {
476  h_n_TChit_injHER->Fill(NofTCHit);
477  } else if (cond_injLER) {
478  h_n_TChit_injLER->Fill(NofTCHit);
479  }
480 
481  for (int iclk = 0; iclk < 8 - (m_grpclknum - 1); iclk++) {
482  int group_tcnum_total = 0;
483  int group_tcnum_part[3] = {0};
484  for (int igrp = 0; igrp < m_grpclknum; igrp++) {
485  group_tcnum_total += nTCHitPerClk_total[iclk + igrp];
486  for (int ipart = 0; ipart < 3; ipart++) {
487  group_tcnum_part[ipart] += nTCHitPerClk_part[ipart][iclk + igrp];
488  }
489  }
490 
491  h_n_TChit_event_clkgrp->Fill(group_tcnum_total);
492  h_nTChit_injtime_clkgrp->Fill(group_tcnum_total, diff);
493 
494  if (cond_clean) {
495  h_n_TChit_clean_clkgrp->Fill(group_tcnum_total);
496  } else if (cond_injHER) {
497  h_n_TChit_injHER_clkgrp->Fill(group_tcnum_total);
498  } else if (cond_injLER) {
499  h_n_TChit_injLER_clkgrp->Fill(group_tcnum_total);
500  }
501 
502  for (int ipart = 0; ipart < 3; ipart++) {
503  h_n_TChit_part_event_clkgrp[ipart]->Fill(group_tcnum_part[ipart]);
504  h_nTChit_part_injtime_clkgrp[ipart]->Fill(group_tcnum_part[ipart], diff);
505 
506  if (cond_clean) {
507  h_n_TChit_part_clean_clkgrp[ipart]->Fill(group_tcnum_part[ipart]);
508  } else if (cond_injHER) {
509  h_n_TChit_part_injHER_clkgrp[ipart]->Fill(group_tcnum_part[ipart]);
510  } else if (cond_injLER) {
511  h_n_TChit_part_injLER_clkgrp[ipart]->Fill(group_tcnum_part[ipart]);
512  }
513  }
514  }
515 
516  double trgtiming = 8 * HitRevoTrg - (128 * HitRevoEvtTiming + HitFineTiming);
517 
518  if (trgtiming < 0) {trgtiming = trgtiming + 10240;}
519  h_TRGTiming -> Fill(trgtiming);
520  h_Cal_TRGTiming -> Fill(caltrgtiming);
521  h_TotalEnergy -> Fill(totalEnergy);
522  h_Narrow_TotalEnergy -> Fill(totalEnergy);
523 
524  // usleep(100);
525 }
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:208
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:80
Example Detector.
Definition: TRGECLCluster.h:24
double getEnergyDep() const
The method to get deposited energy.
int getMaxTCId() const
The method to get the Maximum(center) TC id.
TH1 * h_n_TChit_part_injHER_clkgrp[3]
N of TC Hit / events in the injection BG clean region vs. time since injection per two ETM clocks for...
TH1 * h_n_TChit_clean
N of TC Hit / event vs. time since injection.
TH1 * h_n_TChit_injLER_clkgrp
N of TC Hit / events in the HER injection BG region vs. time since injection per two ETM clocks.
std::vector< double > TCEnergy
Hit TC Energy.
TH1 * h_TCId
TCId histogram.
TH1 * h_TCphiId_BWD
TCphiId histogram.
TH1 * h_Narrow_TotalEnergy
Total Energy on narrow range.
TH1 * h_TotalEnergy
Total Energy.
TH1 * h_n_TChit_clean_clkgrp
N of TC Hit / event vs. time since injection per two ETM clocks.
TH1 * h_Cluster
N of Cluster / event.
virtual void initialize() override
initialize
TH1 * h_n_TChit_part_event_clkgrp[3]
N of TC Hit / event per two ETM clocks.
StoreArray< TRGECLUnpackerStore > trgeclHitArray
Trg ECL Unpakcer TC output.
virtual void event() override
Event.
TH1 * h_Cal_TRGTiming
Event Timing / event.
TH1 * h_n_TChit_injHER
N of TC Hit / events in the injection BG clean region vs. time since injection.
TH1 * h_Narrow_TCEnergy
TC Energy histogram on narrow range.
TH2 * h_nTChit_injtime_clkgrp
N of TC Hit / events in the LER injection BG region vs. time since injection per two ETM clocks.
virtual void endRun() override
End Run.
DBObjPtr< HardwareClockSettings > m_hwclkdb
DB pointerto access the hardware clock information.
TH1 * h_Cluster_Energy_Sum
Energy sum of 2 Top energetic clusters when 3D bhabnha bit on.
virtual void terminate() override
terminate
TH1 * h_n_TChit_event_clkgrp
N of TC Hit / event per two ETM clocks.
virtual ~TRGECLDQMModule()
Destrunctor.
TH1 * h_n_TChit_injHER_clkgrp
N of TC Hit / events in the injection BG clean region vs. time since injection per two ETM clocks.
std::vector< int > TCHitWin
Hit TCHitWin.
TH1 * h_TCphiId_FWD
TCphiId histogram.
TH1 * h_TRGTiming
Event Timing / event.
TH2 * h_nTChit_injtime
N of TC Hit / events in the LER injection BG region vs. time since injection.
TH1 * h_ECL_TriggerBit
ECL Trigger Bit.
std::vector< double > TCTiming
Hit TC Timing.
StoreArray< TRGECLUnpackerEvtStore > trgeclEvtArray
Trg ECL Unpakcer Event output.
TH1 * h_TCthetaId
TCthetaId histogram.
virtual void beginRun() override
begin Run
StoreArray< TRGECLUnpackerSumStore > trgeclSumArray
Trg Ecl Unpacker Summary output.
TH1 * h_n_TChit_part_clean_clkgrp[3]
N of TC Hit / event vs. time since injection per two ETM clocks for each part (FWD,...
TH2 * h_nTChit_part_injtime_clkgrp[3]
N of TC Hit / events in the LER injection BG region vs. time since injection per two ETM clocks for e...
TH1 * h_TCphiId_BR
TCphiId histogram.
TH1 * h_Cal_TCTiming
TC Timing / event.
TH1 * h_TCEnergy
TC Energy.
std::vector< double > FineTiming
Event Timing.
int m_grpclknum
The number of clocks to group the number of TCs of ECLTRG.
std::vector< double > RevoTrg
GDL Revolution Clk.
TH1 * h_TCTiming
TC Timing / event.
TH1 * h_n_TChit_event
N of TC Hit / event.
std::vector< double > RevoFAM
FAM Revolution Clk.
TH1 * h_n_TChit_part_injLER_clkgrp[3]
N of TC Hit / events in the HER injection BG region vs. time since injection per two ETM clocks for e...
TH1 * h_n_TChit_injLER
N of TC Hit / events in the HER injection BG region vs. time since injection.
StoreArray< TRGECLCluster > trgeclCluster
Trg ECL Cluster output.
StoreObjPtr< EventLevelTriggerTimeInfo > m_trgTime
Array to access the FTSW information.
std::vector< int > TCId
Hit TCId.
virtual void defineHisto() override
Define Histogram.
int getTCCALTime() const
The method to get cal timing.
int getTCId() const
The method to get cell id.
int getTCTime() const
The method to get hit average time.
int get2DBhabha() const
The mothod to get 2D Bhabha bit.
int getPhysics() const
The mothod to get Physics bit.
int getTimeType() const
The mothod to get Timing Type.
A Class of ECL Trigger clustering
Definition: TrgEclCluster.h:30
void setICN(const std::vector< int > &)
set ICN for each part(Fw,Br,Bw)
int getNofCluster()
0 : center , 1; upper , 2: right , 3: lower , 4: lower right
Definition: TrgEclCluster.h:95
class TrgEclDataBase;
int Get3DBhabhaLUT(int)
TC CM Phi
A class of TC Mapping.
Definition: TrgEclMapping.h:26
REG_MODULE(arichBtest)
Register the Module.
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
Abstract base class for different kinds of events.