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, 1, 17, 2, 0, 2);
284 m_skipProcFlag->GetXaxis()->SetNdivisions(16);
285 m_skipProcFlag->GetXaxis()->CenterLabels();
286 m_skipProcFlag->GetYaxis()->SetNdivisions(2);
287 m_skipProcFlag->GetYaxis()->SetBinLabel(1, "0");
288 m_skipProcFlag->GetYaxis()->SetBinLabel(2, "1");
289 m_skipProcFlag->GetYaxis()->SetLabelSize(0.05);
290 m_skipProcFlag->GetYaxis()->SetAlphanumeric();
291 m_skipProcFlag->SetMinimum(0);
292
293 m_injVetoFlag = new TH2F("injVetoFlag", "Injection veto flag; slot number; flag value", 64, 1, 17, 2, 0, 2);
294 m_injVetoFlag->GetXaxis()->SetNdivisions(16);
295 m_injVetoFlag->GetXaxis()->CenterLabels();
296 m_injVetoFlag->GetYaxis()->SetNdivisions(2);
297 m_injVetoFlag->GetYaxis()->SetBinLabel(1, "0");
298 m_injVetoFlag->GetYaxis()->SetBinLabel(2, "1");
299 m_injVetoFlag->GetYaxis()->SetLabelSize(0.05);
300 m_injVetoFlag->GetYaxis()->SetAlphanumeric();
301 m_injVetoFlag->SetMinimum(0);
302
303 m_injVetoFlagDiff = new TH1F("injVetoFlagDiff", "Injection veto flags check; ; fraction of events", 2, -0.5, 1.5);
304 m_injVetoFlagDiff->GetXaxis()->SetNdivisions(2);
305 m_injVetoFlagDiff->GetXaxis()->SetBinLabel(1, "flags the same");
306 m_injVetoFlagDiff->GetXaxis()->SetBinLabel(2, "flags differ");
307 m_injVetoFlagDiff->GetXaxis()->SetLabelSize(0.05);
308 m_injVetoFlagDiff->GetXaxis()->SetAlphanumeric();
309 m_injVetoFlagDiff->SetMinimum(0);
310
311 m_PSBypassMode = new TH2F("PSBypassMode", "PS-bypass mode; slot number; mode", 64, 1, 17, 8, 0, 8);
312 m_PSBypassMode->GetXaxis()->SetNdivisions(16);
313 m_PSBypassMode->GetYaxis()->SetNdivisions(8);
314 m_PSBypassMode->GetXaxis()->CenterLabels();
315 m_PSBypassMode->GetYaxis()->CenterLabels();
316 m_PSBypassMode->SetMinimum(0);
317
318 m_unpackErr = new TProfile("unpackErr", "Unpacker errors (per boardstack); slot number; fraction of events", 64, 1, 17, 0, 2);
319 m_unpackErr->GetXaxis()->SetNdivisions(16);
320 m_unpackErr->GetXaxis()->CenterLabels();
321 m_unpackErr->SetMinimum(0);
322 m_unpackErr->SetMarkerStyle(24);
323 m_unpackErr->SetFillColor(2);
324 m_unpackErr->SetDrawOption("hist");
325
326 // cd back to root directory
327 oldDir->cd();
328 }
329
331 {
332 // Register histograms (calls back defineHisto)
333
334 REG_HISTOGRAM;
335
336 // register dataobjects
337
338 m_rawFTSWs.isOptional();
339 m_productionEventDebugs.isOptional(); // not in sim
340 m_unpackerErrors.isOptional(); // not in sim
341 m_digits.isRequired();
342 m_recBunch.isOptional();
343 m_timeZeros.isOptional();
344 m_tracks.isOptional();
345 }
346
348 {
349 m_BoolEvtMonitor->Reset();
350 m_window_vs_slot->Reset();
351 m_eventT0->Reset();
352 m_bunchOffset->Reset();
353 m_time->Reset();
354 m_timeBG->Reset();
355 m_signalHits->Reset();
356 m_backgroundHits->Reset();
357 m_trackHits->Reset();
358 m_goodTDCAll->Reset();
359 m_badTDCAll->Reset();
360 m_goodHitsPerEventAll->Reset();
361 m_badHitsPerEventAll->Reset();
362 m_TOPOccAfterInjLER->Reset();
363 m_TOPOccAfterInjHER->Reset();
364 m_TOPEOccAfterInjLER->Reset();
365 m_TOPEOccAfterInjHER->Reset();
366 m_nhitInjLER->Reset();
367 m_nhitInjHER->Reset();
368 m_nhitInjLERcut->Reset();
369 m_nhitInjHERcut->Reset();
370 m_eventInjLER->Reset();
371 m_eventInjHER->Reset();
372 m_eventInjLERcut->Reset();
373 m_eventInjHERcut->Reset();
374
375 for (int i = 0; i < m_numModules; i++) {
376 m_window_vs_asic[i]->Reset();
377 m_goodHitsXY[i]->Reset();
378 m_badHitsXY[i]->Reset();
379 m_goodHitsAsics[i]->Reset();
380 m_badHitsAsics[i]->Reset();
381 m_goodTDC[i]->Reset();
382 m_badTDC[i]->Reset();
383 m_goodTiming[i]->Reset();
384 m_goodTimingBG[i]->Reset();
385 m_goodChannelHits[i]->Reset();
386 m_badChannelHits[i]->Reset();
387 m_pulseHeights[i]->Reset();
388 }
389
390 m_skipProcFlag->Reset();
391 m_injVetoFlag->Reset();
392 m_injVetoFlagDiff->Reset();
393 m_PSBypassMode->Reset();
394 m_unpackErr->Reset();
395 }
396
398 {
399
400 // check if event time is reconstructed; distinguish collision data and cosmics
401
402 bool recBunchValid = false;
403 bool cosmics = false;
404 if (m_recBunch.isValid()) { // collision data
405 recBunchValid = m_recBunch->isReconstructed(); // event time is reconstructed
406 } else if (m_timeZeros.getEntries() == 1) { // cosmics w/ reconstructed event time
407 cosmics = true;
408 m_eventT0->Fill(m_timeZeros[0]->getTime());
409 }
410
411 // fill bunch offset
412
413 if (recBunchValid) {
414 double t0 = m_commonT0->isRoughlyCalibrated() ? m_commonT0->getT0() : 0;
415 double offset = m_recBunch->getCurrentOffset() - t0;
416 offset -= m_bunchTimeSep * lround(offset / m_bunchTimeSep); // wrap around
417 m_bunchOffset->Fill(offset);
418 m_eventT0->Fill(m_recBunch->getTime());
419 }
420
421 // fill event desynchronization
422
423 if (m_digits.getEntries() > 0) {
424 for (const auto& digit : m_digits) {
425 int x = digit.getFirstWindow() != m_digits[0]->getFirstWindow() ? 1 : 0 ;
426 m_BoolEvtMonitor->Fill(x);
427 }
428 }
429
430 // count tracks in the modules and store the momenta
431
432 std::vector<int> numTracks(m_numModules, 0);
433 std::vector<double> trackMomenta(m_numModules, 0);
434 for (const auto& track : m_tracks) {
435 const auto* fitResult = track.getTrackFitResultWithClosestMass(Const::pion);
436 if (not fitResult) continue;
437 int slot = getModuleID(track);
438 if (slot == 0) continue;
439 numTracks[slot - 1]++;
440 trackMomenta[slot - 1] = std::max(trackMomenta[slot - 1], fitResult->getMomentum().R());
441 }
442
443 // count events w/ and w/o track in the slot
444
445 if (recBunchValid or cosmics) {
446 for (size_t i = 0; i < numTracks.size(); i++) {
447 bool hit = numTracks[i] > 0;
448 m_trackHits->Fill(i + 1, hit);
449 }
450 }
451
452 // select modules for counting hits in signal and background time windows
453
454 std::vector<bool> selectedSlots(m_numModules, false);
455 for (size_t i = 0; i < selectedSlots.size(); i++) {
456 selectedSlots[i] = (recBunchValid or cosmics) and numTracks[i] == 1 and trackMomenta[i] > m_momentumCut;
457 }
458
459 // prepare counters
460
461 int nHits_good = 0;
462 int nHits_bad = 0;
463 std::vector<int> numSignalHits(m_numModules, 0);
464 std::vector<int> numBackgroundHits(m_numModules, 0);
465
466 // loop over digits, fill histograms and increment counters
467
468 for (const auto& digit : m_digits) {
469 int slot = digit.getModuleID();
470 if (slot < 1 or slot > m_numModules) {
471 B2ERROR("Invalid slot ID found in TOPDigits: ID = " << slot);
472 continue;
473 }
474 int asic_no = digit.getChannel() / 8;
475 int asic_ch = digit.getChannel() % 8;
476
477 m_window_vs_slot->Fill(digit.getModuleID(), digit.getRawTime() / 64 + 220);
478 m_window_vs_asic[slot - 1]->Fill(asic_no, digit.getRawTime() / 64 + 220);
479
480 if (digit.getHitQuality() != TOPDigit::c_Junk) { // good hits
481 m_goodHitsXY[slot - 1]->Fill(digit.getPixelCol(), digit.getPixelRow());
482 m_goodHitsAsics[slot - 1]->Fill(asic_no, asic_ch);
483 m_goodTDC[slot - 1]->Fill(digit.getRawTime());
484 m_goodTDCAll->Fill(digit.getRawTime());
485 m_goodChannelHits[slot - 1]->Fill(digit.getChannel());
486 m_pulseHeights[slot - 1]->Fill(digit.getPMTNumber(), digit.getPulseHeight());
487 nHits_good++;
488 if (digit.hasStatus(TOPDigit::c_EventT0Subtracted)) {
489 double time = digit.getTime();
490 if (cosmics) time += 10; // move for 10 ns in order to get the complete signal at positive times
491 if (numTracks[slot - 1] > 0) {
492 m_goodTiming[slot - 1]->Fill(time);
493 m_time->Fill(time);
494 } else {
495 m_goodTimingBG[slot - 1]->Fill(time);
496 m_timeBG->Fill(time);
497 }
498 if (selectedSlots[slot - 1] and abs(time) < 50) {
499 if (time > 0) numSignalHits[slot - 1]++;
500 else numBackgroundHits[slot - 1]++;
501 }
502 }
503 } else { // bad hits: FE not valid, pedestal jump, too short or too wide pulses
504 m_badHitsXY[slot - 1]->Fill(digit.getPixelCol(), digit.getPixelRow());
505 m_badHitsAsics[slot - 1]->Fill(asic_no, asic_ch);
506 m_badTDC[slot - 1]->Fill(digit.getRawTime());
507 m_badTDCAll->Fill(digit.getRawTime());
508 m_badChannelHits[slot - 1]->Fill(digit.getChannel());
509 nHits_bad++;
510 }
511 }
512
513 // histogram counters
514
515 m_goodHitsPerEventAll->Fill(nHits_good);
516 m_badHitsPerEventAll->Fill(nHits_bad);
517 for (int slot = 1; slot <= m_numModules; slot++) {
518 if (selectedSlots[slot - 1]) {
519 m_signalHits->Fill(slot, numSignalHits[slot - 1]);
520 m_backgroundHits->Fill(slot, numBackgroundHits[slot - 1]);
521 }
522 }
523
524 // fill injection histograms
525
526 for (auto& it : m_rawFTSWs) {
527 B2DEBUG(29, "TTD FTSW : " << hex << it.GetTTUtime(0) << " " << it.GetTTCtime(0) << " EvtNr " << it.GetEveNo(0) << " Type " <<
528 (it.GetTTCtimeTRGType(0) & 0xF) << " TimeSincePrev " << it.GetTimeSincePrevTrigger(0) << " TimeSinceInj " <<
529 it.GetTimeSinceLastInjection(0) << " IsHER " << it.GetIsHER(0) << " Bunch " << it.GetBunchNumber(0));
530 auto time_clk = it.GetTimeSinceLastInjection(0); // expressed in system clock
531 if (time_clk != 0x7FFFFFFF) {
532 unsigned int nentries = m_digits.getEntries();
533 double time_us = time_clk / 127.0; // 127MHz clock ticks to us, inexact rounding
534 double time_ms = time_us / 1000;
535 double time_1280 = time_clk % 1280; // within beam cycle, expressed in system clock
536 if (it.GetIsHER(0)) {
537 m_TOPOccAfterInjHER->Fill(time_us, nentries);
538 m_TOPEOccAfterInjHER->Fill(time_us);
539 m_nhitInjHER->Fill(time_ms, time_1280, nHits_good);
540 m_eventInjHER->Fill(time_ms, time_1280);
541 if (nHits_good > 1000) {
542 m_nhitInjHERcut->Fill(time_ms, time_1280, nHits_good);
543 m_eventInjHERcut->Fill(time_ms, time_1280);
544 }
545 } else {
546 m_TOPOccAfterInjLER->Fill(time_us, nentries);
547 m_TOPEOccAfterInjLER->Fill(time_us);
548 m_nhitInjLER->Fill(time_ms, time_1280, nHits_good);
549 m_eventInjLER->Fill(time_ms, time_1280);
550 if (nHits_good > 1000) {
551 m_nhitInjLERcut->Fill(time_ms, time_1280, nHits_good);
552 m_eventInjLERcut->Fill(time_ms, time_1280);
553 }
554 }
555 }
556 }
557
558 // fill raw data header histograms
559
560 const auto& feMapper = TOPGeometryPar::Instance()->getFrontEndMapper();
561 double differ = 0;
562 for (const auto& dbg : m_productionEventDebugs) {
563 auto scrodID = dbg.getScrodID();
564 const auto* femap = feMapper.getMap(scrodID);
565 if (not femap) {
566 B2ERROR("No front-end map available for scrodID " << scrodID);
567 continue;
568 }
569 auto slot = femap->getModuleID();
570 auto bs = femap->getBoardstackNumber();
571 double x = slot + bs / 4.0;
572 m_skipProcFlag->Fill(x, dbg.getSkipProcessingFlag());
573 m_injVetoFlag->Fill(x, dbg.getInjectionVetoFlag());
574 m_PSBypassMode->Fill(x, dbg.getPSBypassMode());
575 if (dbg.getInjectionVetoFlag() != m_productionEventDebugs[0]->getInjectionVetoFlag()) differ = 1;
576 }
577 m_injVetoFlagDiff->Fill(differ);
578
579 if (m_unpackerErrors.isValid()) {
580 const auto& flags = m_unpackerErrors->getErrorFlags();
581 for (size_t bit = 0; bit < flags.size(); bit++) {
582 int slot = bit / 4 + 1;
583 int bs = bit % 4;
584 double x = slot + bs / 4.0;
585 m_unpackErr->Fill(x, flags[bit]);
586 }
587 }
588
589 }
590
591
592 int TOPDQMModule::getModuleID(const Track& track) const
593 {
594 Const::EDetector myDetID = Const::EDetector::TOP;
595 int pdgCode = std::abs(Const::pion.getPDGCode());
596
597 RelationVector<ExtHit> extHits = track.getRelationsWith<ExtHit>();
598 for (const auto& extHit : extHits) {
599 if (std::abs(extHit.getPdgCode()) != pdgCode) continue;
600 if (extHit.getDetectorID() != myDetID) continue;
601 if (extHit.getCopyID() < 1 or extHit.getCopyID() > m_numModules) continue;
602 return extHit.getCopyID();
603 }
604
605 return 0;
606 }
607
608
610} // end Belle2 namespace
611
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
TProfile * m_unpackErr
fraction of events with boardstack unpacker errors
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.
StoreObjPtr< TOPUnpackerErrors > m_unpackerErrors
unpacker error flags
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.