Belle II Software development
TOPDQMModule.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#include <top/modules/TOPDQM/TOPDQMModule.h>
10#include <top/geometry/TOPGeometryPar.h>
11#include <tracking/dataobjects/ExtHit.h>
12#include <framework/gearbox/Unit.h>
13#include <framework/gearbox/Const.h>
14#include <framework/logging/Logger.h>
15#include <TDirectory.h>
16#include <boost/format.hpp>
17#include <algorithm>
18#include <cmath>
19
20using namespace std;
21using boost::format;
22
23namespace Belle2 {
28
29 using namespace TOP;
30
31 //-----------------------------------------------------------------
33 //-----------------------------------------------------------------
34
35 REG_MODULE(TOPDQM);
36
37 //-----------------------------------------------------------------
38 // Implementation
39 //-----------------------------------------------------------------
40
42 {
43 // set module description (e.g. insert text)
44 setDescription("TOP DQM histogrammer");
46
47 // Add parameters
48 addParam("histogramDirectoryName", m_histogramDirectoryName,
49 "histogram directory in ROOT file", string("TOP"));
50 addParam("momentumCut", m_momentumCut,
51 "momentum cut used to histogram number of photons per track", 0.5);
52 }
53
54
56 {
57 // Create a separate histogram directory and cd into it.
58
59 TDirectory* oldDir = gDirectory;
60 oldDir->mkdir(m_histogramDirectoryName.c_str())->cd();
61
62 // Variables needed for booking
63
64 const auto* geo = TOPGeometryPar::Instance()->getGeometry();
65 m_numModules = geo->getNumModules();
66 m_bunchTimeSep = geo->getNominalTDC().getSyncTimeBase() / 24;
67
68 // Histograms
69
70 m_BoolEvtMonitor = new TH1D("BoolEvtMonitor", "Event synchronization", 2, -0.5, 1.5);
71 m_BoolEvtMonitor->GetYaxis()->SetTitle("number of digits");
72 m_BoolEvtMonitor->GetXaxis()->SetBinLabel(1, "synchronized");
73 m_BoolEvtMonitor->GetXaxis()->SetBinLabel(2, "de-synchronized");
74 m_BoolEvtMonitor->GetXaxis()->SetLabelSize(0.05);
75 m_BoolEvtMonitor->GetXaxis()->SetAlphanumeric();
76 m_BoolEvtMonitor->SetMinimum(0);
77
78 m_window_vs_slot = new TH2F("window_vs_slot", "Asic windows", 16, 0.5, 16.5, 512, 0, 512);
79 m_window_vs_slot->SetXTitle("slot number");
80 m_window_vs_slot->SetYTitle("window number w.r.t reference window");
81 m_window_vs_slot->SetStats(kFALSE);
82 m_window_vs_slot->SetMinimum(0);
83 m_window_vs_slot->GetXaxis()->SetNdivisions(16);
84
85 int nbinsT0 = 75;
86 double rangeT0 = nbinsT0 * m_bunchTimeSep;
87 m_eventT0 = new TH1F("eventT0", "Event T0; event T0 [ns]; events per bin", nbinsT0, -rangeT0 / 2, rangeT0 / 2);
88
89 m_bunchOffset = new TH1F("bunchOffset", "Bunch offset", 100, -m_bunchTimeSep / 2, m_bunchTimeSep / 2);
90 m_bunchOffset->SetXTitle("offset [ns]");
91 m_bunchOffset->SetYTitle("events per bin");
92 m_bunchOffset->SetMinimum(0);
93
94 m_time = new TH1F("goodHitTimes", "Time distribution of good hits", 1000, -20, 80);
95 m_time->SetXTitle("time [ns]");
96 m_time->SetYTitle("hits per bin");
97
98 m_timeBG = new TH1F("goodHitTimesBG", "Time distribution of good hits (background)", 1000, -20, 80);
99 m_timeBG->SetXTitle("time [ns]");
100 m_timeBG->SetYTitle("hits per bin");
101
102 m_signalHits = new TProfile("signalHits", "Number of good hits per track in [0, 50] ns", 16, 0.5, 16.5, 0, 1000);
103 m_signalHits->SetXTitle("slot number");
104 m_signalHits->SetYTitle("hits per track");
105 m_signalHits->SetMinimum(0);
106 m_signalHits->GetXaxis()->SetNdivisions(16);
107
108 m_backgroundHits = new TProfile("backgroundHits", "Number of good hits pet track in [-50, 0] ns", 16, 0.5, 16.5, 0, 1000);
109 m_backgroundHits->SetXTitle("slot number");
110 m_backgroundHits->SetYTitle("hits per track");
111 m_backgroundHits->SetMinimum(0);
112 m_backgroundHits->GetXaxis()->SetNdivisions(16);
113
114 m_trackHits = new TH2F("trackHits", "Number of events w/ and w/o track in the slot", 16, 0.5, 16.5, 2, 0, 2);
115 m_trackHits->SetXTitle("slot number");
116 m_trackHits->SetYTitle("numTracks > 0");
117
118 int nbinsHits = 1000;
119 double xmaxHits = 10000;
120 m_goodHitsPerEventAll = new TH1F("goodHitsPerEventAll", "Number of good hits per event", nbinsHits, 0, xmaxHits);
121 m_goodHitsPerEventAll->GetXaxis()->SetTitle("hits per event");
122 m_goodHitsPerEventAll->GetYaxis()->SetTitle("entries per bin");
123
124 m_badHitsPerEventAll = new TH1F("badHitsPerEventAll", "Number of junk hits per event", nbinsHits, 0, xmaxHits);
125 m_badHitsPerEventAll->GetXaxis()->SetTitle("hits per event");
126 m_badHitsPerEventAll->GetYaxis()->SetTitle("entries per bin");
127
128 int nbinsTDC = 400;
129 double xminTDC = -100;
130 double xmaxTDC = 700;
131 m_goodTDCAll = new TH1F("goodTDCAll", "Raw time distribution of good hits", nbinsTDC, xminTDC, xmaxTDC);
132 m_goodTDCAll->SetXTitle("raw time [samples]");
133 m_goodTDCAll->SetYTitle("hits per sample");
134
135 m_badTDCAll = new TH1F("badTDCAll", "Raw time distribution of junk hits", nbinsTDC, xminTDC, xmaxTDC);
136 m_badTDCAll->SetXTitle("raw time [samples]");
137 m_badTDCAll->SetYTitle("hits per sample");
138
139 m_TOPOccAfterInjLER = new TH1F("TOPOccInjLER", "TOPOccInjLER/Time;Time in #mus;Nhits/Time (#mus bins)", 4000, 0, 20000);
140 m_TOPOccAfterInjHER = new TH1F("TOPOccInjHER", "TOPOccInjHER/Time;Time in #mus;Nhits/Time (#mus bins)", 4000, 0, 20000);
141 m_TOPEOccAfterInjLER = new TH1F("TOPEOccInjLER", "TOPEOccInjLER/Time;Time in #mus;Triggers/Time (#mus bins)", 4000, 0, 20000);
142 m_TOPEOccAfterInjHER = new TH1F("TOPEOccInjHER", "TOPEOccInjHER/Time;Time in #mus;Triggers/Time (#mus bins)", 4000, 0, 20000);
143
144 m_nhitInjLER = new TProfile2D("nhitInjLER",
145 "LER: average Nhits; "
146 "time since injection [msec]; time within beam cycle [syst. clk]; Nhits average",
147 160, 0, 80, 256, 0, 1280, 0, 100000);
148 m_nhitInjHER = new TProfile2D("nhitInjHER",
149 "HER: average Nhits; "
150 "time since injection [msec]; time within beam cycle [syst. clk]; Nhits average",
151 160, 0, 80, 256, 0, 1280, 0, 100000);
152 m_nhitInjLERcut = new TProfile2D("nhitInjLERcut",
153 "LER: average Nhits (Nhits>1000); "
154 "time since injection [msec]; time within beam cycle [syst. clk]; Nhits average",
155 160, 0, 80, 256, 0, 1280, 0, 100000);
156 m_nhitInjHERcut = new TProfile2D("nhitInjHERcut",
157 "HER: average Nhits (Nhits>1000); "
158 "time since injection [msec]; time within beam cycle [syst. clk]; Nhits average",
159 160, 0, 80, 256, 0, 1280, 0, 100000);
160 m_eventInjLER = new TH2F("eventInjLER",
161 "LER: event distribution; "
162 "time since injection [msec]; time within beam cycle [syst. clk]; events per bin",
163 160, 0, 80, 256, 0, 1280);
164 m_eventInjHER = new TH2F("eventInjHER",
165 "HER: event distribution; "
166 "time since injection [msec]; time within beam cycle [syst. clk]; events per bin",
167 160, 0, 80, 256, 0, 1280);
168 m_eventInjLERcut = new TH2F("eventInjLERcut",
169 "LER: event distribution (Nhits>1000); "
170 "time since injection [msec]; time within beam cycle [syst. clk]; events per bin",
171 160, 0, 80, 256, 0, 1280);
172 m_eventInjHERcut = new TH2F("eventInjHERcut",
173 "HER: event distribution (Nhits>1000); "
174 "time since injection [msec]; time within beam cycle [syst. clk]; events per bin",
175 160, 0, 80, 256, 0, 1280);
176
177 for (int i = 0; i < m_numModules; i++) {
178 int module = i + 1;
179 string name, title;
180 TH1F* h1 = 0;
181 TH2F* h2 = 0;
182
183 name = str(format("window_vs_asic_%1%") % (module));
184 title = str(format("Asic windows for slot #%1%") % (module));
185 h2 = new TH2F(name.c_str(), title.c_str(), 64, 0, 64, 512, 0, 512);
186 h2->SetStats(kFALSE);
187 h2->SetXTitle("ASIC number");
188 h2->SetYTitle("window number w.r.t reference window");
189 h2->SetMinimum(0);
190 m_window_vs_asic.push_back(h2);
191
192 name = str(format("good_hits_xy_%1%") % (module));
193 title = str(format("Distribution of good hits for slot #%1%") % (module));
194 h2 = new TH2F(name.c_str(), title.c_str(), 64, 0.5, 64.5, 8, 0.5, 8.5);
195 h2->SetStats(kFALSE);
196 h2->GetXaxis()->SetTitle("pixel column");
197 h2->GetYaxis()->SetTitle("pixel row");
198 h2->SetMinimum(0);
199 m_goodHitsXY.push_back(h2);
200
201 name = str(format("bad_hits_xy_%1%") % (module));
202 title = str(format("Distribution of junk hits for slot #%1%") % (module));
203 h2 = new TH2F(name.c_str(), title.c_str(), 64, 0.5, 64.5, 8, 0.5, 8.5);
204 h2->SetStats(kFALSE);
205 h2->GetXaxis()->SetTitle("pixel column");
206 h2->GetYaxis()->SetTitle("pixel row");
207 h2->SetMinimum(0);
208 m_badHitsXY.push_back(h2);
209
210 name = str(format("good_hits_asics_%1%") % (module));
211 title = str(format("Distribution of good hits for slot #%1%") % (module));
212 h2 = new TH2F(name.c_str(), title.c_str(), 64, 0, 64, 8, 0, 8);
213 h2->SetStats(kFALSE);
214 h2->GetXaxis()->SetTitle("ASIC number");
215 h2->GetYaxis()->SetTitle("ASIC channel");
216 h2->SetMinimum(0);
217 m_goodHitsAsics.push_back(h2);
218
219 name = str(format("bad_hits_asics_%1%") % (module));
220 title = str(format("Distribution of junk hits for slot #%1%") % (module));
221 h2 = new TH2F(name.c_str(), title.c_str(), 64, 0, 64, 8, 0, 8);
222 h2->SetStats(kFALSE);
223 h2->GetXaxis()->SetTitle("ASIC number");
224 h2->GetYaxis()->SetTitle("ASIC channel");
225 h2->SetMinimum(0);
226 m_badHitsAsics.push_back(h2);
227
228 name = str(format("good_TDC_%1%") % (module));
229 title = str(format("Raw time distribution of good hits for slot #%1%") % (module));
230 h1 = new TH1F(name.c_str(), title.c_str(), nbinsTDC, xminTDC, xmaxTDC);
231 h1->GetXaxis()->SetTitle("raw time [samples]");
232 h1->GetYaxis()->SetTitle("hits per sample");
233 m_goodTDC.push_back(h1);
234
235 name = str(format("bad_TDC_%1%") % (module));
236 title = str(format("Raw time distribution of junk hits for slot #%1%") % (module));
237 h1 = new TH1F(name.c_str(), title.c_str(), nbinsTDC, xminTDC, xmaxTDC);
238 h1->GetXaxis()->SetTitle("raw time [samples]");
239 h1->GetYaxis()->SetTitle("hits per sample");
240 m_badTDC.push_back(h1);
241
242 name = str(format("good_timing_%1%") % (module));
243 title = str(format("Time distribution of good hits for slot #%1%") % (module));
244 h1 = new TH1F(name.c_str(), title.c_str(), 200, -20, 80);
245 h1->GetXaxis()->SetTitle("time [ns]");
246 h1->GetYaxis()->SetTitle("hits per time bin");
247 m_goodTiming.push_back(h1);
248
249 name = str(format("good_timing_%1%BG") % (module));
250 title = str(format("Time distribution of good hits (background) for slot #%1%") % (module));
251 h1 = new TH1F(name.c_str(), title.c_str(), 200, -20, 80);
252 h1->GetXaxis()->SetTitle("time [ns]");
253 h1->GetYaxis()->SetTitle("hits per time bin");
254 m_goodTimingBG.push_back(h1);
255
256 name = str(format("good_channel_hits_%1%") % (module));
257 title = str(format("Distribution of good hits for slot #%1%") % (module));
258 int numPixels = geo->getModule(i + 1).getPMTArray().getNumPixels();
259 h1 = new TH1F(name.c_str(), title.c_str(), numPixels, 0, numPixels);
260 h1->GetXaxis()->SetTitle("channel number");
261 h1->GetYaxis()->SetTitle("hits per channel");
262 h1->SetMinimum(0);
263 m_goodChannelHits.push_back(h1);
264
265 name = str(format("bad_channel_hits_%1%") % (module));
266 title = str(format("Distribution of junk hits for slot #%1%") % (module));
267 h1 = new TH1F(name.c_str(), title.c_str(), numPixels, 0, numPixels);
268 h1->GetXaxis()->SetTitle("channel number");
269 h1->GetYaxis()->SetTitle("hits per channel");
270 h1->SetMinimum(0);
271 m_badChannelHits.push_back(h1);
272
273 name = str(format("pulseHeights_%1%") % (module));
274 title = str(format("Average pulse heights for slot #%1%") % (module));
275 auto* prof = new TProfile(name.c_str(), title.c_str(), 32, 0.5, 32.5, 0, 2000);
276 prof->GetXaxis()->SetTitle("PMT number");
277 prof->GetYaxis()->SetTitle("pulse height [ADC counts]");
278 prof->SetMarkerStyle(20);
279 prof->SetMinimum(0);
280 m_pulseHeights.push_back(prof);
281 }
282
283 m_skipProcFlag = new TH2F("skipProcFlag", "Skip processing flag; slot number; flag value", 64, 0.5, 16.5, 2, -0.5, 1.5);
284 m_skipProcFlag->GetXaxis()->SetNdivisions(16);
285 m_skipProcFlag->GetYaxis()->SetNdivisions(2);
286 m_skipProcFlag->SetMinimum(0);
287
288 m_injVetoFlag = new TH2F("injVetoFlag", "Injection veto flag; slot number; flag value", 64, 0.5, 16.5, 2, -0.5, 1.5);
289 m_injVetoFlag->GetXaxis()->SetNdivisions(16);
290 m_injVetoFlag->GetYaxis()->SetNdivisions(2);
291 m_injVetoFlag->SetMinimum(0);
292
293 m_injVetoFlagDiff = new TH1F("injVetoFlagDiff", "Injection veto flags check; ; fraction of events", 2, -0.5, 1.5);
294 m_injVetoFlagDiff->GetXaxis()->SetNdivisions(2);
295 m_injVetoFlagDiff->GetXaxis()->SetBinLabel(1, "flags the same");
296 m_injVetoFlagDiff->GetXaxis()->SetBinLabel(2, "flags differ");
297 m_injVetoFlagDiff->GetXaxis()->SetLabelSize(0.05);
298 m_injVetoFlagDiff->GetXaxis()->SetAlphanumeric();
299 m_injVetoFlagDiff->SetMinimum(0);
300
301 m_PSBypassMode = new TH2F("PSBypassMode", "PS-bypass mode; slot number; mode", 64, 0.5, 16.5, 8, -0.5, 7.5);
302 m_PSBypassMode->GetXaxis()->SetNdivisions(16);
303 m_PSBypassMode->GetYaxis()->SetNdivisions(8);
304 m_PSBypassMode->SetMinimum(0);
305
306 // cd back to root directory
307 oldDir->cd();
308 }
309
311 {
312 // Register histograms (calls back defineHisto)
313
314 REG_HISTOGRAM;
315
316 // register dataobjects
317
318 m_rawFTSWs.isOptional();
319 m_productionEventDebugs.isOptional(); // not in sim
320 m_digits.isRequired();
321 m_recBunch.isOptional();
322 m_timeZeros.isOptional();
323 m_tracks.isOptional();
324 }
325
327 {
328 m_BoolEvtMonitor->Reset();
329 m_window_vs_slot->Reset();
330 m_eventT0->Reset();
331 m_bunchOffset->Reset();
332 m_time->Reset();
333 m_timeBG->Reset();
334 m_signalHits->Reset();
335 m_backgroundHits->Reset();
336 m_trackHits->Reset();
337 m_goodTDCAll->Reset();
338 m_badTDCAll->Reset();
339 m_goodHitsPerEventAll->Reset();
340 m_badHitsPerEventAll->Reset();
341 m_TOPOccAfterInjLER->Reset();
342 m_TOPOccAfterInjHER->Reset();
343 m_TOPEOccAfterInjLER->Reset();
344 m_TOPEOccAfterInjHER->Reset();
345 m_nhitInjLER->Reset();
346 m_nhitInjHER->Reset();
347 m_nhitInjLERcut->Reset();
348 m_nhitInjHERcut->Reset();
349 m_eventInjLER->Reset();
350 m_eventInjHER->Reset();
351 m_eventInjLERcut->Reset();
352 m_eventInjHERcut->Reset();
353
354 for (int i = 0; i < m_numModules; i++) {
355 m_window_vs_asic[i]->Reset();
356 m_goodHitsXY[i]->Reset();
357 m_badHitsXY[i]->Reset();
358 m_goodHitsAsics[i]->Reset();
359 m_badHitsAsics[i]->Reset();
360 m_goodTDC[i]->Reset();
361 m_badTDC[i]->Reset();
362 m_goodTiming[i]->Reset();
363 m_goodTimingBG[i]->Reset();
364 m_goodChannelHits[i]->Reset();
365 m_badChannelHits[i]->Reset();
366 m_pulseHeights[i]->Reset();
367 }
368
369 m_skipProcFlag->Reset();
370 m_injVetoFlag->Reset();
371 m_injVetoFlagDiff->Reset();
372 m_PSBypassMode->Reset();
373 }
374
376 {
377
378 // check if event time is reconstructed; distinguish collision data and cosmics
379
380 bool recBunchValid = false;
381 bool cosmics = false;
382 if (m_recBunch.isValid()) { // collision data
383 recBunchValid = m_recBunch->isReconstructed(); // event time is reconstructed
384 } else if (m_timeZeros.getEntries() == 1) { // cosmics w/ reconstructed event time
385 cosmics = true;
386 m_eventT0->Fill(m_timeZeros[0]->getTime());
387 }
388
389 // fill bunch offset
390
391 if (recBunchValid) {
392 double t0 = m_commonT0->isRoughlyCalibrated() ? m_commonT0->getT0() : 0;
393 double offset = m_recBunch->getCurrentOffset() - t0;
394 offset -= m_bunchTimeSep * lround(offset / m_bunchTimeSep); // wrap around
395 m_bunchOffset->Fill(offset);
396 m_eventT0->Fill(m_recBunch->getTime());
397 }
398
399 // fill event desynchronization
400
401 if (m_digits.getEntries() > 0) {
402 for (const auto& digit : m_digits) {
403 int x = digit.getFirstWindow() != m_digits[0]->getFirstWindow() ? 1 : 0 ;
404 m_BoolEvtMonitor->Fill(x);
405 }
406 }
407
408 // count tracks in the modules and store the momenta
409
410 std::vector<int> numTracks(m_numModules, 0);
411 std::vector<double> trackMomenta(m_numModules, 0);
412 for (const auto& track : m_tracks) {
413 const auto* fitResult = track.getTrackFitResultWithClosestMass(Const::pion);
414 if (not fitResult) continue;
415 int slot = getModuleID(track);
416 if (slot == 0) continue;
417 numTracks[slot - 1]++;
418 trackMomenta[slot - 1] = std::max(trackMomenta[slot - 1], fitResult->getMomentum().R());
419 }
420
421 // count events w/ and w/o track in the slot
422
423 if (recBunchValid or cosmics) {
424 for (size_t i = 0; i < numTracks.size(); i++) {
425 bool hit = numTracks[i] > 0;
426 m_trackHits->Fill(i + 1, hit);
427 }
428 }
429
430 // select modules for counting hits in signal and background time windows
431
432 std::vector<bool> selectedSlots(m_numModules, false);
433 for (size_t i = 0; i < selectedSlots.size(); i++) {
434 selectedSlots[i] = (recBunchValid or cosmics) and numTracks[i] == 1 and trackMomenta[i] > m_momentumCut;
435 }
436
437 // prepare counters
438
439 int nHits_good = 0;
440 int nHits_bad = 0;
441 std::vector<int> numSignalHits(m_numModules, 0);
442 std::vector<int> numBackgroundHits(m_numModules, 0);
443
444 // loop over digits, fill histograms and increment counters
445
446 for (const auto& digit : m_digits) {
447 int slot = digit.getModuleID();
448 if (slot < 1 or slot > m_numModules) {
449 B2ERROR("Invalid slot ID found in TOPDigits: ID = " << slot);
450 continue;
451 }
452 int asic_no = digit.getChannel() / 8;
453 int asic_ch = digit.getChannel() % 8;
454
455 m_window_vs_slot->Fill(digit.getModuleID(), digit.getRawTime() / 64 + 220);
456 m_window_vs_asic[slot - 1]->Fill(asic_no, digit.getRawTime() / 64 + 220);
457
458 if (digit.getHitQuality() != TOPDigit::c_Junk) { // good hits
459 m_goodHitsXY[slot - 1]->Fill(digit.getPixelCol(), digit.getPixelRow());
460 m_goodHitsAsics[slot - 1]->Fill(asic_no, asic_ch);
461 m_goodTDC[slot - 1]->Fill(digit.getRawTime());
462 m_goodTDCAll->Fill(digit.getRawTime());
463 m_goodChannelHits[slot - 1]->Fill(digit.getChannel());
464 m_pulseHeights[slot - 1]->Fill(digit.getPMTNumber(), digit.getPulseHeight());
465 nHits_good++;
466 if (digit.hasStatus(TOPDigit::c_EventT0Subtracted)) {
467 double time = digit.getTime();
468 if (cosmics) time += 10; // move for 10 ns in order to get the complete signal at positive times
469 if (numTracks[slot - 1] > 0) {
470 m_goodTiming[slot - 1]->Fill(time);
471 m_time->Fill(time);
472 } else {
473 m_goodTimingBG[slot - 1]->Fill(time);
474 m_timeBG->Fill(time);
475 }
476 if (selectedSlots[slot - 1] and abs(time) < 50) {
477 if (time > 0) numSignalHits[slot - 1]++;
478 else numBackgroundHits[slot - 1]++;
479 }
480 }
481 } else { // bad hits: FE not valid, pedestal jump, too short or too wide pulses
482 m_badHitsXY[slot - 1]->Fill(digit.getPixelCol(), digit.getPixelRow());
483 m_badHitsAsics[slot - 1]->Fill(asic_no, asic_ch);
484 m_badTDC[slot - 1]->Fill(digit.getRawTime());
485 m_badTDCAll->Fill(digit.getRawTime());
486 m_badChannelHits[slot - 1]->Fill(digit.getChannel());
487 nHits_bad++;
488 }
489 }
490
491 // histogram counters
492
493 m_goodHitsPerEventAll->Fill(nHits_good);
494 m_badHitsPerEventAll->Fill(nHits_bad);
495 for (int slot = 1; slot <= m_numModules; slot++) {
496 if (selectedSlots[slot - 1]) {
497 m_signalHits->Fill(slot, numSignalHits[slot - 1]);
498 m_backgroundHits->Fill(slot, numBackgroundHits[slot - 1]);
499 }
500 }
501
502 // fill injection histograms
503
504 for (auto& it : m_rawFTSWs) {
505 B2DEBUG(29, "TTD FTSW : " << hex << it.GetTTUtime(0) << " " << it.GetTTCtime(0) << " EvtNr " << it.GetEveNo(0) << " Type " <<
506 (it.GetTTCtimeTRGType(0) & 0xF) << " TimeSincePrev " << it.GetTimeSincePrevTrigger(0) << " TimeSinceInj " <<
507 it.GetTimeSinceLastInjection(0) << " IsHER " << it.GetIsHER(0) << " Bunch " << it.GetBunchNumber(0));
508 auto time_clk = it.GetTimeSinceLastInjection(0); // expressed in system clock
509 if (time_clk != 0x7FFFFFFF) {
510 unsigned int nentries = m_digits.getEntries();
511 double time_us = time_clk / 127.0; // 127MHz clock ticks to us, inexact rounding
512 double time_ms = time_us / 1000;
513 double time_1280 = time_clk % 1280; // within beam cycle, expressed in system clock
514 if (it.GetIsHER(0)) {
515 m_TOPOccAfterInjHER->Fill(time_us, nentries);
516 m_TOPEOccAfterInjHER->Fill(time_us);
517 m_nhitInjHER->Fill(time_ms, time_1280, nHits_good);
518 m_eventInjHER->Fill(time_ms, time_1280);
519 if (nHits_good > 1000) {
520 m_nhitInjHERcut->Fill(time_ms, time_1280, nHits_good);
521 m_eventInjHERcut->Fill(time_ms, time_1280);
522 }
523 } else {
524 m_TOPOccAfterInjLER->Fill(time_us, nentries);
525 m_TOPEOccAfterInjLER->Fill(time_us);
526 m_nhitInjLER->Fill(time_ms, time_1280, nHits_good);
527 m_eventInjLER->Fill(time_ms, time_1280);
528 if (nHits_good > 1000) {
529 m_nhitInjLERcut->Fill(time_ms, time_1280, nHits_good);
530 m_eventInjLERcut->Fill(time_ms, time_1280);
531 }
532 }
533 }
534 }
535
536 // fill raw data header histograms
537
538 const auto& feMapper = TOPGeometryPar::Instance()->getFrontEndMapper();
539 double differ = 0;
540 for (const auto& dbg : m_productionEventDebugs) {
541 auto scrodID = dbg.getScrodID();
542 const auto* femap = feMapper.getMap(scrodID);
543 if (not femap) {
544 B2ERROR("No front-end map available for scrodID " << scrodID);
545 continue;
546 }
547 auto slot = femap->getModuleID();
548 auto bs = femap->getBoardstackNumber();
549 double x = slot + bs / 4. - 0.5;
550 m_skipProcFlag->Fill(x, dbg.getSkipProcessingFlag());
551 m_injVetoFlag->Fill(x, dbg.getInjectionVetoFlag());
552 m_PSBypassMode->Fill(x, dbg.getPSBypassMode());
553 if (dbg.getInjectionVetoFlag() != m_productionEventDebugs[0]->getInjectionVetoFlag()) differ = 1;
554 }
555 m_injVetoFlagDiff->Fill(differ);
556
557 }
558
559
560 int TOPDQMModule::getModuleID(const Track& track) const
561 {
562 Const::EDetector myDetID = Const::EDetector::TOP;
563 int pdgCode = std::abs(Const::pion.getPDGCode());
564
565 RelationVector<ExtHit> extHits = track.getRelationsWith<ExtHit>();
566 for (const auto& extHit : extHits) {
567 if (std::abs(extHit.getPdgCode()) != pdgCode) continue;
568 if (extHit.getDetectorID() != myDetID) continue;
569 if (extHit.getCopyID() < 1 or extHit.getCopyID() > m_numModules) continue;
570 return extHit.getCopyID();
571 }
572
573 return 0;
574 }
575
576
578} // end Belle2 namespace
579
EDetector
Enum for identifying the detector components (detector and subdetector).
Definition Const.h:42
static const ChargedStable pion
charged pion particle
Definition Const.h:661
Store one Ext hit as a ROOT object.
Definition ExtHit.h:31
HistoModule()
Constructor.
Definition HistoModule.h:32
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
Class for type safe access to objects that are referred to in relations.
std::vector< TH2F * > m_window_vs_asic
Histograms window w.r.t reference vs.
std::vector< TH1F * > m_badChannelHits
Histograms for bad channel hits.
int m_numModules
number of TOP modules
std::vector< TH2F * > m_goodHitsXY
Histograms (2D) for good hits in pixels.
std::vector< TH2F * > m_goodHitsAsics
Histograms (2D) for good hits in asic channels.
StoreObjPtr< TOPRecBunch > m_recBunch
reconstructed bunch and event T0
TH2F * m_eventInjHER
event distribution (HER injection)
DBObjPtr< TOPCalCommonT0 > m_commonT0
common T0 calibration constants
double m_momentumCut
momentum cut
TH2F * m_injVetoFlag
injection veto flag vs.
std::vector< TH2F * > m_badHitsXY
Histograms (2D) for bad hits in pixels.
TH1F * m_bunchOffset
reconstructed bunch: current offset
TProfile2D * m_nhitInjLER
average number of good digits (LER injection)
TProfile * m_backgroundHits
number of hits in the background time window vs.
StoreArray< TOPProductionEventDebug > m_productionEventDebugs
collection of event debug data
std::vector< TH1F * > m_badTDC
Histograms for TDC distribution of bad hits.
std::vector< TH1F * > m_goodChannelHits
Histograms for good channel hits.
TH1F * m_eventT0
reconstructed event T0
TH2F * m_eventInjLERcut
event distribution after cut (LER injection)
TH2F * m_window_vs_slot
Histogram window w.r.t reference vs.
TH1F * m_injVetoFlagDiff
check if injection veto flags differ in the event
std::vector< TH1F * > m_goodTDC
Histograms for TDC distribution of good hits.
std::vector< TH1F * > m_goodTimingBG
Histograms for timing distribution of good hits (background)
TH1F * m_timeBG
time distribution of good hits (background)
TH2F * m_eventInjLER
event distribution (LER injection)
StoreArray< RawFTSW > m_rawFTSWs
Input array for DAQ Status.
TH1F * m_TOPEOccAfterInjHER
Histogram for Nr Entries (=Triggrs) for normalization after HER injection.
std::string m_histogramDirectoryName
histogram directory in ROOT file
TProfile2D * m_nhitInjLERcut
average number of good digits after cut (LER injection)
StoreArray< Track > m_tracks
collection of tracks
TH1F * m_badHitsPerEventAll
Number of bad hits per event (all slots)
TH1F * m_TOPEOccAfterInjLER
Histogram for Nr Entries (=Triggrs) for normalization after LER injection.
std::vector< TH1F * > m_goodTiming
Histograms for timing distribution of good hits.
TH1F * m_time
time distribution of good hits
TH1F * m_goodTDCAll
TDC distribution of good hits (all slots)
TProfile2D * m_nhitInjHER
average number of good digits (HER injection)
TH1F * m_TOPOccAfterInjLER
Histogram Ndigits after LER injection.
StoreArray< TOPDigit > m_digits
collection of digits
std::vector< TProfile * > m_pulseHeights
Pulse heights of good hits.
TH2F * m_eventInjHERcut
event distribution after cut (HER injection)
StoreArray< TOPTimeZero > m_timeZeros
reconstructed event T0 in case of cosmics
TProfile * m_signalHits
number of hits in the signal time window vs.
TH2F * m_trackHits
counting events w/ and w/o track in the slot vs.
TH2F * m_PSBypassMode
PS-bypass mode vs.
double m_bunchTimeSep
bunch separation time
TH1D * m_BoolEvtMonitor
Event desynchronization monitoring.
TProfile2D * m_nhitInjHERcut
average number of good digits after cut (HER injection)
std::vector< TH2F * > m_badHitsAsics
Histograms (2D) for bad hits in asic channels.
TH1F * m_TOPOccAfterInjHER
Histogram Ndigits after HER injection.
TH1F * m_badTDCAll
TDC distribution of bad hits (all slots)
TH1F * m_goodHitsPerEventAll
Number of good hits per event (all slots)
TH2F * m_skipProcFlag
skip processing flag vs.
const TOPGeometry * getGeometry() const
Returns pointer to geometry object using basf2 units.
static TOPGeometryPar * Instance()
Static method to obtain the pointer to its instance.
const FrontEndMapper & getFrontEndMapper() const
Returns front-end mapper (mapping of SCROD's to positions within TOP modules)
Class that bundles various TrackFitResults.
Definition Track.h:25
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:559
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition Module.h:649
TOPDQMModule()
Constructor.
virtual void initialize() override
Initialize the Module.
virtual void event() override
Event processor.
virtual void beginRun() override
Called when entering a new run.
int getModuleID(const Track &track) const
Returns slot ID of the module that is hit by the track.
virtual void defineHisto() override
Histogram definitions such as TH1(), TH2(), TNtuple(), TTree()....
Abstract base class for different kinds of events.
STL namespace.