Belle II Software  release-06-01-15
TOPCommonT0CalibratorModule.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/TOPCommonT0Calibrator/TOPCommonT0CalibratorModule.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/logging/Logger.h>
16 #include <TRandom.h>
17 
18 using namespace std;
19 
20 namespace Belle2 {
25  using namespace TOP;
26 
27  //-----------------------------------------------------------------
28  // Register module
29  //-----------------------------------------------------------------
30 
31  REG_MODULE(TOPCommonT0Calibrator)
32 
33  //-----------------------------------------------------------------
34  // Implementation
35  //-----------------------------------------------------------------
36 
38 
39  {
40  // set module description
41  setDescription("Common T0 calibration with dimuons or bhabhas.");
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("commonT0_r*.root"));
67  addParam("pdfOption", m_pdfOption,
68  "PDF option, one of 'rough', 'fine', 'optimal'", std::string("rough"));
69 
70  }
71 
72 
73  void TOPCommonT0CalibratorModule::initialize()
74  {
75  // input collections
76  m_digits.isRequired();
77  m_tracks.isRequired();
78  m_extHits.isRequired();
79  m_recBunch.isOptional();
80 
81  // bunch separation in time
82 
83  const auto* geo = TOPGeometryPar::Instance()->getGeometry();
84  m_bunchTimeSep = geo->getNominalTDC().getSyncTimeBase() / 24;
85 
86  // Parse PDF option
87  if (m_pdfOption == "rough") {
88  m_PDFOption = PDFConstructor::c_Rough;
89  } else if (m_pdfOption == "fine") {
90  m_PDFOption = PDFConstructor::c_Fine;
91  } else if (m_pdfOption == "optimal") {
92  m_PDFOption = PDFConstructor::c_Optimal;
93  } else {
94  B2ERROR("TOPPDFDebuggerModule: unknown PDF option '" << m_pdfOption << "'");
95  }
96 
97  // set track selector
98  m_selector = TrackSelector(m_sample);
99  m_selector.setMinMomentum(m_minMomentum);
100  m_selector.setDeltaEcms(m_deltaEcms);
101  m_selector.setCutOnPOCA(m_dr, m_dz);
102  m_selector.setCutOnLocalZ(m_minZ, m_maxZ);
103 
104  // Chi2 minimum finders
105  double tmin = -m_timeRange / 2;
106  double tmax = m_timeRange / 2;
107  for (unsigned i = 0; i < c_numSets; i++) {
108  m_finders[i] = Chi2MinimumFinder1D(m_numBins, tmin, tmax);
109  }
110 
111  // if file name includes *'s replace them with a run number
112  auto pos = m_outFileName.find("*");
113  if (pos != std::string::npos) {
114  StoreObjPtr<EventMetaData> evtMetaData;
115  auto run = std::to_string(evtMetaData->getRun());
116  while (run.size() < 5) run = "0" + run;
117  while (pos != std::string::npos) {
118  m_outFileName.replace(pos, 1, run);
119  pos = m_outFileName.find("*");
120  }
121  }
122 
123  // open root file for ntuple and histogram output
124  m_file = TFile::Open(m_outFileName.c_str(), "RECREATE");
125  if (not m_file) {
126  B2ERROR("Cannot open output file '" << m_outFileName << "'");
127  return;
128  }
129 
130  // control histograms
131  m_hits1D = TH1F("numHits", "Number of photons per slot",
132  c_numModules, 0.5, c_numModules + 0.5);
133  m_hits1D.SetXTitle("slot number");
134  m_hits1D.SetYTitle("hits per slot");
135 
136  m_hits2D = TH2F("timeHits", "Photon times vs. boardstacks",
137  c_numModules * 4, 0.5, c_numModules + 0.5, 200, 0.0, 20.0);
138  m_hits2D.SetXTitle("slot number");
139  m_hits2D.SetYTitle("time [ns]");
140 
141  // create output tree
142  m_tree = new TTree("tree", "Channel T0 calibration results");
143  m_tree->Branch("slot", &m_moduleID);
144  m_tree->Branch("numPhotons", &m_numPhotons);
145  m_tree->Branch("x", &m_x);
146  m_tree->Branch("y", &m_y);
147  m_tree->Branch("z", &m_z);
148  m_tree->Branch("p", &m_p);
149  m_tree->Branch("theta", &m_theta);
150  m_tree->Branch("phi", &m_phi);
151  m_tree->Branch("r_poca", &m_pocaR);
152  m_tree->Branch("z_poca", &m_pocaZ);
153  m_tree->Branch("x_poca", &m_pocaX);
154  m_tree->Branch("y_poca", &m_pocaY);
155  m_tree->Branch("Ecms", &m_cmsE);
156  m_tree->Branch("charge", &m_charge);
157  m_tree->Branch("PDG", &m_PDG);
158 
159  }
160 
161 
162  void TOPCommonT0CalibratorModule::event()
163  {
164  /* check bunch reconstruction status and run alignment:
165  - if object exists and bunch is found (collision data w/ bunch finder in the path)
166  - if object doesn't exist (cosmic data and other cases w/o bunch finder)
167  */
168 
169  if (m_recBunch.isValid()) {
170  if (not m_recBunch->isReconstructed()) return;
171  }
172 
173  TOPRecoManager::setDefaultTimeWindow();
174  double timeMin = TOPRecoManager::getMinTime();
175  double timeMax = TOPRecoManager::getMaxTime();
176 
177  // running offset must not be subtracted in TOPDigits: issue an error if it is
178 
179  if (isRunningOffsetSubtracted()) {
180  B2ERROR("Running offset subtracted in TOPDigits: common T0 will not be correct");
181  }
182 
183  // loop over reconstructed tracks, make a selection and accumulate log likelihoods
184 
185  for (const auto& track : m_tracks) {
186 
187  // track selection
188  TOPTrack trk(track);
189  if (not trk.isValid()) continue;
190 
191  if (not m_selector.isSelected(trk)) continue;
192 
193  // construct PDF
194  PDFConstructor pdfConstructor(trk, m_selector.getChargedStable(), m_PDFOption);
195  if (not pdfConstructor.isValid()) continue;
196 
197  // minimization procedure: accumulate
198  int sub = gRandom->Integer(c_numSets); // generate sub-sample number
199  auto& finder = m_finders[sub];
200  const auto& binCenters = finder.getBinCenters();
201  for (unsigned ibin = 0; ibin < binCenters.size(); ibin++) {
202  double t0 = binCenters[ibin];
203  finder.add(ibin, -2 * pdfConstructor.getLogL(t0, m_sigmaSmear).logL);
204  }
205 
206  // fill histograms of hits
207  m_numPhotons = 0;
208  for (const auto& digit : m_digits) {
209  if (digit.getHitQuality() != TOPDigit::c_Good) continue;
210  if (digit.getModuleID() != trk.getModuleID()) continue;
211  if (digit.getTime() < timeMin) continue;
212  if (digit.getTime() > timeMax) continue;
213  m_numPhotons++;
214  m_hits1D.Fill(digit.getModuleID());
215  int bs = digit.getBoardstackNumber();
216  m_hits2D.Fill((digit.getModuleID() * 4 + bs - 1.5) / 4.0 , digit.getTime());
217  }
218 
219  // fill output tree
220  m_moduleID = trk.getModuleID();
221  const auto& localPosition = m_selector.getLocalPosition();
222  m_x = localPosition.X();
223  m_y = localPosition.Y();
224  m_z = localPosition.Z();
225  const auto& localMomentum = m_selector.getLocalMomentum();
226  m_p = localMomentum.Mag();
227  m_theta = localMomentum.Theta();
228  m_phi = localMomentum.Phi();
229  const auto& pocaPosition = m_selector.getPOCAPosition();
230  m_pocaR = pocaPosition.Perp();
231  m_pocaZ = pocaPosition.Z();
232  m_pocaX = pocaPosition.X();
233  m_pocaY = pocaPosition.Y();
234  m_cmsE = m_selector.getCMSEnergy();
235  m_charge = trk.getCharge();
236  m_PDG = trk.getPDGCode();
237  m_tree->Fill();
238  }
239 
240  }
241 
242 
243  void TOPCommonT0CalibratorModule::terminate()
244  {
245 
246  // determine scaling factor for errors from statistically independent results
247 
248  TH1F h_pulls("pulls", "Pulls of statistically independent results",
249  200, -15.0, 15.0);
250  h_pulls.SetXTitle("pulls");
251  std::vector<double> pos, err;
252  for (int i = 0; i < c_numSets; i++) {
253  const auto& minimum = m_finders[i].getMinimum();
254  if (not minimum.valid) continue;
255  pos.push_back(minimum.position);
256  err.push_back(minimum.error);
257  }
258  for (unsigned i = 0; i < pos.size(); i++) {
259  for (unsigned j = i + 1; j < pos.size(); j++) {
260  double pull = (pos[i] - pos[j]) / sqrt(err[i] * err[i] + err[j] * err[j]);
261  h_pulls.Fill(pull);
262  }
263  }
264  h_pulls.Write();
265  double scaleError = 1;
266  if (h_pulls.GetEntries() > 1) scaleError = h_pulls.GetRMS();
267 
268  // merge statistically independent finders and store results into histograms
269 
270  auto finder = m_finders[0];
271  for (int i = 1; i < c_numSets; i++) {
272  finder.add(m_finders[i]);
273  }
274 
275  TH1F h_relCommonT0("relCommonT0", "relative common T0", 1, 0, 1);
276  h_relCommonT0.SetYTitle("common T0 residual [ns]");
277  TH1F h_commonT0("commonT0", "Common T0", 1, 0, 1);
278  h_commonT0.SetYTitle("common T0 [ns]");
279 
280  const auto& minimum = finder.getMinimum();
281  auto h = finder.getHistogram("chi2", "chi2");
282  h.Write();
283  if (minimum.valid) {
284  h_relCommonT0.SetBinContent(1, minimum.position);
285  h_relCommonT0.SetBinError(1, minimum.error * scaleError);
286  double T0 = minimum.position;
287  if (m_commonT0->isCalibrated()) T0 += m_commonT0->getT0();
288  T0 -= round(T0 / m_bunchTimeSep) * m_bunchTimeSep; // wrap around
289  h_commonT0.SetBinContent(1, T0);
290  h_commonT0.SetBinError(1, minimum.error * scaleError);
291  }
292  h_relCommonT0.Write();
293  h_commonT0.Write();
294 
295  // write other histograms and ntuple; close the file
296 
297  m_hits1D.Write();
298  m_hits2D.Write();
299  m_tree->Write();
300  m_file->Close();
301 
302  B2RESULT("Results available in " << m_outFileName);
303  }
304 
305  bool TOPCommonT0CalibratorModule::isRunningOffsetSubtracted()
306  {
307  for (const auto& digit : m_digits) {
308  if (digit.hasStatus(TOPDigit::c_BunchOffsetSubtracted)) return true;
309  }
310  return false;
311  }
312 
314 } // end Belle2 namespace
315 
Base class for Modules.
Definition: Module.h:72
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:95
A module for common T0 calibration with collision data (dimuons or bhabhas)
Minimum finder using tabulated chi^2 values in one dimension.
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.
Reconstructed track at TOP.
Definition: TOPTrack.h:40
int getPDGCode() const
Returns PDG code of associated MCParticle (returns 0 if none)
Definition: TOPTrack.h:230
double getCharge() const
Returns charge.
Definition: TOPTrack.h:163
bool isValid() const
Checks if track is successfully constructed.
Definition: TOPTrack.h:139
int getModuleID() const
Returns slot ID.
Definition: TOPTrack.h:145
Utility for the track selection - used in various calibration modules.
Definition: TrackSelector.h:26
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.
double logL
extended log likelihood