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