Belle II Software  release-08-01-10
TOPValidationCollectorModule.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/TOPValidationCollectorModule.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 #include <framework/dataobjects/EventMetaData.h>
19 
20 // root
21 #include <TRandom.h>
22 #include <TH1F.h>
23 #include <TH2F.h>
24 #include <TTree.h>
25 
26 #include <limits>
27 
28 using namespace std;
29 
30 namespace Belle2 {
36  using namespace TOP;
37 
38  //-----------------------------------------------------------------
40  //-----------------------------------------------------------------
41 
42  REG_MODULE(TOPValidationCollector);
43 
44  //-----------------------------------------------------------------
45  // Implementation
46  //-----------------------------------------------------------------
47 
48  TOPValidationCollectorModule::TOPValidationCollectorModule()
49  {
50  // set module description and processing properties
51  setDescription("A collector for automatic validation of TOP calibration");
52  setPropertyFlags(c_ParallelProcessingCertified);
53 
54  // module parameters
55  addParam("numBins", m_numBins, "number of bins of the search region", 100);
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 
76  void TOPValidationCollectorModule::prepare()
77  {
78  // input collections
79 
80  m_digits.isRequired();
81  m_tracks.isRequired();
82  m_extHits.isRequired();
83  m_recBunch.isRequired();
84 
85  // Parse PDF option
86 
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("Unknown PDF option '" << m_pdfOption << "'");
95  }
96 
97  // set track selector
98 
99  if (m_sample == "dimuon" or m_sample == "bhabha") {
100  m_selector = TrackSelector(m_sample);
101  m_selector.setDeltaEcms(m_deltaEcms);
102  m_selector.setCutOnPOCA(m_dr, m_dz);
103  m_selector.setCutOnLocalZ(m_minZ, m_maxZ);
104  } else {
105  B2ERROR("Invalid sample type '" << m_sample << "'");
106  }
107 
108  // create chi2 minimum finders
109 
110  double tmin = -m_timeRange / 2;
111  double tmax = m_timeRange / 2;
112  for (int set = 0; set < c_numSets; set++) {
113  for (int slot = 1; slot <= c_numModules; slot++) {
114  m_finders[set].push_back(Chi2MinimumFinder1D(m_numBins, tmin, tmax));
115  }
116  }
117 
118  // create and register histograms and tree
119 
120  for (int slot = 1; slot <= c_numModules; slot++) {
121  string slotName = to_string(slot);
122  if (slot < 10) slotName.insert(0, "0");
123  string name = "chi2_slot" + slotName;
124  string title = "Chi2 scan, slot" + slotName + "; channel; t0 [ns]";
125  auto h = new TH2F(name.c_str(), title.c_str(), c_numChannels, 0, c_numChannels, m_numBins, tmin, tmax);
126  registerObject<TH2F>(name, h);
127  m_namesChi.push_back(name);
128  }
129 
130  for (int slot = 1; slot <= c_numModules; slot++) {
131  string slotName = to_string(slot);
132  if (slot < 10) slotName.insert(0, "0");
133  string name = "hits_slot" + slotName;
134  string title = "Photon hits, slot" + slotName + "; channel; time [ns]";
135  auto h = new TH2F(name.c_str(), title.c_str(), c_numChannels, 0, c_numChannels, 100, 0., 20.);
136  registerObject<TH2F>(name, h);
137  m_namesHit.push_back(name);
138  }
139 
140  auto h = new TH1F("moduleT0_pulls", "Module T0 pulls; pulls", 200, -15.0, 15.0);
141  registerObject<TH1F>("moduleT0_pulls", h);
142 
143  auto tree = new TTree("tree", "TOP calibration validation tree");
144  tree->Branch("expNo", &m_treeEntry.expNo, "expNo/I");
145  tree->Branch("runNo", &m_treeEntry.runNo, "runNo/I");
146  tree->Branch("numTracks", &m_treeEntry.numTracks, "numTracks/I");
147  tree->Branch("commonT0", &m_treeEntry.commonT0, "commonT0/F");
148  tree->Branch("commonT0Err", &m_treeEntry.commonT0Err, "commonT0Err/F");
149  tree->Branch("moduleT0", &m_treeEntry.moduleT0, "moduleT0[16]/F");
150  tree->Branch("moduleT0Err", &m_treeEntry.moduleT0Err, "moduleT0Err[16]/F");
151  tree->Branch("numTBCalibrated", &m_treeEntry.numTBCalibrated, "numTBCalibrated[16]/I");
152  tree->Branch("numT0Calibrated", &m_treeEntry.numT0Calibrated, "numT0Calibrated[16]/I");
153  tree->Branch("numActive", &m_treeEntry.numActive, "numActive[16]/I");
154  tree->Branch("numActiveCalibrated", &m_treeEntry.numActiveCalibrated, "numActiveCalibrated[16]/I");
155  tree->Branch("thrEffi", &m_treeEntry.thrEffi, "thrEffi[16]/F");
156  tree->Branch("asicShifts", &m_treeEntry.asicShifts, "asicShifts[4]/F");
157  tree->Branch("svdOffset", &m_treeEntry.svdOffset, "svdOffset/F");
158  tree->Branch("svdSigma", &m_treeEntry.svdSigma, "svdSigma/F");
159  tree->Branch("cdcOffset", &m_treeEntry.cdcOffset, "cdcOffset/F");
160  tree->Branch("cdcSigma", &m_treeEntry.cdcSigma, "cdcSigma/F");
161  tree->Branch("fillPatternOffset", &m_treeEntry.fillPatternOffset, "fillPatternOffset/F");
162  tree->Branch("fillPatternFraction", &m_treeEntry.fillPatternFraction, "fillPatternFraction/F");
163  registerObject<TTree>("tree", tree);
164  }
165 
166 
167  void TOPValidationCollectorModule::startRun()
168  {
169  // initialize tree variables
170 
171  m_treeEntry.clear();
172  StoreObjPtr<EventMetaData> evtMetaData;
173  m_treeEntry.expNo = evtMetaData->getExperiment();
174  m_treeEntry.runNo = evtMetaData->getRun();
175 
176  // clear minimum finders
177 
178  for (auto& finders : m_finders) {
179  for (auto& finder : finders) finder.clear();
180  }
181 
182  // pass payload summaries to tree
183 
184  const auto& fe_mapper = TOPGeometryPar::Instance()->getFrontEndMapper();
185  for (unsigned module = 0; module < c_numModules; module++) {
186  auto& numTBCalibrated = m_treeEntry.numTBCalibrated[module];
187  auto& numT0Calibrated = m_treeEntry.numT0Calibrated[module];
188  auto& numActive = m_treeEntry.numActive[module];
189  auto& numActiveCalibrated = m_treeEntry.numActiveCalibrated[module];
190  auto& thrEffi = m_treeEntry.thrEffi[module];
191  int slot = module + 1;
192  for (unsigned channel = 0; channel < c_numChannels; channel++) {
193  bool tbCalibrated = false;
194  const auto* fe = fe_mapper.getMap(slot, channel / 128);
195  if (fe) {
196  tbCalibrated = m_timebase->isAvailable(fe->getScrodID(), channel);
197  } else {
198  B2ERROR("No front-end map found");
199  }
200  bool t0Calibrated = m_channelT0->isCalibrated(slot, channel);
201  bool active = m_channelMask->isActive(slot, channel);
202  if (tbCalibrated) numTBCalibrated++;
203  if (t0Calibrated) numT0Calibrated++;
204  if (active) numActive++;
205  if (tbCalibrated and t0Calibrated and active) {
206  numActiveCalibrated++;
207  thrEffi += m_thresholdEff->getThrEff(slot, channel);
208  }
209  }
210  if (numActiveCalibrated > 0) thrEffi /= numActiveCalibrated;
211  }
212 
213  for (unsigned carrier = 0; carrier < 4; carrier++) {
214  unsigned asic = (3 * 4 + carrier) * 4;
215  m_treeEntry.asicShifts[carrier] = m_asicShift->isCalibrated(13, asic) ? m_asicShift->getT0(13, asic) :
216  std::numeric_limits<float>::quiet_NaN();
217  }
218 
219  const auto& svd = m_eventT0Offset->get(Const::SVD);
220  m_treeEntry.svdOffset = svd.offset;
221  m_treeEntry.svdSigma = svd.sigma;
222 
223  const auto& cdc = m_eventT0Offset->get(Const::CDC);
224  m_treeEntry.cdcOffset = cdc.offset;
225  m_treeEntry.cdcSigma = cdc.sigma;
226 
227  m_treeEntry.fillPatternOffset = m_fillPatternOffset->isCalibrated() ? m_fillPatternOffset->get() :
228  std::numeric_limits<float>::quiet_NaN();
229  m_treeEntry.fillPatternFraction = m_fillPatternOffset->getFraction();
230  }
231 
232 
233  void TOPValidationCollectorModule::collect()
234  {
235  // bunch must be reconstructed
236 
237  if (not m_recBunch.isValid()) return;
238  if (not m_recBunch->isReconstructed()) return;
239 
240  TOPRecoManager::setDefaultTimeWindow();
241  const auto& chMapper = TOPGeometryPar::Instance()->getChannelMapper();
242 
243  // loop over reconstructed tracks, make a selection and accumulate log likelihoods
244 
245  for (const auto& track : m_tracks) {
246 
247  // track selection
248  TOPTrack trk(track);
249  if (not trk.isValid()) continue;
250 
251  if (not m_selector.isSelected(trk)) continue;
252 
253  // construct PDF
254  PDFConstructor pdfConstructor(trk, m_selector.getChargedStable(), m_PDFOption);
255  if (not pdfConstructor.isValid()) continue;
256 
257  // minimization procedure: accumulate
258  unsigned module = trk.getModuleID() - 1;
259  if (module >= m_namesChi.size()) continue;
260  auto h = getObjectPtr<TH2F>(m_namesChi[module]);
261  int set = gRandom->Integer(c_numSets); // generate sub-sample number
262  auto& finder = m_finders[set][module];
263  for (int ibin = 0; ibin < h->GetNbinsY(); ibin++) {
264  double t0 = h->GetYaxis()->GetBinCenter(ibin + 1);
265  const auto& pixelLogLs = pdfConstructor.getPixelLogLs(t0, m_sigmaSmear);
266  for (unsigned channel = 0; channel < c_numChannels; channel++) {
267  int pix = chMapper.getPixelID(channel) - 1;
268  double chi = h->GetBinContent(channel + 1, ibin + 1);
269  chi += -2 * pixelLogLs[pix].logL;
270  h->SetBinContent(channel + 1, ibin + 1, chi);
271  }
272  double logL = 0;
273  for (auto& LL : pixelLogLs) logL += LL.logL;
274  finder.add(ibin, -2 * logL);
275  }
276  m_treeEntry.numTracks++;
277 
278  auto h1 = getObjectPtr<TH2F>(m_namesHit[module]);
279  for (const auto& digit : m_digits) {
280  if (digit.getHitQuality() != TOPDigit::c_Good) continue;
281  if (digit.getModuleID() != trk.getModuleID()) continue;
282  h1->Fill(digit.getChannel(), digit.getTime());
283  }
284 
285  }
286  }
287 
288  void TOPValidationCollectorModule::closeRun()
289  {
290 
291  // module T0 pulls
292 
293  for (int module = 0; module < c_numModules; module++) {
294  std::vector<double> pos, err;
295  for (int set = 0; set < c_numSets; set++) {
296  const auto& minimum = m_finders[set][module].getMinimum();
297  if (minimum.valid) {
298  pos.push_back(minimum.position);
299  err.push_back(minimum.error);
300  }
301  }
302  auto h_pulls = getObjectPtr<TH1F>("moduleT0_pulls");
303  for (unsigned i = 0; i < pos.size(); i++) {
304  for (unsigned j = i + 1; j < pos.size(); j++) {
305  double pull = (pos[i] - pos[j]) / sqrt(err[i] * err[i] + err[j] * err[j]);
306  h_pulls->Fill(pull);
307  }
308  }
309  }
310 
311  // module T0 residuals
312 
313  for (int module = 0; module < c_numModules; module++) {
314  auto& finder = m_finders[0][module];
315  for (int set = 1; set < c_numSets; set++) {
316  finder.add(m_finders[set][module]);
317  }
318  const auto& minimum = finder.getMinimum();
319  if (minimum.valid) {
320  m_treeEntry.moduleT0[module] = minimum.position;
321  m_treeEntry.moduleT0Err[module] = minimum.error;
322  }
323  }
324 
325  // common T0 residual
326 
327  auto& finder = m_finders[0][0];
328  for (int module = 1; module < c_numModules; module++) {
329  finder.add(m_finders[0][module]);
330  }
331  const auto& minimum = finder.getMinimum();
332  if (minimum.valid) {
333  m_treeEntry.commonT0 = minimum.position;
334  m_treeEntry.commonT0Err = minimum.error;
335  }
336 
337  // fill the tree
338 
339  auto tree = getObjectPtr<TTree>("tree");
340  tree->Fill();
341  }
342 
344 }
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
Minimum finder using tabulated chi^2 values in one dimension.
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.
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
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
Abstract base class for different kinds of events.