Belle II Software  release-08-01-10
TOPCommonT0LLCollectorModule.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/TOPCommonT0LLCollectorModule.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 
15 // framework aux
16 #include <framework/gearbox/Const.h>
17 #include <framework/logging/Logger.h>
18 
19 // root
20 #include <TRandom.h>
21 #include <TH1F.h>
22 #include <TH1D.h>
23 #include <TH2F.h>
24 
25 using namespace std;
26 
27 namespace Belle2 {
33  using namespace TOP;
34 
35  //-----------------------------------------------------------------
37  //-----------------------------------------------------------------
38 
39  REG_MODULE(TOPCommonT0LLCollector);
40 
41  //-----------------------------------------------------------------
42  // Implementation
43  //-----------------------------------------------------------------
44 
45  TOPCommonT0LLCollectorModule::TOPCommonT0LLCollectorModule()
46  {
47  // set module description and processing properties
48  setDescription("A collector for common T0 calibration with dimuons or Bhabha's using "
49  "neg. log likelihood minimization (method LL)");
50  setPropertyFlags(c_ParallelProcessingCertified);
51 
52  // module parameters
53  addParam("bunchesPerSSTclk", m_bunchesPerSSTclk,
54  "number of bunches per SST clock period", 24);
55  addParam("numBins", m_numBins, "number of bins of the search region", 200);
56  addParam("timeRange", m_timeRange,
57  "time range in which to search for the minimum [ns]", 10.0);
58  addParam("sigmaSmear", m_sigmaSmear,
59  "sigma in [ns] for additional smearing of PDF", 0.0);
60  addParam("sample", m_sample,
61  "sample type: one of dimuon or bhabha", std::string("dimuon"));
62  addParam("deltaEcms", m_deltaEcms,
63  "c.m.s energy window (half size) if sample is dimuon or bhabha", 0.1);
64  addParam("dr", m_dr, "cut on POCA in r", 2.0);
65  addParam("dz", m_dz, "cut on POCA in abs(z)", 4.0);
66  addParam("minZ", m_minZ,
67  "minimal local z of extrapolated hit", -130.0);
68  addParam("maxZ", m_maxZ,
69  "maximal local z of extrapolated hit", 130.0);
70  addParam("pdfOption", m_pdfOption,
71  "PDF option, one of 'rough', 'fine', 'optimal'", std::string("rough"));
72  }
73 
74 
75  void TOPCommonT0LLCollectorModule::prepare()
76  {
77  // input collections
78 
79  m_digits.isRequired();
80  m_tracks.isRequired();
81  m_extHits.isRequired();
82  m_recBunch.isRequired();
83 
84  // bunch time separation
85 
86  const auto* geo = TOPGeometryPar::Instance()->getGeometry();
87  m_bunchTimeSep = geo->getNominalTDC().getSyncTimeBase() / m_bunchesPerSSTclk;
88 
89  // Parse PDF option
90 
91  if (m_pdfOption == "rough") {
92  m_PDFOption = PDFConstructor::c_Rough;
93  } else if (m_pdfOption == "fine") {
94  m_PDFOption = PDFConstructor::c_Fine;
95  } else if (m_pdfOption == "optimal") {
96  m_PDFOption = PDFConstructor::c_Optimal;
97  } else {
98  B2ERROR("Unknown PDF option '" << m_pdfOption << "'");
99  }
100 
101  // set track selector
102 
103  if (m_sample == "dimuon" or m_sample == "bhabha") {
104  m_selector = TrackSelector(m_sample);
105  m_selector.setDeltaEcms(m_deltaEcms);
106  m_selector.setCutOnPOCA(m_dr, m_dz);
107  m_selector.setCutOnLocalZ(m_minZ, m_maxZ);
108  } else {
109  B2ERROR("Invalid sample type '" << m_sample << "'");
110  }
111 
112  // create and register histograms
113 
114  double tmin = -m_timeRange / 2;
115  double tmax = m_timeRange / 2;
116  for (unsigned i = 0; i < c_numSets; i++) {
117  string name = "chi2_set" + to_string(i);
118  auto h = new TH1D(name.c_str(), "chi2 scan; t0 [ns]; chi2", m_numBins, tmin, tmax);
119  registerObject<TH1D>(name, h);
120  m_names.push_back(name);
121  }
122 
123  auto h1 = new TH1F("tracks_per_set", "tracks per sample; sample number; num tracks",
124  c_numSets, 0, c_numSets);
125  registerObject<TH1F>("tracks_per_set", h1);
126 
127  auto h2 = new TH1F("numHits", "Number of photons per slot",
128  c_numModules, 0.5, static_cast<float>(c_numModules) + 0.5);
129  h2->SetXTitle("slot number");
130  h2->SetYTitle("hits per slot");
131  registerObject<TH1F>("numHits", h2);
132 
133  auto h3 = new TH2F("timeHits", "Photon times vs. boardstacks",
134  c_numModules * 4, 0.5, static_cast<float>(c_numModules) + 0.5, 200, 0.0, 20.0);
135  h3->SetXTitle("slot number");
136  h3->SetYTitle("time [ns]");
137  registerObject<TH2F>("timeHits", h3);
138 
139  // this one is needed primarely to pass bunch time separation to the algorithm,
140  // since DB interface doesn't work there
141  auto h4 = new TH1F("offset", "current offset from input files; offset [ns]",
142  200, -m_bunchTimeSep / 2, m_bunchTimeSep / 2);
143  registerObject<TH1F>("offset", h4);
144 
145  }
146 
147 
148  void TOPCommonT0LLCollectorModule::collect()
149  {
150  // bunch must be reconstructed
151 
152  if (not m_recBunch.isValid()) return;
153  if (not m_recBunch->isReconstructed()) return;
154 
155  TOPRecoManager::setDefaultTimeWindow();
156  double timeMin = TOPRecoManager::getMinTime();
157  double timeMax = TOPRecoManager::getMaxTime();
158 
159  // correct-back digits for common T0
160 
161  double T0 = 0;
162  if (m_commonT0->isCalibrated()) T0 = m_commonT0->getT0();
163  for (auto& digit : m_digits) {
164  if (digit.isCommonT0Calibrated()) digit.subtractT0(-T0);
165  if (digit.hasStatus(TOPDigit::c_BunchOffsetSubtracted)) digit.subtractT0(-m_recBunch->getAverageOffset());
166  }
167 
168  // loop over reconstructed tracks, make a selection and accumulate log likelihoods
169 
170  int ntra = 0;
171  for (const auto& track : m_tracks) {
172 
173  // track selection
174  TOPTrack trk(track);
175  if (not trk.isValid()) continue;
176 
177  if (not m_selector.isSelected(trk)) continue;
178 
179  // construct PDF
180  PDFConstructor pdfConstructor(trk, m_selector.getChargedStable(), m_PDFOption);
181  if (not pdfConstructor.isValid()) continue;
182 
183  // minimization procedure: accumulate
184  int sub = gRandom->Integer(c_numSets); // generate sub-sample number
185  auto h = getObjectPtr<TH1D>(m_names[sub]);
186  for (int ibin = 0; ibin < h->GetNbinsX(); ibin++) {
187  double t0 = h->GetBinCenter(ibin + 1);
188  double chi = h->GetBinContent(ibin + 1);
189  chi += -2 * pdfConstructor.getLogL(t0, m_sigmaSmear).logL;
190  h->SetBinContent(ibin + 1, chi);
191  }
192  auto h1 = getObjectPtr<TH1F>("tracks_per_set");
193  h1->Fill(sub);
194 
195  // fill histograms of hits
196  auto h2 = getObjectPtr<TH1F>("numHits");
197  auto h3 = getObjectPtr<TH2F>("timeHits");
198  for (const auto& digit : m_digits) {
199  if (digit.getHitQuality() != TOPDigit::c_Good) continue;
200  if (digit.getModuleID() != trk.getModuleID()) continue;
201  if (digit.getTime() < timeMin) continue;
202  if (digit.getTime() > timeMax) continue;
203  h2->Fill(digit.getModuleID());
204  int bs = digit.getBoardstackNumber();
205  h3->Fill((digit.getModuleID() * 4 + bs - 1.5) / 4.0, digit.getTime());
206  }
207  ntra++;
208  }
209 
210  // fill just another control histogram
211 
212  if (ntra > 0) {
213  auto h4 = getObjectPtr<TH1F>("offset");
214  h4->Fill(m_recBunch->getCurrentOffset());
215  }
216 
217  }
218 
220 }
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:39
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
#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