Belle II Software development
TOPPhotonYieldsCollectorModule.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/collectors/TOPPhotonYieldsCollectorModule.h>
10#include <top/reconstruction_cpp/TOPTrack.h>
11#include <top/geometry/TOPGeometryPar.h>
12#include <top/dataobjects/TOPTimeZero.h>
13
14// framework aux
15#include <framework/gearbox/Unit.h>
16#include <framework/logging/Logger.h>
17
18// root
19#include <TH1F.h>
20#include <TH2F.h>
21#include <TProfile.h>
22
23using namespace std;
24
25namespace Belle2 {
30
31 using namespace TOP;
32
33 //-----------------------------------------------------------------
35 //-----------------------------------------------------------------
36
37 REG_MODULE(TOPPhotonYieldsCollector);
38
39 //-----------------------------------------------------------------
40 // Implementation
41 //-----------------------------------------------------------------
42
44 {
45 // set module description and processing properties
46 setDescription("A collector for photon pixel yields aimed for PMT ageing studies and for finding optically decoupled PMT's");
48
49 // module parameters
50 addParam("sample", m_sample, "sample type: one of dimuon or bhabha", std::string("dimuon"));
51 addParam("deltaEcms", m_deltaEcms, "c.m.s energy window (half size) if sample is dimuon or bhabha", 0.1);
52 addParam("dr", m_dr, "cut on POCA in r", 2.0);
53 addParam("dz", m_dz, "cut on POCA in abs(z)", 4.0);
54
55 }
56
57
59 {
60 // input collections
61 m_digits.isRequired();
62 m_tracks.isRequired();
63 m_extHits.isRequired();
64 m_recBunch.isRequired();
65 m_asicMask.isRequired();
66 m_associatedPDFs.isRequired();
67
68 // set track selector
69 if (m_sample == "dimuon" or m_sample == "bhabha") {
71 m_selector.setDeltaEcms(m_deltaEcms);
72 m_selector.setCutOnPOCA(m_dr, m_dz);
73 m_selector.setCutOnLocalZ(m_minZprism, m_maxZ);
74 } else {
75 B2ERROR("Invalid sample type '" << m_sample << "'");
76 }
77
78 // create and register histograms
79
80 const int numModules = 16;
81 const int numPixels = 512;
82
83 // TOF corrections
84 auto* tofCorrections = new TProfile("tofCorrections", "TOF corrections; local z [cm]; t0 [ns]", 1000, m_minZprism, m_maxZ, -1, 1);
85 registerObject<TProfile>("tofCorrections", tofCorrections);
86
87 // time stamp (average unix time and its standard deviation)
88 auto* timeStamp = new TProfile("timeStamp", "Time stamp; ; unix time", 1, 0, 1, 0, 1.0e10, "S");
89 registerObject<TProfile>("timeStamp", timeStamp);
90
91 // number of selected tracks per slot
92 auto* numTracks = new TH1F("numTracks", "Number of tracks per slot; slot number; track count", numModules, 0.5, numModules + 0.5);
93 registerObject<TH1F>("numTracks", numTracks);
94
95 // number of pixel hits in a signal time window
96 for (int slot = 1; slot <= numModules; slot++) {
97 string name = (slot < 10) ? "signalHits_0" + to_string(slot) : "signalHits_" + to_string(slot);
98 string title = "Hits in signal window for slot " + to_string(slot);
99 auto h = new TH1F(name.c_str(), title.c_str(), numPixels, 0.5, numPixels + 0.5);
100 h->SetXTitle("pixel number");
101 h->SetYTitle("hit count");
102 registerObject<TH1F>(name, h);
103 m_signalNames.push_back(name);
104 }
105
106 // number of pixel hits in a background time window
107 for (int slot = 1; slot <= numModules; slot++) {
108 string name = (slot < 10) ? "bkgHits_0" + to_string(slot) : "bkgHits_" + to_string(slot);
109 string title = "Hits in background window for slot " + to_string(slot);
110 auto h = new TH1F(name.c_str(), title.c_str(), numPixels, 0.5, numPixels + 0.5);
111 h->SetXTitle("pixel number");
112 h->SetYTitle("hit count");
113 registerObject<TH1F>(name, h);
114 m_bkgNames.push_back(name);
115 }
116
117 // active pixels
118 for (int slot = 1; slot <= numModules; slot++) {
119 string name = (slot < 10) ? "activePixels_0" + to_string(slot) : "activePixels_" + to_string(slot);
120 string title = "Active pixels for slot " + to_string(slot);
121 auto h = new TH1F(name.c_str(), title.c_str(), numPixels, 0.5, numPixels + 0.5);
122 h->SetXTitle("pixel number");
123 h->SetYTitle("track count");
124 registerObject<TH1F>(name, h);
125 m_activeNames.push_back(name);
126 }
127
128 // number of effective signal hits in pixels
129 for (int slot = 1; slot <= numModules; slot++) {
130 string name = (slot < 10) ? "effectiveSignalHits_0" + to_string(slot) : "effectiveSignalHits_" + to_string(slot);
131 string title = "Effective signal hits for slot " + to_string(slot);
132 auto h = new TH1F(name.c_str(), title.c_str(), numPixels, 0.5, numPixels + 0.5);
133 h->SetXTitle("pixel number");
134 h->SetYTitle("hit count");
135 registerObject<TH1F>(name, h);
136 m_effectiveSignalNames.push_back(name);
137 }
138
139 // number of pixel hits with low impact angle on photo cathode
140 for (int slot = 1; slot <= numModules; slot++) {
141 string name = (slot < 10) ? "alphaLow_0" + to_string(slot) : "alphaLow_" + to_string(slot);
142 string title = "Hits w/ low alpha for slot " + to_string(slot);
143 auto h = new TH1F(name.c_str(), title.c_str(), numPixels, 0.5, numPixels + 0.5);
144 h->SetXTitle("pixel number");
145 h->SetYTitle("hit count");
146 registerObject<TH1F>(name, h);
147 m_alphaLowNames.push_back(name);
148 }
149
150 // number of pixel hits with high impact angle on photo cathode
151 for (int slot = 1; slot <= numModules; slot++) {
152 string name = (slot < 10) ? "alphaHigh_0" + to_string(slot) : "alphaHigh_" + to_string(slot);
153 string title = "Hits w/ high alpha for slot " + to_string(slot);
154 auto h = new TH1F(name.c_str(), title.c_str(), numPixels, 0.5, numPixels + 0.5);
155 h->SetXTitle("pixel number");
156 h->SetYTitle("hit count");
157 registerObject<TH1F>(name, h);
158 m_alphaHighNames.push_back(name);
159 }
160
161 // pixel pulse-height distributions
162 for (int slot = 1; slot <= numModules; slot++) {
163 string name = (slot < 10) ? "pulseHeights_0" + to_string(slot) : "pulseHeights_" + to_string(slot);
164 string title = "Pulse height distributions for slot " + to_string(slot);
165 auto h = new TH2F(name.c_str(), title.c_str(), numPixels, 0.5, numPixels + 0.5, 200, 0, 2000);
166 h->SetXTitle("pixel number");
167 h->SetYTitle("pulse height");
168 registerObject<TH2F>(name, h);
169 m_pulseHeightNames.push_back(name);
170 }
171
172 // local z-distribution of tracks
173 for (int slot = 1; slot <= numModules; slot++) {
174 string name = (slot < 10) ? "muonZ_0" + to_string(slot) : "muonZ_" + to_string(slot);
175 string title = "Track z-distribution for slot " + to_string(slot);
176 auto h = new TH1F(name.c_str(), title.c_str(), 100, m_minZ, m_maxZ);
177 h->SetXTitle("local z [cm]");
178 h->SetYTitle("track count");
179 registerObject<TH1F>(name, h);
180 m_muonZNames.push_back(name);
181 }
182
183 }
184
185
187 {
188 // bunch must be reconstructed
189
190 if (not m_recBunch->isReconstructed()) return;
191
192 // loop over reconstructed tracks, make a selection and fill histograms
193
194 for (const auto& track : m_tracks) {
195 // track selection
196 TOPTrack trk(track);
197 if (not trk.isValid()) continue;
198 if (not m_selector.isSelected(trk)) continue;
199
200 // fill TOF corrections
201 double localZ = m_selector.getLocalPosition().Z();
202 const auto* timeZero = trk.getExtHit()->getRelated<TOPTimeZero>();
203 if (timeZero and timeZero->isValid()) {
204 auto tofCorrections = getObjectPtr<TProfile>("tofCorrections");
205 tofCorrections->Fill(localZ, timeZero->getTime());
206 }
207 if (localZ < m_minZ) continue;
208
209 // fill all the other histograms
210 auto timeStamp = getObjectPtr<TProfile>("timeStamp");
211 timeStamp->Fill(0.5, m_eventMetaData->getTime() / 1000000000);
212
213 int slot = trk.getModuleID();
214 auto numTracks = getObjectPtr<TH1F>("numTracks");
215 numTracks->Fill(slot);
216
217 auto muonZ = getObjectPtr<TH1F>(m_muonZNames[slot - 1]);
218 muonZ->Fill(m_selector.getLocalPosition().Z());
219
220 auto signalHits = getObjectPtr<TH1F>(m_signalNames[slot - 1]);
221 auto bkgHits = getObjectPtr<TH1F>(m_bkgNames[slot - 1]);
222 auto pulseHeight = getObjectPtr<TH2F>(m_pulseHeightNames[slot - 1]);
223 for (const auto& digit : m_digits) {
224 if (digit.getModuleID() != slot) continue;
225 if (digit.getHitQuality() != TOPDigit::c_Good) continue; // junk hit or pixel masked-out
226 if (std::abs(digit.getTime()) > m_timeWindow) continue;
227 if (digit.getTime() > 0) {
228 signalHits->Fill(digit.getPixelID());
229 pulseHeight->Fill(digit.getPixelID(), digit.getPulseHeight());
230 } else {
231 bkgHits->Fill(digit.getPixelID());
232 }
233 }
234
235 auto activePixels = getObjectPtr<TH1F>(m_activeNames[slot - 1]);
236 const auto& chMapper = TOPGeometryPar::Instance()->getChannelMapper();
237 for (int pixel = 1; pixel <= activePixels->GetNbinsX(); pixel++) {
238 unsigned channel = chMapper.getChannel(pixel);
239 if (m_channelMask->isActive(slot, channel) and m_asicMask->isActive(slot, channel)) activePixels->Fill(pixel);
240 }
241
242 auto alphaLow = getObjectPtr<TH1F>(m_alphaLowNames[slot - 1]);
243 auto alphaHigh = getObjectPtr<TH1F>(m_alphaHighNames[slot - 1]);
244 auto effectiveSignalHits = getObjectPtr<TH1F>(m_effectiveSignalNames[slot - 1]);
245 for (const auto& digit : m_digits) {
246 if (digit.getModuleID() != slot) continue;
247 if (digit.getHitQuality() != TOPDigit::c_Good) continue; // junk hit or pixel masked-out
248 const auto* pdf = digit.getRelated<TOPAssociatedPDF>();
249 if (not pdf) continue;
250 const auto* peak = pdf->getSinglePeak();
251 if (not peak) continue; // hit associated with background
252 effectiveSignalHits->Fill(digit.getPixelID());
253 if (std::abs(m_selector.getLocalPosition().Z()) > m_excludedZ) {
254 double alpha = acos(std::abs(peak->kzd)) / Unit::deg;
255 if (alpha > 60) continue;
256 if (alpha < 30) alphaLow->Fill(digit.getPixelID());
257 else alphaHigh->Fill(digit.getPixelID());
258 }
259 }
260 } // loop over tracks
261
262 }
263
264
266} // end namespace Belle2
void registerObject(const std::string &name, T *obj)
Register object with a name, takes ownership, do not access the pointer beyond prepare()
T * getObjectPtr(const std::string &name)
Calls the CalibObjManager to get the requested stored collector data.
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
T * getRelated(const std::string &name="", const std::string &namedRelation="") const
Get the object to or from which this object has a relation.
Class to store analytic PDF associated with a photon.
DBObjPtr< TOPCalChannelMask > m_channelMask
masked channels
StoreObjPtr< TOPRecBunch > m_recBunch
reconstructed bunch
TOP::TrackSelector m_selector
track selection utility
std::vector< std::string > m_activeNames
histogram names for active pixels count
std::vector< std::string > m_alphaHighNames
histogram names for counting hits w/ high impact angle on photo cathode
std::vector< std::string > m_pulseHeightNames
histogram names for pulse heights
StoreObjPtr< EventMetaData > m_eventMetaData
event meta data object
std::vector< std::string > m_effectiveSignalNames
histogram names for effective signal hits in pixels
std::vector< std::string > m_muonZNames
histogram names for track z-distribution
std::vector< std::string > m_alphaLowNames
histogram names for counting hits w/ low impact angle on photo cathode
const double m_timeWindow
time window for counting photon hits (half size)
StoreArray< Track > m_tracks
collection of tracks
StoreArray< TOPAssociatedPDF > m_associatedPDFs
collection of PDF's associated to TOP digits
const double m_minZ
minimal local z of extrapolated track
const double m_excludedZ
excluded central region of extrapolated track for photon impact angle counting
const double m_maxZ
maximal local z of extrapolated track
std::vector< std::string > m_signalNames
histogram names for signal window hit counts
std::vector< std::string > m_bkgNames
histogram names for background window hit counts
const double m_minZprism
most minimal local z of extrapolated track (this includes prism)
StoreArray< TOPDigit > m_digits
collection of TOP digits
StoreArray< ExtHit > m_extHits
collection of extrapolated hits
StoreObjPtr< TOPAsicMask > m_asicMask
online masked Asics
Class to store T0 information.
Definition TOPTimeZero.h:24
const ChannelMapper & getChannelMapper() const
Returns default channel mapper (mapping of channels to pixels)
static TOPGeometryPar * Instance()
Static method to obtain the pointer to its instance.
Reconstructed track at TOP.
Definition TOPTrack.h:40
bool isValid() const
Checks if track is successfully constructed.
Definition TOPTrack.h:144
const ExtHit * getExtHit() const
Returns extrapolated hit (track entrance to the bar)
Definition TOPTrack.h:223
int getModuleID() const
Returns slot ID.
Definition TOPTrack.h:150
Utility for the track selection - used in various calibration modules.
static const double deg
degree to radians
Definition Unit.h:109
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
virtual void collect() final
Replacement for event().
virtual void prepare() final
Replacement for initialize().
Abstract base class for different kinds of events.
STL namespace.