Belle II Software  release-05-02-19
TOPChannelMaskCollectorModule.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2021 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Marko Staric *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #include <top/modules/collectors/TOPChannelMaskCollectorModule.h>
12 #include <TH1F.h>
13 #include <TH2F.h>
14 
15 using namespace std;
16 
17 namespace Belle2 {
23  //-----------------------------------------------------------------
24  // Register module
25  //-----------------------------------------------------------------
26 
27  REG_MODULE(TOPChannelMaskCollector)
28 
29  //-----------------------------------------------------------------
30  // Implementation
31  //-----------------------------------------------------------------
32 
34  {
35  // set module description and processing properties
36  setDescription("A collector for preparing masks of dead and hot channels");
37  setPropertyFlags(c_ParallelProcessingCertified);
38  }
39 
40 
41  void TOPChannelMaskCollectorModule::prepare()
42  {
43 
44  m_digits.isRequired();
45 
46  auto nhits = new TH1F("nhits", "Number of good hits per event; hits per event; entries per bin", 200, 0, 2000);
47  registerObject<TH1F>("nhits", nhits);
48 
49  for (int slot = 1; slot <= 16; slot++) {
50  string name = "hits_" + to_string(slot);
51  string title = "Channel occupancies for slot " + to_string(slot);
52  auto h = new TH1F(name.c_str(), title.c_str(), 512, 0, 512);
53  h->SetXTitle("channel number");
54  h->SetYTitle("number of hits per channel");
55  registerObject<TH1F>(name, h);
56  m_names.push_back(name);
57  }
58 
59  for (int slot = 1; slot <= 16; slot++) {
60  string name = "window_vs_asic_" + to_string(slot);
61  string title = "Window vs. asic for slot " + to_string(slot);
62  auto h = new TH2F(name.c_str(), title.c_str(), 64, 0, 64, 512, 0, 512);
63  h->SetXTitle("ASIC number");
64  h->SetYTitle("window number w.r.t reference window");
65  registerObject<TH2F>(name, h);
66  m_asicNames.push_back(name);
67  }
68 
69  }
70 
71 
72  void TOPChannelMaskCollectorModule::collect()
73  {
74 
75  std::vector<TH1F*> slots;
76  for (const auto& name : m_names) {
77  auto h = getObjectPtr<TH1F>(name);
78  slots.push_back(h);
79  }
80 
81  std::vector<TH2F*> asics;
82  for (const auto& name : m_asicNames) {
83  auto h = getObjectPtr<TH2F>(name);
84  asics.push_back(h);
85  }
86 
87  int NGood = 0;
88  for (const auto& digit : m_digits) {
89  if (digit.getHitQuality() == TOPDigit::c_Junk) continue;
90  NGood++;
91  unsigned m = digit.getModuleID() - 1;
92  if (m < slots.size()) slots[m]->Fill(digit.getChannel());
93  if (m < asics.size()) asics[m]->Fill(digit.getChannel() / 8, digit.getRawTime() / 64 + 220);
94  }
95 
96  auto h = getObjectPtr<TH1F>("nhits");
97  h->Fill(NGood);
98 
99  }
100 
102 }
103 
Belle2::TOPChannelMaskCollectorModule
Collector for preparing masks of hot and dead channels.
Definition: TOPChannelMaskCollectorModule.h:39
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19