Belle II Software development
TOPChannelT0CalibratorModule.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/TOPChannelT0Calibrator/TOPChannelT0CalibratorModule.h>
10#include <top/geometry/TOPGeometryPar.h>
11#include <top/reconstruction_cpp/TOPTrack.h>
12#include <top/reconstruction_cpp/PDFConstructor.h>
13#include <top/reconstruction_cpp/TOPRecoManager.h>
14#include <framework/dataobjects/EventMetaData.h>
15#include <framework/gearbox/Const.h>
16#include <framework/logging/Logger.h>
17#include <TRandom.h>
18
19using namespace std;
20
21namespace Belle2 {
26 using namespace TOP;
27
28 //-----------------------------------------------------------------
30 //-----------------------------------------------------------------
31
32 REG_MODULE(TOPChannelT0Calibrator);
33
34 //-----------------------------------------------------------------
35 // Implementation
36 //-----------------------------------------------------------------
37
39
40 {
41 // set module description
42 setDescription("Alternative channel T0 calibration with dimuons or bhabhas. "
43 "This module can also be used to check the calibration "
44 "(channel, module and common T0).\n"
45 "Note: after this kind of calibration one cannot do the alignment.");
46 // setPropertyFlags(c_ParallelProcessingCertified);
47
48 // Add parameters
49 addParam("numBins", m_numBins, "number of bins of the search region", 200);
50 addParam("timeRange", m_timeRange,
51 "time range in which to search for the minimum [ns]", 10.0);
52 addParam("sigmaSmear", m_sigmaSmear,
53 "sigma in [ns] for additional smearing of PDF", 0.0);
54 addParam("sample", m_sample,
55 "sample type: one of dimuon, bhabha or cosmics", std::string("dimuon"));
56 addParam("minMomentum", m_minMomentum,
57 "minimal track momentum if sample is cosmics", 1.0);
58 addParam("deltaEcms", m_deltaEcms,
59 "c.m.s energy window (half size) if sample is dimuon or bhabha", 0.1);
60 addParam("dr", m_dr, "cut on POCA in r", 2.0);
61 addParam("dz", m_dz, "cut on POCA in abs(z)", 4.0);
62 addParam("minZ", m_minZ,
63 "minimal local z of extrapolated hit", -130.0);
64 addParam("maxZ", m_maxZ,
65 "maximal local z of extrapolated hit", 130.0);
66 addParam("outputFileName", m_outFileName,
67 "Root output file name containing calibration results. "
68 "File name can include *'s; "
69 "they will be replaced with a run number from the first input file",
70 std::string("calibrationT0_r*.root"));
71 addParam("pdfOption", m_pdfOption,
72 "PDF option, one of 'rough', 'fine', 'optimal'", std::string("rough"));
73
74 }
75
77 {
78 // input collections
79 m_digits.isRequired();
80 m_tracks.isRequired();
81 m_extHits.isRequired();
82 m_recBunch.isOptional();
83
84 // Parse PDF option
85 if (m_pdfOption == "rough") {
87 } else if (m_pdfOption == "fine") {
89 } else if (m_pdfOption == "optimal") {
91 } else {
92 B2ERROR("TOPPDFDebuggerModule: unknown PDF option '" << m_pdfOption << "'");
93 }
94
95 // set track selector
101
102 // minimum finders
103 double tmin = -m_timeRange / 2;
104 double tmax = m_timeRange / 2;
105 for (unsigned i = 0; i < 2; i++) {
106 for (unsigned m = 0; m < c_numModules; m++) {
107 for (unsigned c = 0; c < c_numChannels; c++) {
108 m_finders[i][m][c] = Chi2MinimumFinder1D(m_numBins, tmin, tmax);
109 }
110 }
111 }
112
113 // if file name includes *'s replace them with a run number
114 auto pos = m_outFileName.find("*");
115 if (pos != std::string::npos) {
116 StoreObjPtr<EventMetaData> evtMetaData;
117 auto run = std::to_string(evtMetaData->getRun());
118 while (run.size() < 5) run = "0" + run;
119 while (pos != std::string::npos) {
120 m_outFileName.replace(pos, 1, run);
121 pos = m_outFileName.find("*");
122 }
123 }
124
125 // open root file for ntuple and histogram output
126 m_file = TFile::Open(m_outFileName.c_str(), "RECREATE");
127 if (not m_file) {
128 B2ERROR("Cannot open output file '" << m_outFileName << "'");
129 return;
130 }
131
132 // histograms
133 for (unsigned module = 0; module < c_numModules; module++) {
134 int moduleID = module + 1;
135
136 std::string slotNum = std::to_string(moduleID);
137 if (moduleID < 10) slotNum = "0" + slotNum;
138
139 std::string name = "numHits_slot" + slotNum;
140 std::string title = "Number of hits per channel for slot " + slotNum;
141 TH1F h1(name.c_str(), title.c_str(), c_numChannels, 0, c_numChannels);
142 h1.SetXTitle("channel number");
143 h1.SetYTitle("hits per channel");
144 m_hits1D.push_back(h1);
145
146 name = "timeHits_slot" + slotNum;
147 title = "hit time vs. channel for slot " + slotNum;
148 TH2F h2(name.c_str(), title.c_str(), c_numChannels, 0, c_numChannels,
149 200, 0.0, 20.0);
150 h2.SetXTitle("channel number");
151 h2.SetYTitle("time [ns]");
152 m_hits2D.push_back(h2);
153 }
154
155 // create output tree
156
157 m_tree = new TTree("tree", "Channel T0 calibration results");
158 m_tree->Branch("slot", &m_moduleID);
159 m_tree->Branch("numPhotons", &m_numPhotons);
160 m_tree->Branch("x", &m_x);
161 m_tree->Branch("y", &m_y);
162 m_tree->Branch("z", &m_z);
163 m_tree->Branch("p", &m_p);
164 m_tree->Branch("theta", &m_theta);
165 m_tree->Branch("phi", &m_phi);
166 m_tree->Branch("r_poca", &m_pocaR);
167 m_tree->Branch("z_poca", &m_pocaZ);
168 m_tree->Branch("x_poca", &m_pocaX);
169 m_tree->Branch("y_poca", &m_pocaY);
170 m_tree->Branch("Ecms", &m_cmsE);
171 m_tree->Branch("charge", &m_charge);
172 m_tree->Branch("PDG", &m_PDG);
173
174 }
175
177 {
178 /* check bunch reconstruction status and run alignment:
179 - if object exists and bunch is found (collision data w/ bunch finder in the path)
180 - if object doesn't exist (cosmic data and other cases w/o bunch finder)
181 */
182
183 if (m_recBunch.isValid()) {
184 if (not m_recBunch->isReconstructed()) return;
185 }
186
188 double timeMin = TOPRecoManager::getMinTime();
189 double timeMax = TOPRecoManager::getMaxTime();
190 const auto& chMapper = TOPGeometryPar::Instance()->getChannelMapper();
191
192 // loop over reconstructed tracks, make a selection and accumulate log likelihoods
193
194 for (const auto& track : m_tracks) {
195
196 // track selection
197 TOPTrack trk(track);
198 if (not trk.isValid()) continue;
199
200 if (not m_selector.isSelected(trk)) continue;
201
202 // construct PDF
204 if (not pdfConstructor.isValid()) continue;
205
206 // minimization procedure: accumulate
207 const unsigned module = trk.getModuleID() - 1;
208 if (module >= c_numModules) continue;
209 int sub = gRandom->Integer(2); // generate sub-sample number
210 auto& finders = m_finders[sub][module];
211 const auto& binCenters = finders[0].getBinCenters();
212 for (unsigned ibin = 0; ibin < binCenters.size(); ibin++) {
213 double t0 = binCenters[ibin];
214 const auto& pixelLogLs = pdfConstructor.getPixelLogLs(t0, m_sigmaSmear);
215 for (unsigned channel = 0; channel < c_numChannels; channel++) {
216 int pix = chMapper.getPixelID(channel) - 1;
217 finders[channel].add(ibin, -2 * pixelLogLs[pix].logL);
218 }
219 }
220
221 // fill histograms of hits
222 m_numPhotons = 0;
223 for (const auto& digit : m_digits) {
224 if (digit.getHitQuality() != TOPDigit::c_Good) continue;
225 if (digit.getModuleID() != trk.getModuleID()) continue;
226 if (digit.getTime() < timeMin) continue;
227 if (digit.getTime() > timeMax) continue;
228 m_numPhotons++;
229 auto& h1 = m_hits1D[module];
230 h1.Fill(digit.getChannel());
231 auto& h2 = m_hits2D[module];
232 h2.Fill(digit.getChannel(), digit.getTime());
233 }
234
235 // fill output tree
236 m_moduleID = trk.getModuleID();
237 const auto& localPosition = m_selector.getLocalPosition();
238 m_x = localPosition.X();
239 m_y = localPosition.Y();
240 m_z = localPosition.Z();
241 const auto& localMomentum = m_selector.getLocalMomentum();
242 m_p = localMomentum.R();
243 m_theta = localMomentum.Theta();
244 m_phi = localMomentum.Phi();
245 const auto& pocaPosition = m_selector.getPOCAPosition();
246 m_pocaR = pocaPosition.Rho();
247 m_pocaZ = pocaPosition.Z();
248 m_pocaX = pocaPosition.X();
249 m_pocaY = pocaPosition.Y();
251 m_charge = trk.getCharge();
252 m_PDG = trk.getPDGCode();
253 m_tree->Fill();
254 }
255
256 }
257
258
260 {
261
262 // determine scaling factor for errors from two statistically independent results
263
264 TH1F h_pulls("pulls", "Pulls of the two statistically independent results",
265 200, -15.0, 15.0);
266 h_pulls.SetXTitle("pulls");
267 for (unsigned module = 0; module < c_numModules; module++) {
268 for (unsigned channel = 0; channel < c_numChannels; channel++) {
269 std::vector<double> pos, err;
270 for (int i = 0; i < 2; i++) {
271 const auto& minimum = m_finders[i][module][channel].getMinimum();
272 if (not minimum.valid) continue;
273 pos.push_back(minimum.position);
274 err.push_back(minimum.error);
275 }
276 if (pos.size() < 2) continue;
277 double pull = (pos[0] - pos[1]) / sqrt(err[0] * err[0] + err[1] * err[1]);
278 h_pulls.Fill(pull);
279 }
280 }
281 h_pulls.Write();
282 double scaleError = h_pulls.GetRMS();
283
284 // merge two statistically independent finders and store results into histograms
285
286 for (unsigned module = 0; module < c_numModules; module++) {
287 int moduleID = module + 1;
288
289 std::string slotNum = std::to_string(moduleID);
290 if (moduleID < 10) slotNum = "0" + slotNum;
291
292 std::string name = "channelT0_slot" + slotNum;
293 std::string title = "Relative channel T0 for slot " + slotNum;
294 TH1F h_channelT0(name.c_str(), title.c_str(), c_numChannels, 0, c_numChannels);
295 h_channelT0.SetXTitle("channel number");
296 h_channelT0.SetYTitle("relative channel T0 [ns]");
297
298 for (unsigned channel = 0; channel < c_numChannels; channel++) {
299 auto& finder = m_finders[0][module][channel].add(m_finders[1][module][channel]);
300 const auto& minimum = finder.getMinimum();
301 if (minimum.valid) {
302 h_channelT0.SetBinContent(channel + 1, minimum.position);
303 h_channelT0.SetBinError(channel + 1, minimum.error * scaleError);
304 }
305 }
306
307 h_channelT0.Write();
308
309 }
310
311 // merge all finders of a slot to find module T0
312
313 TH1F h_moduleT0("moduleT0", "Relative module T0",
314 c_numModules, 0.5, static_cast<float>(c_numModules) + 0.5);
315 h_moduleT0.SetXTitle("slot number");
316 h_moduleT0.SetYTitle("relative module T0 [ns]");
317 auto finderCommon = m_finders[0][0][0];
318 finderCommon.clear();
319 for (unsigned module = 0; module < c_numModules; module++) {
320 auto finder = m_finders[0][module][0];
321 finder.clear();
322 for (unsigned channel = 0; channel < c_numChannels; channel++) {
323 finder.add(m_finders[0][module][channel]);
324 }
325 finderCommon.add(finder);
326 const auto& minimum = finder.getMinimum();
327 if (minimum.valid) {
328 h_moduleT0.SetBinContent(module + 1, minimum.position);
329 h_moduleT0.SetBinError(module + 1, minimum.error * scaleError);
330 }
331 }
332 h_moduleT0.Write();
333
334 // find common T0
335
336 TH1F h_commonT0("commonT0", "Relative common T0", 1, 0, 1);
337 h_commonT0.SetYTitle("relative common T0 [ns]");
338 const auto& minimum = finderCommon.getMinimum();
339 if (minimum.valid) {
340 h_commonT0.SetBinContent(1, minimum.position);
341 h_commonT0.SetBinError(1, minimum.error * scaleError);
342 }
343 h_commonT0.Write();
344
345 // write other histograms and ntuple; close the file
346
347 for (auto& h : m_hits1D) h.Write();
348 for (auto& h : m_hits2D) h.Write();
349
350 m_tree->Write();
351 m_file->Close();
352
353 B2RESULT("Results available in " << m_outFileName);
354 }
355
356
358} // end Belle2 namespace
359
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
StoreObjPtr< TOPRecBunch > m_recBunch
reconstructed bunch
TTree * m_tree
TTree containing selected track parameters etc.
TOP::TrackSelector m_selector
track selection utility
double m_sigmaSmear
additional smearing of PDF in [ns]
int m_PDG
track MC truth (simulated data only)
float m_phi
track: extrapolated hit momentum in local (module) frame
double m_minMomentum
minimal track momentum if sample is "cosmics"
int m_numBins
number of bins to which search region is divided
int m_numPhotons
number of photons in this slot
double m_maxZ
maximal local z of extrapolated hit
int m_moduleID
slot to which the track is extrapolated to
std::vector< TH1F > m_hits1D
number photon hits in a channel
TOP::PDFConstructor::EPDFOption m_PDFOption
PDF option.
double m_minZ
minimal local z of extrapolated hit
@ c_numChannels
number of channels per module
float m_p
track: extrapolated hit momentum in local (module) frame
StoreArray< Track > m_tracks
collection of tracks
double m_deltaEcms
c.m.s energy window if sample is "dimuon" or "bhabha"
double m_timeRange
time range in which to search for the minimum [ns]
float m_y
track: extrapolated hit coordinate in local (module) frame
StoreArray< TOPDigit > m_digits
collection of digits
StoreArray< ExtHit > m_extHits
collection of extrapolated hits
std::string m_outFileName
Root output file name containing results.
float m_z
track: extrapolated hit coordinate in local (module) frame
float m_x
track: extrapolated hit coordinate in local (module) frame
TOP::Chi2MinimumFinder1D m_finders[2][c_numModules][c_numChannels]
finders
float m_theta
track: extrapolated hit momentum in local (module) frame
Minimum finder using tabulated chi^2 values in one dimension.
const std::vector< double > & getBinCenters() const
Returns vector of bin centers.
void add(unsigned i, double chi2)
Add chi^2 value to bin i.
const Minimum & getMinimum()
Returns parabolic minimum.
void clear()
Set chi^2 values to zero.
PDF construction and log likelihood determination for a given track and particle hypothesis.
bool isValid() const
Checks the object status.
const std::vector< LogL > & getPixelLogLs(double t0, double sigt=0) const
Returns extended log likelihoods in pixels for PDF shifted in time.
@ c_Optimal
y dependent only where necessary
@ c_Fine
y dependent everywhere
@ c_Rough
no dependence on y
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.
static void setDefaultTimeWindow()
Sets default time window (functions getMinTime(), getMaxTime() will then return default values from D...
static double getMaxTime()
Returns time window upper edge.
static double getMinTime()
Returns time window lower edge.
Reconstructed track at TOP.
Definition: TOPTrack.h:39
int getPDGCode() const
Returns PDG code of associated MCParticle (returns 0 if none)
Definition: TOPTrack.h:228
double getCharge() const
Returns charge.
Definition: TOPTrack.h:161
bool isValid() const
Checks if track is successfully constructed.
Definition: TOPTrack.h:137
int getModuleID() const
Returns slot ID.
Definition: TOPTrack.h:143
Utility for the track selection - used in various calibration modules.
Definition: TrackSelector.h:27
double getCMSEnergy() const
Returns c.m.s energy of the track in last isSelected call.
const ROOT::Math::XYZPoint & getLocalPosition() const
Returns position at TOP in local frame of the track in last isSelected call.
void setDeltaEcms(double deltaEcms)
Sets cut on c.m.s.
Definition: TrackSelector.h:63
void setMinMomentum(double minMomentum)
Sets momentum cut (used for "cosmics" only)
Definition: TrackSelector.h:57
void setCutOnPOCA(double dr, double dz)
Sets cut on point of closest approach to (0, 0, 0)
Definition: TrackSelector.h:70
const ROOT::Math::XYZVector & getPOCAPosition() const
Returns position of POCA of the track in last isSelected call.
bool isSelected(const TOPTrack &track) const
Returns selection status.
const ROOT::Math::XYZVector & getLocalMomentum() const
Returns momentum at TOP in local frame of the track in last isSelected call.
void setCutOnLocalZ(double minZ, double maxZ)
Sets cut on local z coordinate (module frame) of the track extrapolated to TOP.
Definition: TrackSelector.h:81
const Const::ChargedStable & getChargedStable() const
Returns track hypothesis.
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
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
virtual void initialize() override
Initialize the Module.
virtual void event() override
Event processor.
virtual void terminate() override
Termination action.
Abstract base class for different kinds of events.
STL namespace.