Belle II Software  release-05-02-19
TOPRingPlotterModule.cc
1 /**************************************************************************
2 * BASF2 (Belle Analysis Framework 2) *
3 * Copyright(C) 2020 - Belle II Collaboration *
4 * *
5 * Author: The Belle II Collaboration *
6 * Contributors: Umberto Tamponi *
7 * *
8 * This software is provided "as is" without any warranty. *
9 **************************************************************************/
10 
11 #include <top/modules/TOPRingPlotter/TOPRingPlotterModule.h>
12 
13 #include <analysis/dataobjects/ParticleList.h>
14 #include <analysis/dataobjects/Particle.h>
15 
16 #include <framework/datastore/StoreArray.h>
17 #include <framework/datastore/StoreObjPtr.h>
18 
19 #include <framework/logging/Logger.h>
20 #include <framework/utilities/MakeROOTCompatible.h>
21 
22 #include <iostream>
23 #include <TTree.h>
24 #include <TH2F.h>
25 #include <TFile.h>
26 #include <TRandom.h>
27 
28 // framework - DataStore
29 #include <framework/datastore/DataStore.h>
30 #include <framework/datastore/StoreArray.h>
31 #include <framework/datastore/StoreObjPtr.h>
32 
33 // framework aux
34 #include <framework/gearbox/Unit.h>
35 #include <framework/gearbox/Const.h>
36 #include <framework/logging/Logger.h>
37 
38 // DataStore classes
39 #include <framework/dataobjects/EventMetaData.h>
40 #include <mdst/dataobjects/Track.h>
41 #include <mdst/dataobjects/MCParticle.h>
42 #include <top/dataobjects/TOPLikelihood.h>
43 #include <top/dataobjects/TOPBarHit.h>
44 #include <tracking/dataobjects/ExtHit.h>
45 #include <top/dataobjects/TOPDigit.h>
46 #include <mdst/dataobjects/ECLCluster.h>
47 
48 #include <analysis/VariableManager/Manager.h>
49 
50 #include <top/geometry/TOPGeometryPar.h>
51 #include <top/reconstruction/TOPreco.h>
52 #include <top/reconstruction/TOPtrack.h>
53 #include <top/reconstruction/TOPconfigure.h>
54 
55 using namespace Belle2;
56 using namespace TOP;
57 using namespace Belle2::Variable;
58 
59 
60 REG_MODULE(TOPRingPlotter)
61 
63 {
64  // Set module properties
65  setDescription(
66  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");
67 
68  // Parameter definitions
69  addParam("particleList", m_particleList, "List of particles to be used for plotting", std::string("pi+:all"));
70  addParam("variables", m_variables, "List of variables to be saved", {});
71  addParam("pdgHyp", m_pdgHyp, "List of pdg codes for which the PDF is sampled ansd saved. Valid values are 11, 13, 211, 321, 2212.",
72  m_pdgHyp);
73  addParam("outputName", m_outputName, "Name of the output file", std::string("TOPRings.root"));
74  addParam("nToys", m_toyNumber,
75  "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) ",
76  1);
77  addParam("saveHistograms", m_saveHistograms,
78  "Set true to save the histograms of the maps in addition to the list of sampled pseudo-digits", false);
79 }
80 
81 
83 {
84  m_pdfSamplesP = 0;
85  m_pdfToysP = 0;
86  m_pdfSamplesK = 0;
87  m_pdfToysK = 0;
88  m_pdfSamplesPI = 0;
89  m_pdfToysPI = 0;
90  m_pdfSamplesMU = 0;
91  m_pdfToysMU = 0;
92  m_pdfSamplesE = 0;
93  m_pdfToysE = 0;
94  m_hitMapMCK->Reset();
95  m_hitMapMCPI->Reset();
96  m_hitMapMCP->Reset();
97  m_hitMapMCE->Reset();
98  m_hitMapMCMU->Reset();
99 
100  for (int i = 0; i < m_MaxPDFPhotons; i++) {
101  m_pdfPixelE[i] = 0;
102  m_pdfTimeE[i] = 0;
103  m_pdfPixelMU[i] = 0;
104  m_pdfTimeMU[i] = 0;
105  m_pdfPixelPI[i] = 0;
106  m_pdfTimePI[i] = 0;
107  m_pdfPixelK[i] = 0;
108  m_pdfTimeK[i] = 0;
109  m_pdfPixelP[i] = 0;
110  m_pdfTimeP[i] = 0;
111  }
112 
113 
114  m_nDigits = 0;
115 
116  for (int i = 0; i < m_MaxPhotons; i++) {
117  m_digitTime[i] = 0;
118  m_digitChannel[i] = 0;
119  m_digitPixel[i] = 0;
120  m_digitPixelCol[i] = 0;
121  m_digitPixelRow[i] = 0;
122  m_digitAmplitude[i] = 0;
123  m_digitWidth[i] = 0;
124  m_digitASICChannel[i] = 0;
125  m_digitPMTNumber[i] = 0;
126  m_digitSlot[i] = {0};
127  m_digitBoardstack[i] = {0};
128  m_digitCarrier[i] = {0};
129  m_digitAsic[i] = {0};
130  m_digitQuality[i] = {0};
131  }
132 
133 
134  return;
135 }
136 
137 
138 void TOPRingPlotterModule::fillPDF(Belle2::Const::ChargedStable ch, const Track* track, TH2F* histo, short* pixelArray,
139  float* timeArray, int& arrayGuard, int& toyCounter)
140 {
141  double mass = ch.getMass();
142  int pdgCode = ch.getPDGCode();
143 
144  TOPtrack trkPDF(track, ch);
145 
146  if (!trkPDF.isValid())
147  return;
148 
149  TOPreco reco(1, &mass, &pdgCode);
150  reco.setPDFoption(TOPreco::c_Fine);
151  reco.setTimeWindow(0, 100);
152  reco.setMass(mass, pdgCode);
153  reco.reconstruct(trkPDF, pdgCode);
154 
155  if (reco.getFlag() != 1)
156  return ; // track is not in the acceptance of TOP
157 
158  arrayGuard = 0; // this keeps track of the total number of PDF samples we have saved so far
159  m_pdfAsHisto->Reset(); // clan the histogram that will contain the PDf to be sampled
160 
161  // loops of pixels and PDF peaks to make a full sampling of the PDF
162  for (int pixelID = 1; pixelID <= 512; pixelID++) {
163  for (int iBinY = 1; iBinY < 500 + 1; iBinY++) {
164  // MARKO'S Trick HERE
165  auto time = m_pdfAsHisto->GetYaxis()->GetBinCenter(iBinY);
166  auto pdfVal = reco.getPDF(pixelID, time, mass, pdgCode); // jitter = 0 understood
167  m_pdfAsHisto->Fill(pixelID - 0.5, time, pdfVal * 0.2); // bins are 1 px X 0.2 ns wide
168  short pixelCol = (pixelID - 1) % 64 + 1;
169  histo->Fill(pixelCol, time, pdfVal * 1.6); // bins in the map are 8 px X 0.2 ns wide
170  }
171  }
172 
173  auto nExpPhot = reco.getExpectedPhotons();
174 
175  m_pdfAsHisto->Scale(nExpPhot);
176  histo->Scale(nExpPhot);
177 
178  // loop over the toys
179  for (int iSim = 0; iSim < m_toyNumber; iSim++) {
180  // tmp arrays. Before appending the PDF-based digits to the list, we want to be sure
181  // that their total number does not exceeds teh size of pixelArray and timeArray.
182  // We run the ful toy for an event storing the results here, and we copy them ro timeArray and pixelArray only
183  // if there's enough room.
184  std::vector<float> tmpTime;
185  std::vector<short> tmpPixel;
186 
187  short numPhot = gRandom->Poisson(nExpPhot);
188 
189  for (short rn = 0; rn < numPhot; rn++) {
190  double time = 0;
191  double pxID = 0;
192  m_pdfAsHisto->GetRandom2(pxID, time);
193  tmpTime.push_back(time);
194  tmpPixel.push_back(short(pxID));
195  }
196 
197 
198  // before transferring the results of the PDF sampling into pixelArray and timeArray, check if there's
199  // enough room.
200  if (arrayGuard + tmpTime.size() < m_MaxPDFPhotons) {
201  for (unsigned int iPh = 0; iPh < tmpTime.size(); iPh++) {
202  pixelArray[arrayGuard] = tmpPixel[iPh];
203  timeArray[arrayGuard] = tmpTime[iPh];
204  arrayGuard++;
205  }
206  toyCounter++;
207  } else {
208  return; // we already reached the max size of the arrays for sampling, no point in keep doing this
209  }
210  }
211 
212  return;
213 }
214 
215 
216 
217 
218 
220 {
221 
222  // Check the list of pdg hypotheses
223  for (auto pdg : m_pdgHyp) {
224  if ((pdg != 11) & (pdg != 211) & (pdg != 321) & (pdg != 13) & (pdg != 2212))
225  B2FATAL("Invalid PDG hypothesis for the PDF evaluation: " << pdg);
226  short duplicateCount = 0;
227  for (auto pdg2 : m_pdgHyp) {
228  if (pdg == pdg2)
229  duplicateCount++;
230  }
231  if (duplicateCount > 1)
232  B2FATAL("Duplicate PDG hypothesis found");
233  }
234 
235  m_hitMapMCK = new TH2F("hitMapMCK", "x-t map of the kaon PDF", 64, 0, 64, 500, 0, 100);
236  m_hitMapMCPI = new TH2F("hitMapMCPI", "x-t map of the pion PDF", 64, 0, 64, 500, 0, 100);
237  m_hitMapMCP = new TH2F("hitMapMCP", "x-t map of the proton PDF", 64, 0, 64, 500, 0, 100);
238  m_hitMapMCE = new TH2F("hitMapMCE", "x-t map of the electron PDF", 64, 0, 64, 500, 0, 100);
239  m_hitMapMCMU = new TH2F("hitMapMCMU", "x-t map of the muon PDF", 64, 0, 64, 500, 0, 100);
240 
241  m_pdfAsHisto = new TH2F("pdfAsHisto", "tms holder for the pdf to be sampled", 512, 0, 512, 500, 0, 100);
242 
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 == 11)
250  m_tree->Branch("hitMapMCE", "TH2F", &m_hitMapMCE, 32000, 0);
251  else if (pdg == 13)
252  m_tree->Branch("hitMapMCMU", "TH2F", &m_hitMapMCMU, 32000, 0);
253  else if (pdg == 211)
254  m_tree->Branch("hitMapMCPI", "TH2F", &m_hitMapMCPI, 32000, 0);
255  else if (pdg == 321)
256  m_tree->Branch("hitMapMCK", "TH2F", &m_hitMapMCK, 32000, 0);
257  else if (pdg == 2212)
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  // This parts parses the list of variables and creates the branches in the output tree.
304  // This is copy-pasted from VariablesToNtupleModule
305  m_variables = Variable::Manager::Instance().resolveCollections(m_variables);
306  m_branchAddresses.resize(m_variables.size() + 1);
307  size_t enumerate = 0;
308  for (const std::string& varStr : m_variables) {
309  std::string branchName = makeROOTCompatible(varStr);
310  m_tree->Branch(branchName.c_str(), &m_branchAddresses[enumerate], (branchName + "/D").c_str());
311 
313  if (!var) {
314  B2ERROR("Variable '" << varStr << "' is not available in Variable::Manager!");
315  } else {
316  if (m_particleList.empty() && var->description.find("[Eventbased]") == std::string::npos) {
317  B2ERROR("Variable '" << varStr << "' is not an event-based variable, "
318  "but you are using VariablesToNtuple without a decay string, i.e. in the event-wise mode.\n"
319  "If you have created an event-based alias you can wrap your alias with `eventCached` to "
320  "declare it as event based, which avoids this error.\n\n"
321  "vm.addAlias('myAliasName', 'eventCached(myAlias)')");
322  continue;
323  }
324  m_functions.push_back(var->function);
325  }
326  enumerate++;
327  }
328 
329  StoreArray<Track> tracks;
330  tracks.isRequired();
331  StoreArray<ExtHit> extHits;
332  extHits.isRequired();
333  StoreArray<TOPLikelihood> likelihoods;
334  likelihoods.isRequired();
335  StoreArray<TOPDigit> digits;
336  digits.isRequired();
337  StoreArray<TOPBarHit> barHits;
338  barHits.isOptional();
339 
340  TOPconfigure config;
341 }
342 
343 
344 
345 
347 {
348 
349  StoreObjPtr<EventMetaData> evtMetaData;
350  StoreArray<Track> tracks;
351  StoreArray<TOPDigit> digits;
352  StoreObjPtr<ParticleList> particles(m_particleList);
353 
354  // skip on invalid lists
355  if (! particles.isValid())
356  return;
357 
358  // loops over all particles
359  int n_part = particles->getListSize();
360  for (int i = 0; i < n_part; i++) {
361 
362  // Reset the tree variables to the default values
363  resetTree();
364 
365  // Get the track associated to the particle
366  const Particle* particle = particles->getParticle(i);
367  const Track* track = particle->getTrack();
368 
369  bool isFromTrack = true;
370  short moduleID = 0;
371  short deltaModule = 0;
372 
373  // If there's no track, we will use the cluster position
374  // to pin-point the moduleID. This is meant to allow to use
375  // the module to study ECL backsplashes and conversions in the
376  // TOP volume
377  if (!track) {
378  const ECLCluster* clstr = particle->getECLCluster();
379  if (! clstr)
380  continue; // if there's not even a cluster, we can't do much
381  else {
382  if (clstr->getDetectorRegion() != 2)
383  continue; // the cluster is not in the barrel
384  isFromTrack = false;
385  // find the two slots closer to the cluster
386  auto phi = particle->getECLCluster()->getPhi();
387  if (phi < 0)
388  phi = 2 * TMath::Pi() - phi;
389  auto rawModuleID = phi / (TMath::Pi() / 8) + 1;
390  moduleID = int(rawModuleID);
391  if (rawModuleID - moduleID > 0.5) {
392  deltaModule = 1;
393  if (moduleID == 16)
394  deltaModule = -15; // ad-hoc treatment for the 16-1 boundary
395  } else {
396  deltaModule = -1;
397  if (moduleID == 1)
398  deltaModule = 15; // ad-hoc treatment for the 1-16 boundary
399  }
400  }
401  } else {
402  // Check if one can make a TOP track
403  TOPtrack trk(track, Const::pion);
404  if (!trk.isValid())
405  continue;
406  moduleID = trk.getModuleID();
407  }
408 
409  if (moduleID < 1)
410  continue;
411 
412  // Save the digit info for all the topDigits in the module where the track was extrapolated
413  for (const auto& digi : digits) {
414  // if we have a cluster, save the digits from the two modules closer to the cluster
415  if ((digi.getModuleID() == moduleID) || (!isFromTrack && digi.getModuleID() == moduleID + deltaModule)) {
416  m_digitTime[m_nDigits] = digi.getTime();
417  m_digitChannel[m_nDigits] = digi.getChannel();
418  m_digitPixel[m_nDigits] = digi.getPixelID();
419  m_digitPixelCol[m_nDigits] = digi.getPixelCol();
420  m_digitPixelRow[m_nDigits] = digi.getPixelRow();
421  m_digitAmplitude[m_nDigits] = digi.getPulseHeight();
422  m_digitWidth[m_nDigits] = digi.getPulseWidth();
423  m_digitASICChannel[m_nDigits] = digi.getASICChannel();
424  m_digitPMTNumber[m_nDigits] = digi.getPMTNumber();
425  m_digitSlot[m_nDigits] = digi.getModuleID();
426  m_digitBoardstack[m_nDigits] = digi.getBoardstackNumber();
427  m_digitCarrier[m_nDigits] = digi.getCarrierNumber();
428  m_digitAsic[m_nDigits] = digi.getASICNumber();
429  m_digitQuality[m_nDigits] = digi.getHitQuality();
430  m_nDigits++;
431  if (m_nDigits > m_MaxPhotons) break;
432  }
433  }
434 
435  // Save the track variables from the VM
436  for (unsigned int iVar = 0; iVar < m_variables.size(); iVar++) {
437  m_branchAddresses[iVar] = m_functions[iVar](particle);
438  }
439 
440  // Save the PDF samplings
441  if (isFromTrack) {
442  for (auto pdg : m_pdgHyp) {
443  if (pdg == 11)
444  fillPDF(Belle2::Const::electron, track, m_hitMapMCE, m_pdfPixelE, m_pdfTimeE, m_pdfSamplesE, m_pdfToysE);
445  else if (pdg == 13)
446  fillPDF(Belle2::Const::muon, track, m_hitMapMCMU, m_pdfPixelMU, m_pdfTimeMU, m_pdfSamplesMU, m_pdfToysMU);
447  else if (pdg == 211)
448  fillPDF(Belle2::Const::pion, track, m_hitMapMCPI, m_pdfPixelPI, m_pdfTimePI, m_pdfSamplesPI, m_pdfToysPI);
449  else if (pdg == 321)
450  fillPDF(Belle2::Const::kaon, track, m_hitMapMCK, m_pdfPixelK, m_pdfTimeK, m_pdfSamplesK, m_pdfToysK);
451  else if (pdg == 2212)
452  fillPDF(Belle2::Const::proton, track, m_hitMapMCP, m_pdfPixelP, m_pdfTimeP, m_pdfSamplesP, m_pdfToysP);
453  }
454  }
455 
456  m_tree->Fill();
457  }
458  return;
459 }
460 
461 
463 {
464  m_tree->Write();
465  m_outputFile->Close();
466 }
467 
468 
Belle2::Variable::Manager::Var
A variable returning a floating-point value for a given Particle.
Definition: Manager.h:137
Belle2::TOP::TOPreco::getFlag
int getFlag()
Return status.
Definition: TOPreco.cc:278
Belle2::TOP::TOPreco::getExpectedPhotons
double getExpectedPhotons(int i=0)
Return expected number of photons (signal + background) for i-th mass hypothesis.
Definition: TOPreco.cc:284
Belle2::TOP::TOPreco::reconstruct
void reconstruct(const TOPtrack &trk, int pdg=0)
Run reconstruction for a given track.
Definition: TOPreco.cc:268
Belle2::Variable::Manager::getVariable
const Var * getVariable(std::string name)
Get the variable belonging to the given key.
Definition: Manager.cc:33
Belle2::ECLCluster
ECL cluster data.
Definition: ECLCluster.h:39
Belle2::Const::electron
static const ChargedStable electron
electron particle
Definition: Const.h:533
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::TOPRingPlotterModule::fillPDF
void fillPDF(Belle2::Const::ChargedStable, const Track *, TH2F *, short *, float *, int &, int &)
Fills the pdf-related branches.
Definition: TOPRingPlotterModule.cc:138
Belle2::TOPRingPlotterModule::initialize
void initialize() override
Prepares the tree.
Definition: TOPRingPlotterModule.cc:219
Belle2::TOP::TOPtrack::isValid
bool isValid() const
Check if track is properly defined.
Definition: TOPtrack.h:87
Belle2::Const::ParticleType::getPDGCode
int getPDGCode() const
PDG code.
Definition: Const.h:349
Belle2::TOP::TOPtrack::getModuleID
int getModuleID() const
Return module ID.
Definition: TOPtrack.h:217
Belle2::Const::kaon
static const ChargedStable kaon
charged kaon particle
Definition: Const.h:536
Belle2::Module
Base class for Modules.
Definition: Module.h:74
Belle2::Const::pion
static const ChargedStable pion
charged pion particle
Definition: Const.h:535
Belle2::TOP::TOPreco
TOP reconstruction: this class provides interface to fortran code.
Definition: TOPreco.h:36
Belle2::makeROOTCompatible
std::string makeROOTCompatible(std::string str)
Remove special characters that ROOT dislikes in branch names, e.g.
Definition: MakeROOTCompatible.cc:74
Belle2::TOP::TOPreco::setMass
void setMass(double mass, int pdg)
Set mass of the particle hypothesis (overrides settings in the constructor)
Definition: TOPreco.cc:172
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::StoreObjPtr
Type-safe access to single objects in the data store.
Definition: ParticleList.h:33
Belle2::TOPRingPlotterModule
A module to plot the x-t images from TOP, and in general study the distribution of the digits associa...
Definition: TOPRingPlotterModule.h:33
Belle2::TOP::TOPconfigure
Configure TOP geometry for reconstruction: provides interface to fortran code.
Definition: TOPconfigure.h:36
Belle2::TOPRingPlotterModule::terminate
void terminate() override
Writes the tree.
Definition: TOPRingPlotterModule.cc:462
Belle2::TOPRingPlotterModule::event
void event() override
Fills the tree.
Definition: TOPRingPlotterModule.cc:346
Belle2::Variable::Manager::resolveCollections
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:139
Belle2::Particle
Class to store reconstructed particles.
Definition: Particle.h:77
Belle2::Const::proton
static const ChargedStable proton
proton particle
Definition: Const.h:537
Belle2::rn
TString rn()
Get random string.
Definition: tools.h:41
Belle2::Const::ChargedStable
Provides a type-safe way to pass members of the chargedStableSet set.
Definition: Const.h:465
Belle2::Const::muon
static const ChargedStable muon
muon particle
Definition: Const.h:534
Belle2::Track
Class that bundles various TrackFitResults.
Definition: Track.h:35
Belle2::StoreArray
Accessor to arrays stored in the data store.
Definition: ECLMatchingPerformanceExpertModule.h:33
Belle2::ECLCluster::getDetectorRegion
int getDetectorRegion() const
Return detector region: 0: below acceptance, 1: FWD, 2: BRL, 3: BWD, 11: FWDGAP, 13: BWDGAP.
Definition: ECLCluster.cc:68
Belle2::TOP::TOPreco::getPDF
double getPDF(int pixelID, double t, double mass, int PDG, double jitter=0)
Return PDF for pixel pixelID at time t for mass hypothesis mass.
Definition: TOPreco.cc:428
Belle2::TOP::TOPreco::setPDFoption
void setPDFoption(PDFoption opt, int NP=0, int NC=0)
Sets PDF option.
Definition: TOPreco.cc:196
Belle2::TOP::TOPreco::setTimeWindow
void setTimeWindow(double Tmin, double Tmax)
Set time window for photons.
Definition: TOPreco.cc:180
Belle2::TOPRingPlotterModule::resetTree
void resetTree()
reset the tree variables.
Definition: TOPRingPlotterModule.cc:82
Belle2::Variable::Manager::Instance
static Manager & Instance()
get singleton instance.
Definition: Manager.cc:27
Belle2::Const::ParticleType::getMass
double getMass() const
Particle mass.
Definition: UnitConst.cc:310