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