Belle II Software  release-05-02-19
TOPPulseHeightCollectorModule.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2019 - 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/TOPPulseHeightCollectorModule.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(TOPPulseHeightCollector)
28 
29  //-----------------------------------------------------------------
30  // Implementation
31  //-----------------------------------------------------------------
32 
34  {
35  // set module description and processing properties
36  setDescription("A collector for channel pulse-height distributions");
37  setPropertyFlags(c_ParallelProcessingCertified);
38 
39  // module parameters
40  addParam("nx", m_nx, "number of histogram bins", 200);
41  addParam("xmax", m_xmax, "histogram upper limit [ADC counts]", 2000.0);
42  addParam("pulseWidthWindow", m_widthWindow,
43  "lower and upper bound of pulse-width selection window [ns]. "
44  "Empty list means no selection on the pulse width. "
45  "Note: selection on pulse width will influence pulse-height distribution.",
46  m_widthWindow);
47  addParam("timeWindow", m_timeWindow,
48  "lower and upper bound of time selection window [ns]. "
49  "Empty list means no selection on photon time.", m_timeWindow);
50 
51  }
52 
53 
54  void TOPPulseHeightCollectorModule::prepare()
55  {
56 
57  m_digits.isRequired();
58 
59  auto h1a = new TH1F("time", "time distribution (all hits)", 1000, -100, 250);
60  h1a->SetXTitle("time [ns]");
61  registerObject<TH1F>("time", h1a);
62 
63  auto h1b = new TH1F("time_sel", "time distribution (selected hits)", 1000, -100, 250);
64  h1b->SetXTitle("time [ns]");
65  registerObject<TH1F>("time_sel", h1b);
66 
67  auto h2a = new TH2F("ph_vs_width", "pulse height vs. width (all hits)",
68  200, 0, 10, 200, 0, 2000);
69  h2a->SetXTitle("pulse width [ns]");
70  h2a->SetYTitle("pulse height [ADC counts]");
71  registerObject<TH2F>("ph_vs_width", h2a);
72 
73  auto h2b = new TH2F("ph_vs_width_sel", "pulse height vs. width (selected hits)",
74  200, 0, 10, 200, 0, 2000);
75  h2b->SetXTitle("pulse width [ns]");
76  h2b->SetYTitle("pulse height [ADC counts]");
77  registerObject<TH2F>("ph_vs_width_sel", h2b);
78 
79  for (int slot = 1; slot <= 16; slot++) {
80  string name = "ph_slot_" + to_string(slot);
81  string title = "pulse-height vs. channel for slot " + to_string(slot);
82  auto h = new TH2F(name.c_str(), title.c_str(), 512, 0, 512, m_nx, 0, m_xmax);
83  h->SetXTitle("channel number");
84  h->SetYTitle("pulse height [ADC counts]");
85  registerObject<TH2F>(name, h);
86  m_names.push_back(name);
87  }
88 
89  }
90 
91 
92  void TOPPulseHeightCollectorModule::collect()
93  {
94 
95  std::vector<TH2F*> slots;
96  for (const auto& name : m_names) {
97  auto h = getObjectPtr<TH2F>(name);
98  slots.push_back(h);
99  }
100 
101  auto h1a = getObjectPtr<TH1F>("time");
102  auto h1b = getObjectPtr<TH1F>("time_sel");
103  auto h2a = getObjectPtr<TH2F>("ph_vs_width");
104  auto h2b = getObjectPtr<TH2F>("ph_vs_width_sel");
105 
106  for (const auto& digit : m_digits) {
107  if (digit.getHitQuality() != TOPDigit::c_Good) continue;
108  auto t = digit.getTime();
109  auto w = digit.getPulseWidth();
110  auto ph = digit.getPulseHeight();
111  h1a->Fill(t);
112  h2a->Fill(w, ph);
113  if (m_widthWindow.size() == 2) {
114  if (w < m_widthWindow[0] or w > m_widthWindow[1]) continue;
115  }
116  if (m_timeWindow.size() == 2) {
117  if (t < m_timeWindow[0] or t > m_timeWindow[1]) continue;
118  }
119  h1b->Fill(t);
120  h2b->Fill(w, ph);
121  unsigned m = digit.getModuleID() - 1;
122  if (m < slots.size()) slots[m]->Fill(digit.getChannel(), ph);
123  }
124 
125  }
126 
128 }
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::TOPPulseHeightCollectorModule
Collector for channel pulse-height distributions.
Definition: TOPPulseHeightCollectorModule.h:39
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19