Belle II Software  release-08-01-10
TOPModuleT0CalibratorModule.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/TOPModuleT0Calibrator/TOPModuleT0CalibratorModule.h>
10 #include <top/reconstruction_cpp/TOPTrack.h>
11 #include <top/reconstruction_cpp/PDFConstructor.h>
12 #include <top/reconstruction_cpp/TOPRecoManager.h>
13 #include <framework/dataobjects/EventMetaData.h>
14 #include <framework/logging/Logger.h>
15 #include <TRandom.h>
16 
17 using namespace std;
18 
19 namespace Belle2 {
24  using namespace TOP;
25 
26  //-----------------------------------------------------------------
28  //-----------------------------------------------------------------
29 
30  REG_MODULE(TOPModuleT0Calibrator);
31 
32  //-----------------------------------------------------------------
33  // Implementation
34  //-----------------------------------------------------------------
35 
36  TOPModuleT0CalibratorModule::TOPModuleT0CalibratorModule() : Module()
37 
38  {
39  // set module description
40  setDescription("Module T0 calibration with dimuons or bhabhas. "
41  "Useful when the geometrical alignment is fine.");
42  // setPropertyFlags(c_ParallelProcessingCertified);
43 
44  // Add parameters
45  addParam("numBins", m_numBins, "number of bins of the search region", 200);
46  addParam("timeRange", m_timeRange,
47  "time range in which to search for the minimum [ns]", 10.0);
48  addParam("sigmaSmear", m_sigmaSmear,
49  "sigma in [ns] for additional smearing of PDF", 0.0);
50  addParam("sample", m_sample,
51  "sample type: one of dimuon, bhabha or cosmics", std::string("dimuon"));
52  addParam("minMomentum", m_minMomentum,
53  "minimal track momentum if sample is cosmics", 1.0);
54  addParam("deltaEcms", m_deltaEcms,
55  "c.m.s energy window (half size) if sample is dimuon or bhabha", 0.1);
56  addParam("dr", m_dr, "cut on POCA in r", 2.0);
57  addParam("dz", m_dz, "cut on POCA in abs(z)", 4.0);
58  addParam("minZ", m_minZ,
59  "minimal local z of extrapolated hit", -130.0);
60  addParam("maxZ", m_maxZ,
61  "maximal local z of extrapolated hit", 130.0);
62  addParam("outputFileName", m_outFileName,
63  "Root output file name containing calibration results. "
64  "File name can include *'s; "
65  "they will be replaced with a run number from the first input file",
66  std::string("moduleT0_r*.root"));
67  addParam("pdfOption", m_pdfOption,
68  "PDF option, one of 'rough', 'fine', 'optimal'", std::string("rough"));
69 
70  }
71 
73  {
74  // input collections
75  m_digits.isRequired();
76  m_tracks.isRequired();
77  m_extHits.isRequired();
78  m_recBunch.isOptional();
79 
80  // toggle has changed status to prevent warning in beginRun for the first IOV
81  m_moduleT0.hasChanged();
82 
83  // Parse PDF option
84  if (m_pdfOption == "rough") {
86  } else if (m_pdfOption == "fine") {
88  } else if (m_pdfOption == "optimal") {
90  } else {
91  B2ERROR("TOPPDFDebuggerModule: unknown PDF option '" << m_pdfOption << "'");
92  }
93 
94  // set track selector
100 
101  // Chi2 minimum finders
102  double tmin = -m_timeRange / 2;
103  double tmax = m_timeRange / 2;
104  for (unsigned i = 0; i < 2; i++) {
105  for (unsigned m = 0; m < c_numModules; m++) {
106  m_finders[i][m] = Chi2MinimumFinder1D(m_numBins, tmin, tmax);
107  }
108  }
109 
110  // if file name includes *'s replace them with a run number
111  auto pos = m_outFileName.find("*");
112  if (pos != std::string::npos) {
113  StoreObjPtr<EventMetaData> evtMetaData;
114  auto run = std::to_string(evtMetaData->getRun());
115  while (run.size() < 5) run = "0" + run;
116  while (pos != std::string::npos) {
117  m_outFileName.replace(pos, 1, run);
118  pos = m_outFileName.find("*");
119  }
120  }
121 
122  // open root file for ntuple and histogram output
123  m_file = TFile::Open(m_outFileName.c_str(), "RECREATE");
124  if (not m_file) {
125  B2ERROR("Cannot open output file '" << m_outFileName << "'");
126  return;
127  }
128 
129  // histograms
130  m_hits1D = TH1F("numHits", "Number of photons per slot",
131  c_numModules, 0.5, static_cast<float>(c_numModules) + 0.5);
132  m_hits1D.SetXTitle("slot number");
133  m_hits1D.SetYTitle("hits per slot");
134 
135  m_hits2D = TH2F("timeHits", "Photon times vs. boardstacks",
136  c_numModules * 4, 0.5, static_cast<float>(c_numModules) + 0.5, 200, 0.0, 20.0);
137  m_hits2D.SetXTitle("slot number");
138  m_hits2D.SetYTitle("time [ns]");
139 
140  // create output tree
141  m_tree = new TTree("tree", "Channel T0 calibration results");
142  m_tree->Branch("slot", &m_moduleID);
143  m_tree->Branch("numPhotons", &m_numPhotons);
144  m_tree->Branch("x", &m_x);
145  m_tree->Branch("y", &m_y);
146  m_tree->Branch("z", &m_z);
147  m_tree->Branch("p", &m_p);
148  m_tree->Branch("theta", &m_theta);
149  m_tree->Branch("phi", &m_phi);
150  m_tree->Branch("r_poca", &m_pocaR);
151  m_tree->Branch("z_poca", &m_pocaZ);
152  m_tree->Branch("x_poca", &m_pocaX);
153  m_tree->Branch("y_poca", &m_pocaY);
154  m_tree->Branch("Ecms", &m_cmsE);
155  m_tree->Branch("charge", &m_charge);
156  m_tree->Branch("PDG", &m_PDG);
157 
158  }
159 
161  {
162  if (m_moduleT0.hasChanged()) {
163  StoreObjPtr<EventMetaData> evtMetaData;
164  B2WARNING("ModuleT0 payload has changed. Calibration results may not be correct."
165  << LogVar("run", evtMetaData->getRun()));
166  }
167  }
168 
170  {
171  /* check bunch reconstruction status and run alignment:
172  - if object exists and bunch is found (collision data w/ bunch finder in the path)
173  - if object doesn't exist (cosmic data and other cases w/o bunch finder)
174  */
175 
176  if (m_recBunch.isValid()) {
177  if (not m_recBunch->isReconstructed()) return;
178  }
179 
181  double timeMin = TOPRecoManager::getMinTime();
182  double timeMax = TOPRecoManager::getMaxTime();
183 
184  // running offset must not be subtracted in TOPDigits: issue an error if it is
185 
187  B2ERROR("Running offset subtracted in TOPDigits: module T0 will not be correct");
188  }
189 
190  // loop over reconstructed tracks, make a selection and accumulate log likelihoods
191 
192  for (const auto& track : m_tracks) {
193 
194  // track selection
195  TOPTrack trk(track);
196  if (not trk.isValid()) continue;
197 
198  if (not m_selector.isSelected(trk)) continue;
199 
200  // construct PDF
201  PDFConstructor pdfConstructor(trk, m_selector.getChargedStable(), m_PDFOption);
202  if (not pdfConstructor.isValid()) continue;
203 
204  // minimization procedure: accumulate
205  const unsigned module = trk.getModuleID() - 1;
206  if (module >= c_numModules) continue;
207  int sub = gRandom->Integer(2); // generate sub-sample number
208  auto& finder = m_finders[sub][module];
209  const auto& binCenters = finder.getBinCenters();
210  for (unsigned ibin = 0; ibin < binCenters.size(); ibin++) {
211  double t0 = binCenters[ibin];
212  finder.add(ibin, -2 * pdfConstructor.getLogL(t0, m_sigmaSmear).logL);
213  }
214 
215  // fill histograms of hits
216  m_numPhotons = 0;
217  for (const auto& digit : m_digits) {
218  if (digit.getHitQuality() != TOPDigit::c_Good) continue;
219  if (digit.getModuleID() != trk.getModuleID()) continue;
220  if (digit.getTime() < timeMin) continue;
221  if (digit.getTime() > timeMax) continue;
222  m_numPhotons++;
223  m_hits1D.Fill(digit.getModuleID());
224  int bs = digit.getBoardstackNumber();
225  m_hits2D.Fill((digit.getModuleID() * 4 + bs - 1.5) / 4.0, digit.getTime());
226  }
227 
228  // fill output tree
229  m_moduleID = trk.getModuleID();
230  const auto& localPosition = m_selector.getLocalPosition();
231  m_x = localPosition.X();
232  m_y = localPosition.Y();
233  m_z = localPosition.Z();
234  const auto& localMomentum = m_selector.getLocalMomentum();
235  m_p = localMomentum.R();
236  m_theta = localMomentum.Theta();
237  m_phi = localMomentum.Phi();
238  const auto& pocaPosition = m_selector.getPOCAPosition();
239  m_pocaR = pocaPosition.Rho();
240  m_pocaZ = pocaPosition.Z();
241  m_pocaX = pocaPosition.X();
242  m_pocaY = pocaPosition.Y();
244  m_charge = trk.getCharge();
245  m_PDG = trk.getPDGCode();
246  m_tree->Fill();
247  }
248 
249  }
250 
251 
253  {
254  // determine scaling factor for errors from two statistically independent results
255 
256  TH1F h_pulls("pulls", "Pulls of the two statistically independent results",
257  200, -15.0, 15.0);
258  h_pulls.SetXTitle("pulls");
259  for (unsigned module = 0; module < c_numModules; module++) {
260  std::vector<double> pos, err;
261  for (int i = 0; i < 2; i++) {
262  const auto& minimum = m_finders[i][module].getMinimum();
263  if (not minimum.valid) continue;
264  pos.push_back(minimum.position);
265  err.push_back(minimum.error);
266  }
267  if (pos.size() < 2) continue;
268  double pull = (pos[0] - pos[1]) / sqrt(err[0] * err[0] + err[1] * err[1]);
269  h_pulls.Fill(pull);
270  }
271  h_pulls.Write();
272  double scaleError = h_pulls.GetRMS();
273 
274  // merge two statistically independent finders and store results into histograms
275 
276  TH1F h_relModuleT0("relModuleT0", "Module T0 relative to calibration",
277  c_numModules, 0.5, static_cast<float>(c_numModules) + 0.5);
278  h_relModuleT0.SetXTitle("slot number");
279  h_relModuleT0.SetYTitle("module T0 residual [ns]");
280 
281  for (unsigned module = 0; module < c_numModules; module++) {
282  auto& finder = m_finders[0][module].add(m_finders[1][module]);
283  const auto& minimum = finder.getMinimum();
284  if (minimum.valid) {
285  h_relModuleT0.SetBinContent(module + 1, minimum.position);
286  h_relModuleT0.SetBinError(module + 1, minimum.error * scaleError);
287  }
288  std::string slotNum = to_string(module + 1);
289  auto h = finder.getHistogram("chi2_slot_" + slotNum, "slot " + slotNum);
290  h.Write();
291  }
292  h_relModuleT0.Write();
293 
294  // absolute module T0
295  TH1F h_moduleT0("moduleT0", "Module T0",
296  c_numModules, 0.5, static_cast<float>(c_numModules) + 0.5);
297  h_moduleT0.SetXTitle("slot number");
298  h_moduleT0.SetYTitle("module T0 [ns]");
299 
300  double average = 0; // average to be subtracted
301  int num = 0;
302  for (int slot = 1; slot <= c_numModules; slot++) {
303  if (h_relModuleT0.GetBinError(slot) > 0) {
304  double T0 = 0;
305  if (m_moduleT0->isCalibrated(slot)) T0 = m_moduleT0->getT0(slot);
306  average += h_relModuleT0.GetBinContent(slot) + T0;
307  num++;
308  }
309  }
310  if (num > 0) average /= num;
311 
312  for (int slot = 1; slot <= c_numModules; slot++) {
313  double T0 = 0;
314  if (m_moduleT0->isCalibrated(slot)) T0 = m_moduleT0->getT0(slot);
315  double t0 = h_relModuleT0.GetBinContent(slot) + T0 - average;
316  double err = h_relModuleT0.GetBinError(slot);
317  h_moduleT0.SetBinContent(slot, t0);
318  h_moduleT0.SetBinError(slot, err);
319  }
320  h_moduleT0.Write();
321 
322  // write other histograms and ntuple; close the file
323 
324  m_hits1D.Write();
325  m_hits2D.Write();
326  m_tree->Write();
327  m_file->Close();
328 
329  B2RESULT(num << "/16 calibrated. Results available in " << m_outFileName);
330  }
331 
332 
334  {
335  for (const auto& digit : m_digits) {
336  if (digit.hasStatus(TOPDigit::c_BunchOffsetSubtracted)) return true;
337  }
338  return false;
339  }
340 
342 } // end Belle2 namespace
343 
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
TOP::PDFConstructor::EPDFOption m_PDFOption
PDF option.
double m_minZ
minimal local z of extrapolated hit
float m_p
track: extrapolated hit momentum in local (module) frame
StoreArray< Track > m_tracks
collection of tracks
TOP::Chi2MinimumFinder1D m_finders[2][c_numModules]
finders
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
DBObjPtr< TOPCalModuleT0 > m_moduleT0
module T0 calibration constants
std::string m_outFileName
Root output file name containing results.
TH1F m_hits1D
number photon hits in a slot
float m_z
track: extrapolated hit coordinate in local (module) frame
float m_x
track: extrapolated hit coordinate in local (module) frame
float m_theta
track: extrapolated hit momentum in local (module) frame
Minimum finder using tabulated chi^2 values in one dimension.
void add(unsigned i, double chi2)
Add chi^2 value to bin i.
const Minimum & getMinimum()
Returns parabolic minimum.
PDF construction and log likelihood determination for a given track and particle hypothesis.
LogL getLogL() const
Returns extended log likelihood (using the default time window)
bool isValid() const
Checks the object status.
@ c_Optimal
y dependent only where necessary
@ c_Fine
y dependent everywhere
@ c_Rough
no dependence on y
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.
void setDeltaEcms(double deltaEcms)
Sets cut on c.m.s.
Definition: TrackSelector.h:63
const ROOT::Math::XYZVector & getPOCAPosition() const
Returns position of POCA of the track in last isSelected call.
void setMinMomentum(double minMomentum)
Sets momentum cut (used for "cosmics" only)
Definition: TrackSelector.h:57
const ROOT::Math::XYZVector & getLocalMomentum() const
Returns momentum at TOP in local frame of the track in last isSelected call.
void setCutOnPOCA(double dr, double dz)
Sets cut on point of closest approach to (0, 0, 0)
Definition: TrackSelector.h:70
bool isSelected(const TOPTrack &track) const
Returns selection status.
const ROOT::Math::XYZPoint & getLocalPosition() const
Returns position at TOP in local frame of the track in last isSelected call.
const Const::ChargedStable & getChargedStable() const
Returns track hypothesis.
void setCutOnLocalZ(double minZ, double maxZ)
Sets cut on local z coordinate (module frame) of the track extrapolated to TOP.
Definition: TrackSelector.h:81
Class to store variables with their name which were sent to the logging service.
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.
bool isRunningOffsetSubtracted()
Checks if running offset is subtracted in TOPDigits.
virtual void beginRun() override
Called when entering a new run.
Abstract base class for different kinds of events.
double logL
extended log likelihood