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