Belle II Software  release-05-02-19
TOPModuleT0LLCollectorModule.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2019 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Marko Staric *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #include <top/modules/collectors/TOPModuleT0LLCollectorModule.h>
12 #include <top/geometry/TOPGeometryPar.h>
13 #include <top/reconstruction/TOPreco.h>
14 #include <top/reconstruction/TOPtrack.h>
15 #include <top/reconstruction/TOPconfigure.h>
16 
17 // framework aux
18 #include <framework/gearbox/Const.h>
19 #include <framework/logging/Logger.h>
20 
21 // root
22 #include <TRandom.h>
23 #include <TH1F.h>
24 #include <TH1D.h>
25 #include <TH2F.h>
26 
27 using namespace std;
28 
29 namespace Belle2 {
35  using namespace TOP;
36 
37  //-----------------------------------------------------------------
38  // Register module
39  //-----------------------------------------------------------------
40 
41  REG_MODULE(TOPModuleT0LLCollector)
42 
43  //-----------------------------------------------------------------
44  // Implementation
45  //-----------------------------------------------------------------
46 
48  {
49  // set module description and processing properties
50  setDescription("A collector for module T0 calibration with neg. log likelihood "
51  "minimization (method LL). Aimed for the final (precise) calibration"
52  "after initial one with DeltaT method.");
53  setPropertyFlags(c_ParallelProcessingCertified);
54 
55  // module parameters
56  addParam("numBins", m_numBins, "number of bins of the search region", 200);
57  addParam("timeRange", m_timeRange,
58  "time range in which to search for the minimum [ns]", 10.0);
59  addParam("minBkgPerBar", m_minBkgPerBar,
60  "minimal number of background photons per module", 0.0);
61  addParam("scaleN0", m_scaleN0,
62  "Scale factor for figure-of-merit N0", 1.0);
63  addParam("sigmaSmear", m_sigmaSmear,
64  "sigma in [ns] for additional smearing of PDF", 0.0);
65  addParam("sample", m_sample,
66  "sample type: one of dimuon or bhabha", std::string("dimuon"));
67  addParam("deltaEcms", m_deltaEcms,
68  "c.m.s energy window (half size) if sample is dimuon or bhabha", 0.1);
69  addParam("dr", m_dr, "cut on POCA in r", 2.0);
70  addParam("dz", m_dz, "cut on POCA in abs(z)", 4.0);
71  addParam("minZ", m_minZ,
72  "minimal local z of extrapolated hit", -130.0);
73  addParam("maxZ", m_maxZ,
74  "maximal local z of extrapolated hit", 130.0);
75  addParam("pdfOption", m_pdfOption,
76  "PDF option, one of 'rough', 'fine', 'optimal'", std::string("rough"));
77 
78  }
79 
80 
81  void TOPModuleT0LLCollectorModule::prepare()
82  {
83  // input collections
84 
85  m_digits.isRequired();
86  m_tracks.isRequired();
87  m_extHits.isRequired();
88  m_recBunch.isRequired();
89 
90  // Configure TOP detector for reconstruction
91 
92  TOPconfigure config;
93 
94  // Parse PDF option
95 
96  if (m_pdfOption == "rough") {
97  m_PDFOption = TOPreco::c_Rough;
98  } else if (m_pdfOption == "fine") {
99  m_PDFOption = TOPreco::c_Fine;
100  } else if (m_pdfOption == "optimal") {
101  m_PDFOption = TOPreco::c_Optimal;
102  } else {
103  B2ERROR("Unknown PDF option '" << m_pdfOption << "'");
104  }
105 
106  // set track selector
107 
108  if (m_sample == "dimuon" or m_sample == "bhabha") {
109  m_selector = TrackSelector(m_sample);
110  m_selector.setDeltaEcms(m_deltaEcms);
111  m_selector.setCutOnPOCA(m_dr, m_dz);
112  m_selector.setCutOnLocalZ(m_minZ, m_maxZ);
113  } else {
114  B2ERROR("Invalid sample type '" << m_sample << "'");
115  }
116 
117  // create and register histograms
118 
119  double tmin = -m_timeRange / 2;
120  double tmax = m_timeRange / 2;
121  for (unsigned i = 0; i < c_numSets; i++) {
122  for (int slot = 1; slot <= c_numModules; slot++) {
123  string name = "chi2_set" + to_string(i) + "_slot" + to_string(slot);
124  string title = "chi2 scan, slot" + to_string(slot) + "; t0 [ns]; chi2";
125  auto h = new TH1D(name.c_str(), title.c_str(), m_numBins, tmin, tmax);
126  registerObject<TH1D>(name, h);
127  m_names[i].push_back(name);
128  }
129  }
130 
131  auto h1 = new TH2F("tracks_per_slot", "tracks per slot and sample",
132  c_numModules, 0.5, c_numModules + 0.5, c_numSets, 0, c_numSets);
133  h1->SetXTitle("slot number");
134  h1->SetYTitle("sample number");
135  registerObject<TH2F>("tracks_per_slot", h1);
136 
137  auto h2 = new TH1F("numHits", "Number of photons per slot",
138  c_numModules, 0.5, c_numModules + 0.5);
139  h2->SetXTitle("slot number");
140  h2->SetYTitle("hits per slot");
141  registerObject<TH1F>("numHits", h2);
142 
143  auto h3 = new TH2F("timeHits", "Photon times vs. boardstacks",
144  c_numModules * 4, 0.5, c_numModules + 0.5, 200, 0.0, 20.0);
145  h3->SetXTitle("slot number");
146  h3->SetYTitle("time [ns]");
147  registerObject<TH2F>("timeHits", h3);
148 
149  const auto* geo = TOPGeometryPar::Instance()->getGeometry();
150  double bunchTimeSep = geo->getNominalTDC().getSyncTimeBase() / 24;
151  auto h4 = new TH1F("offset", "current offset from input files; offset [ns]",
152  200, -bunchTimeSep / 2, bunchTimeSep / 2);
153  registerObject<TH1F>("offset", h4);
154 
155  }
156 
157 
158  void TOPModuleT0LLCollectorModule::collect()
159  {
160  // bunch must be reconstructed
161 
162  if (not m_recBunch.isValid()) return;
163  if (not m_recBunch->isReconstructed()) return;
164 
165  // create reconstruction object and set various options
166 
167  int Nhyp = 1;
168  double mass = m_selector.getChargedStable().getMass();
169  int pdg = m_selector.getChargedStable().getPDGCode();
170  TOPreco reco(Nhyp, &mass, &pdg, m_minBkgPerBar, m_scaleN0);
171  reco.setPDFoption(m_PDFOption);
172  const auto& tdc = TOPGeometryPar::Instance()->getGeometry()->getNominalTDC();
173  double timeMin = tdc.getTimeMin();
174  double timeMax = tdc.getTimeMax();
175 
176  // add photon hits to reconstruction object with time corrected-back for module T0
177 
178  for (const auto& digit : m_digits) {
179  if (digit.getHitQuality() == TOPDigit::c_Good) {
180  double t = digit.getTime();
181  auto slot = digit.getModuleID();
182  if (digit.isModuleT0Calibrated()) t += m_moduleT0->getT0(slot);
183  if (digit.hasStatus(TOPDigit::c_BunchOffsetSubtracted)) {
184  t += m_recBunch->getAverageOffset();
185  }
186  reco.addData(digit.getModuleID(), digit.getPixelID(), t, digit.getTimeError());
187  }
188  }
189 
190  // loop over reconstructed tracks, make a selection and accumulate log likelihoods
191 
192  int ntra = 0;
193  for (const auto& track : m_tracks) {
194 
195  // track selection
196  TOPtrack trk(&track);
197  if (not trk.isValid()) continue;
198 
199  if (not m_selector.isSelected(trk)) continue;
200 
201  // run reconstruction
202  reco.reconstruct(trk);
203  if (reco.getFlag() != 1) continue; // track is not in the acceptance of TOP
204 
205  // minimization procedure: accumulate
206  int sub = gRandom->Integer(c_numSets); // generate sub-sample number
207  unsigned module = trk.getModuleID() - 1;
208  if (module >= m_names[sub].size()) continue;
209  auto h = getObjectPtr<TH1D>(m_names[sub][module]);
210  for (int ibin = 0; ibin < h->GetNbinsX(); ibin++) {
211  double t0 = h->GetBinCenter(ibin + 1);
212  double chi = h->GetBinContent(ibin + 1);
213  chi += -2 * reco.getLogL(t0, timeMin, timeMax, m_sigmaSmear);
214  h->SetBinContent(ibin + 1, chi);
215  }
216  auto h1 = getObjectPtr<TH2F>("tracks_per_slot");
217  h1->Fill(trk.getModuleID(), sub);
218 
219  // fill histograms of hits
220  auto h2 = getObjectPtr<TH1F>("numHits");
221  auto h3 = getObjectPtr<TH2F>("timeHits");
222  for (const auto& digit : m_digits) {
223  if (digit.getHitQuality() != TOPDigit::c_Good) continue;
224  if (digit.getModuleID() != trk.getModuleID()) continue;
225  if (digit.getTime() < timeMin) continue;
226  if (digit.getTime() > timeMax) continue;
227  h2->Fill(digit.getModuleID());
228  int bs = digit.getBoardstackNumber();
229  h3->Fill((digit.getModuleID() * 4 + bs - 1.5) / 4.0 , digit.getTime());
230  }
231  ntra++;
232  }
233 
234  // fill just another control histogram
235 
236  if (ntra > 0) {
237  auto h4 = getObjectPtr<TH1F>("offset");
238  h4->Fill(m_recBunch->getCurrentOffset());
239  }
240 
241  }
242 
244 }
Belle2::TOP::TOPreco::getFlag
int getFlag()
Return status.
Definition: TOPreco.cc:278
Belle2::TOP::TOPreco::reconstruct
void reconstruct(const TOPtrack &trk, int pdg=0)
Run reconstruction for a given track.
Definition: TOPreco.cc:268
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::TOP::TOPtrack
Class to hold reconstructed track, interface to fortran.
Definition: TOPtrack.h:42
Belle2::TOP::TOPtrack::isValid
bool isValid() const
Check if track is properly defined.
Definition: TOPtrack.h:87
Belle2::TOP::TOPreco::addData
int addData(int moduleID, int pixelID, double time, double timeError)
Add data.
Definition: TOPreco.cc:214
Belle2::TOP::TOPtrack::getModuleID
int getModuleID() const
Return module ID.
Definition: TOPtrack.h:217
Belle2::TOP::TOPreco
TOP reconstruction: this class provides interface to fortran code.
Definition: TOPreco.h:36
Belle2::TOP::TOPreco::getLogL
double getLogL(int i=0)
Return log likelihood for i-th mass hypothesis.
Definition: TOPreco.cc:303
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::TOPModuleT0LLCollectorModule
Collector for module T0 calibration with neg.
Definition: TOPModuleT0LLCollectorModule.h:51
Belle2::TOP::TOPconfigure
Configure TOP geometry for reconstruction: provides interface to fortran code.
Definition: TOPconfigure.h:36
Belle2::TOP::TrackSelector
Utility for the track selection - used in various calibration modules.
Definition: TrackSelector.h:37
Belle2::TOP::TOPreco::setPDFoption
void setPDFoption(PDFoption opt, int NP=0, int NC=0)
Sets PDF option.
Definition: TOPreco.cc:196