Belle II Software  release-06-01-15
TOPRingPlotterModule.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/TOPRingPlotter/TOPRingPlotterModule.h>
10 
11 #include <analysis/dataobjects/ParticleList.h>
12 #include <analysis/dataobjects/Particle.h>
13 
14 #include <framework/datastore/StoreArray.h>
15 #include <framework/datastore/StoreObjPtr.h>
16 
17 #include <framework/logging/Logger.h>
18 #include <framework/utilities/MakeROOTCompatible.h>
19 
20 #include <iostream>
21 #include <TTree.h>
22 #include <TH2F.h>
23 #include <TFile.h>
24 #include <TRandom.h>
25 #include <TDirectory.h>
26 
27 // framework - DataStore
28 #include <framework/datastore/DataStore.h>
29 #include <framework/datastore/StoreArray.h>
30 #include <framework/datastore/StoreObjPtr.h>
31 
32 // framework aux
33 #include <framework/gearbox/Unit.h>
34 #include <framework/gearbox/Const.h>
35 #include <framework/logging/Logger.h>
36 
37 // DataStore classes
38 #include <mdst/dataobjects/Track.h>
39 #include <mdst/dataobjects/MCParticle.h>
40 #include <top/dataobjects/TOPLikelihoodScanResult.h>
41 #include <top/dataobjects/TOPLikelihood.h>
42 #include <top/dataobjects/TOPBarHit.h>
43 #include <tracking/dataobjects/ExtHit.h>
44 #include <top/dataobjects/TOPDigit.h>
45 #include <mdst/dataobjects/ECLCluster.h>
46 
47 #include <analysis/VariableManager/Manager.h>
48 
49 #include <top/geometry/TOPGeometryPar.h>
50 #include <top/reconstruction_cpp/TOPTrack.h>
51 #include <top/reconstruction_cpp/PDFConstructor.h>
52 #include <top/reconstruction_cpp/TOPRecoManager.h>
53 
54 using namespace Belle2;
55 using namespace TOP;
56 using namespace Belle2::Variable;
57 
58 
59 REG_MODULE(TOPRingPlotter)
60 
62 {
63  // Set module properties
64  setDescription(
65  R"DOC(A module to plot the x-t images from TOP, and in general study the distribution of the digits associated to the particles in a particleList)DOC");
66 
67  // Parameter definitions
68  addParam("particleList", m_particleList, "List of particles to be used for plotting", std::string("pi+:all"));
69  addParam("variables", m_variables, "List of variables to be saved", {});
70  addParam("pdgHyp", m_pdgHyp, "List of pdg codes for which the PDF is sampled ansd saved. Valid values are 11, 13, 211, 321, 2212.",
71  m_pdgHyp);
72  addParam("outputName", m_outputName, "Name of the output file", std::string("TOPRings.root"));
73  addParam("nToys", m_toyNumber,
74  "Number of toys used to sample the PDFs. Set it to a large number for a precise sampling of the PDF in the histograms (suggested for small number of events) ",
75  1);
76  addParam("saveHistograms", m_saveHistograms,
77  "Set true to save the histograms of the maps in addition to the list of sampled pseudo-digits", false);
78  addParam("saveLLScan", m_saveLLScan,
79  "Set true to save the results of the LL scan. Requires to add TOPLLScannerModule to the path", false);
80 
81 }
82 
83 
85 {
86  m_pdfSamplesP = 0;
87  m_pdfToysP = 0;
88  m_pdfSamplesK = 0;
89  m_pdfToysK = 0;
90  m_pdfSamplesPI = 0;
91  m_pdfToysPI = 0;
92  m_pdfSamplesMU = 0;
93  m_pdfToysMU = 0;
94  m_pdfSamplesE = 0;
95  m_pdfToysE = 0;
96  m_hitMapMCK->Reset();
97  m_hitMapMCPI->Reset();
98  m_hitMapMCP->Reset();
99  m_hitMapMCE->Reset();
100  m_hitMapMCMU->Reset();
101  m_nScanPoints = 0;
102 
103  for (int i = 0; i < 10000; i++) {
104  m_scanMass[i] = 0;
105  m_scanLL[i] = 0;
106  }
107 
108  for (int i = 0; i < m_MaxPDFPhotons; i++) {
109  m_pdfPixelE[i] = 0;
110  m_pdfTimeE[i] = 0;
111  m_pdfPixelMU[i] = 0;
112  m_pdfTimeMU[i] = 0;
113  m_pdfPixelPI[i] = 0;
114  m_pdfTimePI[i] = 0;
115  m_pdfPixelK[i] = 0;
116  m_pdfTimeK[i] = 0;
117  m_pdfPixelP[i] = 0;
118  m_pdfTimeP[i] = 0;
119  }
120 
121 
122  m_nDigits = 0;
123 
124  for (int i = 0; i < m_MaxPhotons; i++) {
125  m_digitTime[i] = 0;
126  m_digitChannel[i] = 0;
127  m_digitPixel[i] = 0;
128  m_digitPixelCol[i] = 0;
129  m_digitPixelRow[i] = 0;
130  m_digitAmplitude[i] = 0;
131  m_digitWidth[i] = 0;
132  m_digitASICChannel[i] = 0;
133  m_digitPMTNumber[i] = 0;
134  m_digitSlot[i] = {0};
135  m_digitBoardstack[i] = {0};
136  m_digitCarrier[i] = {0};
137  m_digitAsic[i] = {0};
138  m_digitQuality[i] = {0};
139  }
140 
141 
142  return;
143 }
144 
145 
146 void TOPRingPlotterModule::fillPDF(Belle2::Const::ChargedStable ch, const Track* track, TH2F* histo, short* pixelArray,
147  float* timeArray, int& arrayGuard, int& toyCounter)
148 {
149 
150  if (not track) return;
151 
152  TOPTrack trkPDF(*track);
153 
154  if (!trkPDF.isValid())
155  return;
156 
157  PDFConstructor pdfConstructor(trkPDF, ch, PDFConstructor::c_Fine);
158  if (not pdfConstructor.isValid()) return;
159 
160  arrayGuard = 0; // this keeps track of the total number of PDF samples we have saved so far
161  m_pdfAsHisto->Reset(); // clan the histogram that will contain the PDf to be sampled
162 
163  // loops of pixels and PDF peaks to make a full sampling of the PDF
164  for (int pixelID = 1; pixelID <= 512; pixelID++) {
165  for (int iBinY = 1; iBinY < 500 + 1; iBinY++) {
166  // MARKO'S Trick HERE
167  auto time = m_pdfAsHisto->GetYaxis()->GetBinCenter(iBinY);
168  auto pdfVal = pdfConstructor.getPDFValue(pixelID, time, 0.070); // average electronic time reso = 70 ps
169  m_pdfAsHisto->Fill(pixelID - 0.5, time, pdfVal * 0.2); // bins are 1 px X 0.2 ns wide
170  short pixelCol = (pixelID - 1) % 64 + 1;
171  histo->Fill(pixelCol, time, pdfVal * 1.6); // bins in the map are 8 px X 0.2 ns wide
172  }
173  }
174 
175  auto nExpPhot = pdfConstructor.getExpectedPhotons();
176 
177  m_pdfAsHisto->Scale(nExpPhot);
178  histo->Scale(nExpPhot);
179 
180  // loop over the toys
181  for (int iSim = 0; iSim < m_toyNumber; iSim++) {
182  // tmp arrays. Before appending the PDF-based digits to the list, we want to be sure
183  // that their total number does not exceeds teh size of pixelArray and timeArray.
184  // We run the ful toy for an event storing the results here, and we copy them ro timeArray and pixelArray only
185  // if there's enough room.
186  std::vector<float> tmpTime;
187  std::vector<short> tmpPixel;
188 
189  short numPhot = gRandom->Poisson(nExpPhot);
190 
191  for (short rn = 0; rn < numPhot; rn++) {
192  double time = 0;
193  double pxID = 0;
194  m_pdfAsHisto->GetRandom2(pxID, time);
195  tmpTime.push_back(time);
196  tmpPixel.push_back(short(pxID));
197  }
198 
199  // before transferring the results of the PDF sampling into pixelArray and timeArray, check if there's
200  // enough room.
201  if (arrayGuard + tmpTime.size() < m_MaxPDFPhotons) {
202  for (unsigned int iPh = 0; iPh < tmpTime.size(); iPh++) {
203  pixelArray[arrayGuard] = tmpPixel[iPh];
204  timeArray[arrayGuard] = tmpTime[iPh];
205  arrayGuard++;
206  }
207  toyCounter++;
208  } else {
209  return; // we already reached the max size of the arrays for sampling, no point in keep doing this
210  }
211  }
212 
213  return;
214 }
215 
216 
218 {
219 
220  // Check the list of pdg hypotheses
221  for (auto pdg : m_pdgHyp) {
222  if ((pdg != Const::electron.getPDGCode()) & (pdg != Const::pion.getPDGCode()) & (pdg != Const::kaon.getPDGCode()) &
223  (pdg != Const::muon.getPDGCode()) & (pdg != Const::proton.getPDGCode()))
224  B2FATAL("Invalid PDG hypothesis for the PDF evaluation: " << pdg);
225  short duplicateCount = 0;
226  for (auto pdg2 : m_pdgHyp) {
227  if (pdg == pdg2)
228  duplicateCount++;
229  }
230  if (duplicateCount > 1)
231  B2FATAL("Duplicate PDG hypothesis found");
232  }
233 
234  m_hitMapMCK = new TH2F("hitMapMCK", "x-t map of the kaon PDF", 64, 0, 64, 500, 0, 100);
235  m_hitMapMCPI = new TH2F("hitMapMCPI", "x-t map of the pion PDF", 64, 0, 64, 500, 0, 100);
236  m_hitMapMCP = new TH2F("hitMapMCP", "x-t map of the proton PDF", 64, 0, 64, 500, 0, 100);
237  m_hitMapMCE = new TH2F("hitMapMCE", "x-t map of the electron PDF", 64, 0, 64, 500, 0, 100);
238  m_hitMapMCMU = new TH2F("hitMapMCMU", "x-t map of the muon PDF", 64, 0, 64, 500, 0, 100);
239 
240  m_pdfAsHisto = new TH2F("pdfAsHisto", "tms holder for the pdf to be sampled", 512, 0, 512, 500, 0, 100);
241 
242  TDirectory::TContext context;
243  B2INFO("creating a tree for TOPRingPlotterModule");
244  m_outputFile = new TFile(m_outputName.c_str(), "recreate");
245  m_tree = new TTree("tree", "tree");
246 
247  if (m_saveHistograms) {
248  for (auto pdg : m_pdgHyp) {
249  if (pdg == Const::electron.getPDGCode())
250  m_tree->Branch("hitMapMCE", "TH2F", &m_hitMapMCE, 32000, 0);
251  else if (pdg == Const::muon.getPDGCode())
252  m_tree->Branch("hitMapMCMU", "TH2F", &m_hitMapMCMU, 32000, 0);
253  else if (pdg == Const::pion.getPDGCode())
254  m_tree->Branch("hitMapMCPI", "TH2F", &m_hitMapMCPI, 32000, 0);
255  else if (pdg == Const::kaon.getPDGCode())
256  m_tree->Branch("hitMapMCK", "TH2F", &m_hitMapMCK, 32000, 0);
257  else if (pdg == Const::proton.getPDGCode())
258  m_tree->Branch("hitMapMCP", "TH2F", &m_hitMapMCP, 32000, 0);
259  }
260  }
261 
262  m_tree->Branch("nDigits", &m_nDigits, "m_nDigits/S");
263  m_tree->Branch("digitTime", m_digitTime, "m_digitTime[m_nDigits]/F");
264  m_tree->Branch("digitAmplitude", m_digitAmplitude, "m_digitAmplitude[m_nDigits]/F");
265  m_tree->Branch("digitWidth", m_digitWidth, "m_digitWidth[m_nDigits]/F");
266  m_tree->Branch("digitChannel", m_digitChannel, "m_digitChannel[m_nDigits]/S");
267  m_tree->Branch("digitPixel", m_digitPixel, "m_digitPixel[m_nDigits]/S");
268  m_tree->Branch("digitPixelCol", m_digitPixelCol, "m_digitPixelCol[m_nDigits]/S");
269  m_tree->Branch("digitPixelRow", m_digitPixelRow, "m_digitPixelRow[m_nDigits]/S");
270  m_tree->Branch("digitASICChannel", m_digitASICChannel, "m_digitASICChannel[m_nDigits]/S");
271  m_tree->Branch("digitPMTNumber", m_digitPMTNumber, "m_digitPMTNumber[m_nDigits]/S");
272  m_tree->Branch("digitSlot", m_digitSlot, "m_digitSlot[m_nDigits]/S");
273  m_tree->Branch("digitBoardstack", m_digitBoardstack, "m_digitBoardstack[m_nDigits]/S");
274  m_tree->Branch("digitCarrier", m_digitCarrier, "m_digitCarrier[m_nDigits]/S");
275  m_tree->Branch("digitAsic", m_digitAsic, "m_digitAsic[m_nDigits]/S");
276  m_tree->Branch("digitQuality", m_digitQuality, "m_digitQuality[m_nDigits]/S");
277 
278  m_tree->Branch("pdfSamplesP", &m_pdfSamplesP, "m_pdfSamplesP/I");
279  m_tree->Branch("pdfToysP", &m_pdfToysP, "m_pdfToysP/I");
280  m_tree->Branch("pdfPixelP", m_pdfPixelP, "m_pdfPixelP[m_pdfSamplesP]/S");
281  m_tree->Branch("pdfTimeP", m_pdfTimeP, "m_pdfTimeP[m_pdfSamplesP]/F");
282 
283  m_tree->Branch("pdfSamplesK", &m_pdfSamplesK, "m_pdfSamplesK/I");
284  m_tree->Branch("pdfToysK", &m_pdfToysK, "m_pdfToysK/I");
285  m_tree->Branch("pdfPixelK", m_pdfPixelK, "m_pdfPixelK[m_pdfSamplesK]/S");
286  m_tree->Branch("pdfTimeK", m_pdfTimeK, "m_pdfTimeK[m_pdfSamplesK]/F");
287 
288  m_tree->Branch("pdfSamplesPI", &m_pdfSamplesPI, "m_pdfSamplesPI/I");
289  m_tree->Branch("pdfToysPI", &m_pdfToysPI, "m_pdfToysPI/I");
290  m_tree->Branch("pdfPixelPI", m_pdfPixelPI, "m_pdfPixelPI[m_pdfSamplesPI]/S");
291  m_tree->Branch("pdfTimePI", m_pdfTimePI, "m_pdfTimePI[m_pdfSamplesPI]/F");
292 
293  m_tree->Branch("pdfSamplesMU", &m_pdfSamplesMU, "m_pdfSamplesMU/I");
294  m_tree->Branch("pdfToysMU", &m_pdfToysMU, "m_pdfToysMU/I");
295  m_tree->Branch("pdfPixelMU", m_pdfPixelMU, "m_pdfPixelMU[m_pdfSamplesMU]/S");
296  m_tree->Branch("pdfTimeMU", m_pdfTimeMU, "m_pdfTimeMU[m_pdfSamplesMU]/F");
297 
298  m_tree->Branch("pdfSamplesE", &m_pdfSamplesE, "m_pdfSamplesE/I");
299  m_tree->Branch("pdfToysE", &m_pdfToysE, "m_pdfToysE/I");
300  m_tree->Branch("pdfPixelE", m_pdfPixelE, "m_pdfPixelE[m_pdfSamplesE]/S");
301  m_tree->Branch("pdfTimeE", m_pdfTimeE, "m_pdfTimeE[m_pdfSamplesE]/F");
302 
303  if (m_saveLLScan) {
304  m_tree->Branch("nScanPoints", &m_nScanPoints, "m_nScanPoints/I");
305  m_tree->Branch("scanMass", m_scanMass, "m_scanMass[m_nScanPoints]/F");
306  m_tree->Branch("scanLL", m_scanLL, "m_scanLL[m_nScanPoints]/F");
307  }
308 
309  // This parts parses the list of variables and creates the branches in the output tree.
310  // This is copy-pasted from VariablesToNtupleModule
311  m_variables = Variable::Manager::Instance().resolveCollections(m_variables);
312  m_branchAddresses.resize(m_variables.size() + 1);
313  size_t enumerate = 0;
314  for (const std::string& varStr : m_variables) {
315  std::string branchName = makeROOTCompatible(varStr);
316  m_tree->Branch(branchName.c_str(), &m_branchAddresses[enumerate], (branchName + "/D").c_str());
317 
319  if (!var) {
320  B2ERROR("Variable '" << varStr << "' is not available in Variable::Manager!");
321  } else {
322  if (m_particleList.empty() && var->description.find("[Eventbased]") == std::string::npos) {
323  B2ERROR("Variable '" << varStr << "' is not an event-based variable, "
324  "but you are using VariablesToNtuple without a decay string, i.e. in the event-wise mode.\n"
325  "If you have created an event-based alias you can wrap your alias with `eventCached` to "
326  "declare it as event based, which avoids this error.\n\n"
327  "vm.addAlias('myAliasName', 'eventCached(myAlias)')");
328  continue;
329  }
330  m_functions.push_back(var->function);
331  }
332  enumerate++;
333  }
334 
335  StoreArray<Track> tracks;
336  tracks.isRequired();
337  StoreArray<ExtHit> extHits;
338  extHits.isRequired();
339  StoreArray<TOPLikelihood> likelihoods;
340  likelihoods.isRequired();
341  StoreArray<TOPDigit> digits;
342  digits.isRequired();
343  StoreArray<TOPBarHit> barHits;
344  barHits.isOptional();
345  StoreArray<TOPLikelihoodScanResult> likelihoodScans;
346  likelihoodScans.isOptional();
347 
348 }
349 
350 
352 {
354 
355  StoreArray<TOPDigit> digits;
356  StoreObjPtr<ParticleList> particles(m_particleList);
357 
358  // skip on invalid lists
359  if (! particles.isValid())
360  return;
361 
362  // loops over all particles
363  int n_part = particles->getListSize();
364  for (int i = 0; i < n_part; i++) {
365 
366  // Reset the tree variables to the default values
367  resetTree();
368 
369  // Get the track associated to the particle
370  const Particle* particle = particles->getParticle(i);
371  const Track* track = particle->getTrack();
372 
373  bool isFromTrack = true;
374  short moduleID = 0;
375  short deltaModule = 0;
376 
377  // If there's no track, we will use the cluster position
378  // to pin-point the moduleID. This is meant to allow to use
379  // the module to study ECL backsplashes and conversions in the
380  // TOP volume
381  if (!track) {
382  const ECLCluster* clstr = particle->getECLCluster();
383  if (! clstr)
384  continue; // if there's not even a cluster, we can't do much
385  else {
386  if (clstr->getDetectorRegion() != 2)
387  continue; // the cluster is not in the barrel
388  isFromTrack = false;
389  // find the two slots closer to the cluster
390  auto phi = particle->getECLCluster()->getPhi();
391  if (phi < 0)
392  phi = 2 * TMath::Pi() - phi;
393  auto rawModuleID = phi / (TMath::Pi() / 8) + 1;
394  moduleID = int(rawModuleID);
395  if (rawModuleID - moduleID > 0.5) {
396  deltaModule = 1;
397  if (moduleID == 16)
398  deltaModule = -15; // ad-hoc treatment for the 16-1 boundary
399  } else {
400  deltaModule = -1;
401  if (moduleID == 1)
402  deltaModule = 15; // ad-hoc treatment for the 1-16 boundary
403  }
404  }
405  } else {
406  // Check if one can make a TOP track
407  TOPTrack trk(*track);
408  if (!trk.isValid())
409  continue;
410  moduleID = trk.getModuleID();
411  }
412 
413  if (moduleID < 1)
414  continue;
415 
416  // Save the digit info for all the topDigits in the module where the track was extrapolated
417  for (const auto& digi : digits) {
418  // if we have a cluster, save the digits from the two modules closer to the cluster
419  if (m_nDigits >= m_MaxPhotons) break;
420  if ((digi.getModuleID() == moduleID) || (!isFromTrack && digi.getModuleID() == moduleID + deltaModule)) {
421  m_digitTime[m_nDigits] = digi.getTime();
422  m_digitChannel[m_nDigits] = digi.getChannel();
423  m_digitPixel[m_nDigits] = digi.getPixelID();
424  m_digitPixelCol[m_nDigits] = digi.getPixelCol();
425  m_digitPixelRow[m_nDigits] = digi.getPixelRow();
426  m_digitAmplitude[m_nDigits] = digi.getPulseHeight();
427  m_digitWidth[m_nDigits] = digi.getPulseWidth();
428  m_digitASICChannel[m_nDigits] = digi.getASICChannel();
429  m_digitPMTNumber[m_nDigits] = digi.getPMTNumber();
430  m_digitSlot[m_nDigits] = digi.getModuleID();
431  m_digitBoardstack[m_nDigits] = digi.getBoardstackNumber();
432  m_digitCarrier[m_nDigits] = digi.getCarrierNumber();
433  m_digitAsic[m_nDigits] = digi.getASICNumber();
434  m_digitQuality[m_nDigits] = digi.getHitQuality();
435  m_nDigits++;
436  }
437  }
438 
439  // Save the track variables from the VM
440  for (unsigned int iVar = 0; iVar < m_variables.size(); iVar++) {
441  m_branchAddresses[iVar] = m_functions[iVar](particle);
442  }
443 
444 
445  // Save the PDF samplings
446  if (isFromTrack) {
447 
448  if (m_saveLLScan) {
449  auto llScan = track->getRelated<TOPLikelihoodScanResult>();
450  auto masses = llScan->getCoarseScanMassPoints();
451  auto lls = llScan->getCoarseScanLLValues();
452  m_nScanPoints = lls.size();
453  for (unsigned int j = 0; j < lls.size(); j++) {
454  m_scanMass[j] = masses[j];
455  m_scanLL[j] = lls[j];
456  }
457  }
458 
459  for (auto pdg : m_pdgHyp) {
460  if (pdg == Const::electron.getPDGCode())
461  fillPDF(Belle2::Const::electron, track, m_hitMapMCE, m_pdfPixelE, m_pdfTimeE, m_pdfSamplesE, m_pdfToysE);
462  else if (pdg == Const::muon.getPDGCode())
463  fillPDF(Belle2::Const::muon, track, m_hitMapMCMU, m_pdfPixelMU, m_pdfTimeMU, m_pdfSamplesMU, m_pdfToysMU);
464  else if (pdg == Const::pion.getPDGCode())
465  fillPDF(Belle2::Const::pion, track, m_hitMapMCPI, m_pdfPixelPI, m_pdfTimePI, m_pdfSamplesPI, m_pdfToysPI);
466  else if (pdg == Const::kaon.getPDGCode())
467  fillPDF(Belle2::Const::kaon, track, m_hitMapMCK, m_pdfPixelK, m_pdfTimeK, m_pdfSamplesK, m_pdfToysK);
468  else if (pdg == Const::proton.getPDGCode())
469  fillPDF(Belle2::Const::proton, track, m_hitMapMCP, m_pdfPixelP, m_pdfTimeP, m_pdfSamplesP, m_pdfToysP);
470  }
471  }
472 
473  m_tree->Fill();
474  }
475  return;
476 }
477 
478 
480 {
481  TDirectory::TContext context;
482  m_outputFile->cd();
483  m_tree->Write();
484  m_outputFile->Close();
485 }
486 
487 
Provides a type-safe way to pass members of the chargedStableSet set.
Definition: Const.h:470
static const ChargedStable muon
muon particle
Definition: Const.h:541
static const ChargedStable pion
charged pion particle
Definition: Const.h:542
static const ChargedStable proton
proton particle
Definition: Const.h:544
static const ChargedStable kaon
charged kaon particle
Definition: Const.h:543
static const ChargedStable electron
electron particle
Definition: Const.h:540
ECL cluster data.
Definition: ECLCluster.h:27
int getDetectorRegion() const
Return detector region: 0: below acceptance, 1: FWD, 2: BRL, 3: BWD, 11: FWDGAP, 13: BWDGAP.
Definition: ECLCluster.cc:66
Base class for Modules.
Definition: Module.h:72
Class to store reconstructed particles.
Definition: Particle.h:74
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
bool isOptional(const std::string &name="")
Tell the DataStore about an optional input.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:95
Class to store the result of the TOP LL scan (output of TOPLLScanner).
std::vector< float > getCoarseScanMassPoints() const
Return a std::vector containing the mass points used for the coarse LL scane.
A module to plot the x-t images from TOP, and in general study the distribution of the digits associa...
void initialize() override
Prepares the tree.
void event() override
Fills the tree.
void fillPDF(Belle2::Const::ChargedStable, const Track *, TH2F *, short *, float *, int &, int &)
Fills the pdf-related branches.
void terminate() override
Writes the tree.
void resetTree()
reset the tree variables.
PDF construction and log likelihood determination for a given track and particle hypothesis.
double getPDFValue(int pixelID, double time, double timeErr, double sigt=0) const
Returns PDF value.
bool isValid() const
Checks the object status.
double getExpectedPhotons() const
Returns the expected number of photons within the default time window.
@ c_Fine
y dependent everywhere
static void setTimeWindow(double minTime, double maxTime)
Sets time window.
Reconstructed track at TOP.
Definition: TOPTrack.h:40
bool isValid() const
Checks if track is successfully constructed.
Definition: TOPTrack.h:139
int getModuleID() const
Returns slot ID.
Definition: TOPTrack.h:145
Class that bundles various TrackFitResults.
Definition: Track.h:25
std::vector< std::string > resolveCollections(const std::vector< std::string > &variables)
Resolve Collection Returns variable names corresponding to the given collection or if it is not a col...
Definition: Manager.cc:137
const Var * getVariable(std::string name)
Get the variable belonging to the given key.
Definition: Manager.cc:31
static Manager & Instance()
get singleton instance.
Definition: Manager.cc:25
std::string makeROOTCompatible(std::string str)
Remove special characters that ROOT dislikes in branch names, e.g.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
TString rn()
Get random string.
Definition: tools.h:38
Abstract base class for different kinds of events.
A variable returning a floating-point value for a given Particle.
Definition: Manager.h:133