Belle II Software  release-05-02-19
ECLChargedPIDDataAnalysisValidationModule.cc
1 #include <ecl/modules/eclChargedPIDDataAnalysisExpert/ECLChargedPIDDataAnalysisValidationModule.h>
2 
3 #include <ecl/dataobjects/ECLPidLikelihood.h>
4 #include <framework/datastore/StoreArray.h>
5 #include <mdst/dataobjects/ECLCluster.h>
6 #include <mdst/dataobjects/MCParticle.h>
7 #include <mdst/dataobjects/Track.h>
8 
9 #include <TEfficiency.h>
10 
11 using namespace Belle2;
12 
13 //-----------------------------------------------------------------
14 // Register the Module
15 //-----------------------------------------------------------------
16 
17 REG_MODULE(ECLChargedPIDDataAnalysisValidation)
18 
19 //-----------------------------------------------------------------
20 // Implementation
21 //-----------------------------------------------------------------
22 
24 {
25  // Set module properties
26  setDescription("This module dumps a set of histograms with ECL charged PID-related info used for validation, starting from an input file w/ particle-gun-generated charged stable particles (and antiparticles).");
27 
28  // Default charged stable pdgIds (particles & antiparticles)
29  std::vector<int> defaultChargedPdgIds;
30  for (const auto& hypo : Const::chargedStableSet) {
31  defaultChargedPdgIds.push_back(hypo.getPDGCode());
32  defaultChargedPdgIds.push_back(-hypo.getPDGCode());
33  }
34 
35  addParam("inputPdgIdList", m_inputPdgIdList,
36  "The list of (signed) pdgIds of the charged stable particles for which validation plots should be produced. Default is ALL charged stable particles.",
37  defaultChargedPdgIds);
38  addParam("mergeChargeOfPdgIds", m_mergeChargeOfPdgIds,
39  "The list of (unsigned) pdgIds of the charged stable particles for which particle and antiparticle should be merged together in the plots. Default is no merging, meaning separate plots are generated for +/- charged particles for each input pdgId.",
40  std::vector<unsigned int>());
41  addParam("outputFileName", m_outputFileName,
42  "The base name of the output file. The pdgId of the charged particle is appended to the name.",
43  std::string("ECLChargedPid"));
44  addParam("saveValidationTree", m_saveValidationTree,
45  "If this flag is set to True, save also the validation TTree. Default is False.",
46  bool(false));
47 }
48 
50 {
51 }
52 
53 
55 {
56  B2INFO("Initialising ROOT objects...");
57 
58  // Convert pdgId list to a set to remove any accidental repetitions.
59  m_inputPdgIdSet = std::set<int>(m_inputPdgIdList.begin(), m_inputPdgIdList.end());
60 
61  // By default, do not merge particles and antiparticles together,
62  // unless a particle hypo is in the configured "merge" list.
63  for (const auto& hypo : Const::chargedStableSet) {
64  bool merge = (std::find(m_mergeChargeOfPdgIds.begin(), m_mergeChargeOfPdgIds.end(),
65  hypo.getPDGCode()) == m_mergeChargeOfPdgIds.end()) ? false : true;
66  if (merge) {
67  B2WARNING("For (unsigned) hypothesis " << hypo.getPDGCode() << ", validation plots will be merged for +/- charged particles.");
68  }
69  m_mergeChargeFlagByHypo.insert(std::pair<Const::ChargedStable, bool>(hypo, merge));
70  }
71 
72  std::string chargedPdgIdStr;
73  std::string fname;
74 
75  for (const auto& chargedPdgId : m_inputPdgIdSet) {
76 
77  // Check if this pdgId is that of a legit Const::ChargedStable particle.
78  if (!isValidChargedPdg(std::abs(chargedPdgId))) {
79  B2FATAL("PDG: " << chargedPdgId << " in m_inputPdgIdSet is not that of a valid particle in Const::chargedStableSet! Aborting...");
80  }
81 
82  const auto chargedHypo = Const::chargedStableSet.find(std::abs(chargedPdgId));
83 
84  // If merging particles and antiparticles for this hypo, no need to loop twice:
85  // register one TTree for the '+' charged pdgId only.
86  if (m_mergeChargeFlagByHypo[chargedHypo] and chargedPdgId < 0) continue;
87 
88  // Get the idx of this pdgId in the Const::chargedStableSet
89  auto chargedIdx = chargedHypo.getIndex();
90 
91  if (chargedPdgId > 0) {
92  chargedPdgIdStr = std::to_string(chargedPdgId);
93  } else {
94  chargedPdgIdStr = "anti" + std::to_string(std::abs(chargedPdgId));
95  // Add offset to idx.
96  chargedIdx += Const::ChargedStable::c_SetSize;
97  }
98 
99  fname = m_outputFileName + "_" + chargedPdgIdStr + ".root";
100 
101  m_outputFile[chargedIdx] = new TFile(fname.c_str(), "RECREATE");
102 
103  m_tree[chargedIdx] = new TTree("ECLChargedPid", "ECLChargedPid");
104  m_tree[chargedIdx]->Branch("p", &m_p[chargedIdx], "p/F");
105  m_tree[chargedIdx]->Branch("pt", &m_pt[chargedIdx], "pt/F");
106  m_tree[chargedIdx]->Branch("trkTheta", &m_trkTheta[chargedIdx], "trkTheta/F");
107  m_tree[chargedIdx]->Branch("trkPhi", &m_trkPhi[chargedIdx], "trkPhi/F");
108  m_tree[chargedIdx]->Branch("clusterTheta", &m_clusterTheta[chargedIdx], "clusterTheta/F");
109  m_tree[chargedIdx]->Branch("clusterPhi", &m_clusterPhi[chargedIdx], "clusterPhi/F");
110  m_tree[chargedIdx]->Branch("clusterReg", &m_clusterReg[chargedIdx], "clusterReg/F");
111  m_tree[chargedIdx]->Branch("trackClusterMatch", &m_trackClusterMatch[chargedIdx], "trackClusterMatch/F");
112  m_tree[chargedIdx]->Branch("logl_sig", &m_logl_sig[chargedIdx], "logl_sig/F");
113  m_tree[chargedIdx]->Branch("logl_bkg", &m_logl_bkg[chargedIdx], "logl_bkg/F");
114  m_tree[chargedIdx]->Branch("deltalogl_sig_bkg", &m_deltalogl_sig_bkg[chargedIdx], "deltalogl_sig_bkg/F");
115  m_tree[chargedIdx]->Branch("pids_glob", &m_pids_glob[chargedIdx]);
116 
117  }
118 }
119 
121 {
122 }
123 
125 {
126 
127  std::string chargedPdgIdStr;
128 
129  for (const auto& chargedPdgId : m_inputPdgIdSet) {
130 
131  const auto chargedHypo = Const::chargedStableSet.find(std::abs(chargedPdgId));
132 
133  // If merging particles and antiparticles for this hypo, no need to loop twice:
134  // fill one TTree for the '+' charged pdgId only.
135  if (m_mergeChargeFlagByHypo[chargedHypo] and chargedPdgId < 0) continue;
136 
137  // Get the idx of this pdgId in the Const::chargedStableSet
138  auto chargedIdx = chargedHypo.getIndex();
139 
140  if (chargedPdgId < 0) {
141  // Add offset to idx for antiparticles.
142  chargedIdx += Const::ChargedStable::c_SetSize;
143  }
144 
145  // Initialise branches to unphysical values.
146  m_p[chargedIdx] = std::numeric_limits<float>::quiet_NaN();
147  m_pt[chargedIdx] = std::numeric_limits<float>::quiet_NaN();
148  m_trkTheta[chargedIdx] = std::numeric_limits<float>::quiet_NaN();
149  m_trkPhi[chargedIdx] = std::numeric_limits<float>::quiet_NaN();
150  m_clusterTheta[chargedIdx] = std::numeric_limits<float>::quiet_NaN();
151  m_clusterPhi[chargedIdx] = std::numeric_limits<float>::quiet_NaN();
152  m_clusterReg[chargedIdx] = std::numeric_limits<float>::quiet_NaN();
153  m_trackClusterMatch[chargedIdx] = std::numeric_limits<float>::quiet_NaN();
154  m_logl_sig[chargedIdx] = std::numeric_limits<float>::quiet_NaN();
155  m_logl_bkg[chargedIdx] = std::numeric_limits<float>::quiet_NaN();
156  m_deltalogl_sig_bkg[chargedIdx] = std::numeric_limits<float>::quiet_NaN();
157  for (const auto& chargedStable : Const::chargedStableSet) {
158  m_pids_glob[chargedIdx][chargedStable.getIndex()] = std::numeric_limits<float>::quiet_NaN();
159  }
160 
161  StoreArray<MCParticle> particles;
162 
163  for (const auto& particle : particles) {
164 
165  if (!particle.hasStatus(MCParticle::c_PrimaryParticle)) continue; // Only check primaries.
166  if (particle.hasStatus(MCParticle::c_Initial)) continue; // Ignore initial particles.
167  if (particle.hasStatus(MCParticle::c_IsVirtual)) continue; // Ignore virtual particles.
168 
169  // Skip all particles expect for the one of interest.
170  // If merging particles and antiparticles for this pdgId, use abs so both charges are considered for the MCParticles.
171  if (m_mergeChargeFlagByHypo[chargedHypo]) {
172  // Charge-agnostic check.
173  if (std::abs(particle.getPDG()) != std::abs(chargedPdgId)) continue;
174  } else {
175  // Charge-dependent check.
176  if (particle.getPDG() != chargedPdgId) continue;
177  }
178 
179  // Get the matching track w/ max momentum.
180  int itrack(0);
181  int itrack_max(-1);
182  double p_max(-999.0);
183  for (const auto& track : particle.getRelationsFrom<Track>()) {
184  const auto fitRes = track.getTrackFitResultWithClosestMass(Const::pion);
185  if (!fitRes) continue;
186  if (fitRes->getMomentum().Mag() > p_max) {
187  p_max = fitRes->getMomentum().Mag();
188  itrack_max = itrack;
189  }
190  itrack++;
191  }
192  if (itrack_max < 0) continue; // Go to next particle if no track found.
193 
194  const auto track = particle.getRelationsFrom<Track>()[itrack_max];
195  const auto fitRes = track->getTrackFitResultWithClosestMass(Const::pion);
196 
197  m_p[chargedIdx] = p_max;
198  m_pt[chargedIdx] = fitRes->get4Momentum().Pt();
199  m_trkTheta[chargedIdx] = fitRes->get4Momentum().Theta();
200  m_trkPhi[chargedIdx] = fitRes->get4Momentum().Phi();
201 
202  // Get the index of the ECL cluster matching this track.
203  int icluster_match(-1);
204  auto eclClusters = track->getRelationsTo<ECLCluster>();
205  for (unsigned int icluster(0); icluster < eclClusters.size(); ++icluster) {
206  const auto eclCluster = eclClusters[icluster];
207  if (!eclCluster->hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons)) continue;
208  if (!eclCluster->isTrack()) continue;
209  icluster_match = icluster;
210  break;
211  }
212  // If no cluster match, skip to next particle, but keep track of counter.
213  if (icluster_match < 0) {
214  m_trackClusterMatch[chargedIdx] = 0;
215  continue;
216  }
217 
218  const auto eclCluster = eclClusters[icluster_match];
219 
220  m_clusterTheta[chargedIdx] = eclCluster->getTheta();
221  m_clusterPhi[chargedIdx] = eclCluster->getPhi();
222  m_clusterReg[chargedIdx] = eclCluster->getDetectorRegion();
223 
224  m_trackClusterMatch[chargedIdx] = 1;
225 
226  const auto eclLikelihood = track->getRelated<ECLPidLikelihood>();
227 
228  // The "signal" likelihood corresponds to the current chargedPdgId.
229  const auto chargedStableSig = Const::chargedStableSet.find(std::abs(chargedPdgId));
230  // For deltaLogL, we do a binary comparison sig/bkg.
231  // If sig=pion, use bkg=kaon. Otherwise, bkg=pion.
232  const auto chargedStableBkg = (chargedStableSig == Const::pion) ? Const::kaon : Const::pion;
233 
234  double lh_sig = eclLikelihood->getLikelihood(chargedStableSig);
235  double lh_bkg = eclLikelihood->getLikelihood(chargedStableBkg);
236 
237  m_logl_sig[chargedIdx] = log(lh_sig);
238  m_logl_bkg[chargedIdx] = log(lh_bkg);
239  m_deltalogl_sig_bkg[chargedIdx] = log(lh_bkg) - log(lh_sig);
240 
241  // For the current charged particle candidate, store the global likelihood ratio for all hypotheses.
242  double lh_all(0);
243  for (const auto& chargedStable : Const::chargedStableSet) {
244  lh_all += eclLikelihood->getLikelihood(chargedStable);
245  }
246  for (const auto& chargedStable : Const::chargedStableSet) {
247  m_pids_glob[chargedIdx][chargedStable.getIndex()] = eclLikelihood->getLikelihood(chargedStable) / lh_all;
248  }
249 
250  }
251 
252  m_tree[chargedIdx]->Fill();
253 
254  }
255 }
256 
258 {
259 }
260 
262 {
263 
264  for (const auto& chargedPdgId : m_inputPdgIdSet) {
265 
266  // Define the charged stable particle ("sample") corresponding to the current pdgId.
267  const auto chargedStableSample = Const::chargedStableSet.find(std::abs(chargedPdgId));
268 
269  const auto mergeCharge = m_mergeChargeFlagByHypo[chargedStableSample];
270 
271  // If partcile/antiparticle merging is active for this hypo, the results are stored for the '+' particle only.
272  if (mergeCharge && chargedPdgId < 0) continue;
273 
274  // Extract the sign of the charge.
275  const auto chargeSign = static_cast<int>(chargedPdgId / std::abs(chargedPdgId));
276 
277  // What we call "signal" is equal to the sample under consideration.
278  const auto chargedStableSig = chargedStableSample;
279 
280  // What we call "background" depends on the current "signal" particle hypothesis.
281  const auto chargedStableBkg = (chargedStableSig == Const::pion) ? Const::kaon : Const::pion;
282 
283  // Get the idx of this sample's pdgId to retrieve the correct TTree and output TFile.
284  // Remember to add offset for antiparticles.
285  auto chargedSampleIdx = (chargeSign > 0) ? chargedStableSig.getIndex() : chargedStableSig.getIndex() +
287 
288  m_outputFile[chargedSampleIdx]->cd();
289 
290  auto pdgIdDesc = (mergeCharge) ? std::to_string(chargedPdgId) + " and -" + std::to_string(chargedPdgId) : std::to_string(
291  chargedPdgId);
292 
293  // Add summary description of validation file content.
294  TNamed("Description", TString::Format("ECL Charged PID control plots for charged stable particles/antiparticles ; Sample PDG = %s",
295  pdgIdDesc.c_str()).Data()).Write();
296 
297  // Dump plots of PID variables.
298  dumpPIDVars(m_tree[chargedSampleIdx], chargedStableSig, chargeSign, chargedStableBkg, mergeCharge);
299  // Dump plots of PID "signal" efficiency for this "sample".
300  dumpPIDEfficiencyFakeRate(m_tree[chargedSampleIdx], chargedStableSample, chargeSign, chargedStableSig, mergeCharge);
301  // For pions, dump also the pi->lep fake rate.
302  if (chargedStableSample == Const::pion) {
303  dumpPIDEfficiencyFakeRate(m_tree[chargedSampleIdx], chargedStableSample, chargeSign, Const::electron, mergeCharge);
304  dumpPIDEfficiencyFakeRate(m_tree[chargedSampleIdx], chargedStableSample, chargeSign, Const::muon, mergeCharge);
305  }
306  // Dump plots of matching efficiency for this "sample".
307  dumpTrkClusMatchingEfficiency(m_tree[chargedSampleIdx], chargedStableSample, chargeSign, mergeCharge);
308 
309  // Write the TTree to file if requested.
310  if (m_saveValidationTree) {
311  m_tree[chargedSampleIdx]->Write();
312  }
313 
314  m_outputFile[chargedSampleIdx]->Close();
315 
316  }
317 }
318 
319 void ECLChargedPIDDataAnalysisValidationModule::dumpPIDVars(TTree* sampleTree, const Const::ChargedStable& sigHypo,
320  const int sigCharge, const Const::ChargedStable& bkgHypo, bool mergeSigCharge)
321 {
322 
323  // Get the idx and pdgId of the input sample particle.
324  // This corresponds by construction to the "signal" hypothesis for the likelihood and DeltaLogL.
325  const auto sigHypoIdx = sigHypo.getIndex();
326  const auto sigHypoPdgId = sigHypo.getPDGCode();
327 
328  // Get the pdgId of the "background" hypothesis to test for DeltaLogL.
329  const auto bkgHypoPdgId = bkgHypo.getPDGCode();
330 
331  // Access the "signal" hypothesis's PID component in the sample's
332  // TTree vector branch of global PID values via the idx.
333  TString pidSigBranch = TString::Format("pids_glob[%i]", sigHypoIdx);
334 
335  // Histogram of global PID distribution for the sample particle's signal hypo.
336  TString h_pid_name = TString::Format("h_pid_sig_%i", sigHypoPdgId);
337  TH1F* h_pid = new TH1F(h_pid_name.Data(), h_pid_name.Data(), 50, -0.5, 1.2);
338  h_pid->GetXaxis()->SetTitle(TString::Format("Likelihood ratio (%i/ALL) (ECL)", sigHypoPdgId).Data());
339 
340  // Histogram of deltalogl.
341  TString h_deltalogl_name = TString::Format("h_deltalogl_bkg_%i_sig_%i", bkgHypoPdgId, sigHypoPdgId);
342  double deltalogl_min = -20.0;
343  double deltalogl_max = 20.0;
344  TH1F* h_deltalogl = new TH1F(h_deltalogl_name.Data(), h_deltalogl_name.Data(), 40, deltalogl_min, deltalogl_max);
345  h_deltalogl->GetXaxis()->SetTitle(TString::Format("#Deltaln(L) (%i/%i) (ECL)", bkgHypoPdgId, sigHypoPdgId).Data());
346 
347  // Histogram of track-cluster match flag.
348  TString h_trkclusmatch_name = TString::Format("h_trkclusmatch_sig_%i", sigHypoPdgId);
349  TH1F* h_trkclusmatch = new TH1F(h_trkclusmatch_name.Data(), h_trkclusmatch_name.Data(), 4, -1.5, 2.5);
350  h_trkclusmatch->GetXaxis()->SetTitle(TString::Format("Track-ECLCluster match (%i)", sigHypoPdgId).Data());
351 
352  // Dump histos from TTree.
353  sampleTree->Project(h_pid_name.Data(), pidSigBranch.Data());
354  sampleTree->Project(h_deltalogl_name.Data(), "deltalogl_sig_bkg");
355  sampleTree->Project(h_trkclusmatch_name.Data(), "trackClusterMatch");
356 
357  // Make sure the plots show the u/oflow.
358  paintUnderOverflow(h_pid);
359  paintUnderOverflow(h_deltalogl);
360  paintUnderOverflow(h_trkclusmatch);
361 
362  h_pid->SetOption("HIST");
363  h_deltalogl->SetOption("HIST");
364  h_trkclusmatch->SetOption("HIST");
365 
366  // MetaOptions string.
367  std::string metaopts("pvalue-warn=0.1,pvalue-error=0.01");
368  std::string shifteropt("");
369  // Electron plots should be visible to the shifter by default.
370  if (sigHypo == Const::electron) {
371  shifteropt = "shifter,";
372  }
373 
374  auto pdgIdDesc = (!mergeSigCharge) ? std::to_string(sigHypoPdgId * sigCharge) : std::to_string(
375  sigHypoPdgId) + " and -" + std::to_string(sigHypoPdgId);
376 
377  // Add histogram info.
378  h_pid->GetListOfFunctions()->Add(new TNamed("Description",
379  TString::Format("Sample PDG = %s ; ECL global PID(%i) distribution. U/O flow is added to first (last) bin.",
380  pdgIdDesc.c_str(),
381  sigHypoPdgId).Data()));
382  h_pid->GetListOfFunctions()->Add(new TNamed("Check",
383  "The more peaked at 1, the better. Non-zero O-flow indicates either failure of MC matching for reco tracks (unlikely), or failure of track-ECL-cluster matching (more likely). Both cases result in PID=nan."));
384  h_pid->GetListOfFunctions()->Add(new TNamed("Contact", "Marco Milesi. marco.milesi@desy.de"));
385  h_pid->GetListOfFunctions()->Add(new TNamed("MetaOptions", (shifteropt + metaopts).c_str()));
386 
387  h_deltalogl->GetListOfFunctions()->Add(new TNamed("Description",
388  TString::Format("Sample PDG = %s ; ECL distribution of binary $\\Delta log(L)$ = log(L(%i)) - log(L(%i)). U/O flow is added to first (last) bin.",
389  pdgIdDesc.c_str(),
390  bkgHypoPdgId,
391  sigHypoPdgId).Data()));
392  h_deltalogl->GetListOfFunctions()->Add(new TNamed("Check",
393  "Basic metric for signal/bkg separation. The more negative, the better separation is achieved. Non-zero U-flow indicates a non-normal PDF value (of sig OR bkg) for some p,clusterTheta range, which might be due to a non-optimal definition of the x-axis range of the PDF templates. Non-zero O-flow indicates either failure of MC matching for reco tracks (unlikely), or failure of track-ECL-cluster matching (more likely)."));
394  h_deltalogl->GetListOfFunctions()->Add(new TNamed("Contact", "Marco Milesi. marco.milesi@desy.de"));
395  h_deltalogl->GetListOfFunctions()->Add(new TNamed("MetaOptions", (shifteropt + metaopts).c_str()));
396 
397  h_trkclusmatch->GetListOfFunctions()->Add(new TNamed("Description",
398  TString::Format("Sample PDG = %s ; Track-ECLCluster match flag distribution.",
399  pdgIdDesc.c_str()).Data()));
400  h_trkclusmatch->GetListOfFunctions()->Add(new TNamed("Check",
401  "The more peaked at 1, the better. Non-zero population in the bins w/ flag != 0|1 indicates failure of MC matching for reco tracks. In such cases, flag=nan."));
402  h_trkclusmatch->GetListOfFunctions()->Add(new TNamed("Contact", "Frank Meier. frank.meier@desy.de"));
403  h_trkclusmatch->GetListOfFunctions()->Add(new TNamed("MetaOptions", metaopts.c_str()));
404 
405  h_pid->Write();
406  h_deltalogl->Write();
407  h_trkclusmatch->Write();
408 
409  delete h_pid;
410  delete h_deltalogl;
411  delete h_trkclusmatch;
412 
413 }
414 
415 
416 void ECLChargedPIDDataAnalysisValidationModule::dumpPIDEfficiencyFakeRate(TTree* sampleTree, const Const::ChargedStable& sampleHypo,
417  const int sampleCharge, const Const::ChargedStable& sigHypo, bool mergeSampleCharge)
418 {
419 
420  // The ratio type: EFFICIENCY || FAKE RATE
421  const std::string ratioType = (sampleHypo == sigHypo) ? "Efficiency" : "FakeRate";
422 
423  // Get the *signed* pdgId of the input sample particle.
424  const int sampleHypoPdgId = sampleHypo.getPDGCode() * sampleCharge;
425 
426  // Get the idx and pdgId of the "signal" hypothesis to test.
427  const auto sigHypoIdx = sigHypo.getIndex();
428  const auto sigHypoPdgId = sigHypo.getPDGCode();
429 
430  // Access the "signal" hypothesis's PID component in the sample's
431  // TTree vector branch of global PID values via the idx.
432  TString pidSigCut = TString::Format("pids_glob[%i] > %f", sigHypoIdx, c_PID);
433 
434  // Histograms of p, clusterReg, clusterPhi... for "pass" (N, numerator) and "all" (D, denominator) events.
435  TString h_p_N_name = TString::Format("h_p_N_%i", sigHypoPdgId);
436  TString h_p_D_name = TString::Format("h_p_D_%i", sigHypoPdgId);
437  TH1F* h_p_N = new TH1F(h_p_N_name.Data(), "h_p_N", 10, 0.0, 5.0);
438  TH1F* h_p_D = new TH1F(h_p_D_name.Data(), "h_p_D", 10, 0.0, 5.0);
439 
440  TString h_th_N_name = TString::Format("h_th_N_%i", sigHypoPdgId);
441  TString h_th_D_name = TString::Format("h_th_D_%i", sigHypoPdgId);
442  TH1F* h_th_N = new TH1F(h_th_N_name.Data(), "h_th_N", m_th_binedges.size() - 1, m_th_binedges.data());
443  TH1F* h_th_D = new TH1F(h_th_D_name.Data(), "h_th_D", m_th_binedges.size() - 1, m_th_binedges.data());
444 
445  TString h_eclreg_N_name = TString::Format("h_eclreg_N_%i", sigHypoPdgId);
446  TString h_eclreg_D_name = TString::Format("h_eclreg_D_%i", sigHypoPdgId);
447  TH1F* h_eclreg_N = new TH1F(h_eclreg_N_name.Data(), "h_eclreg_N", 5, -0.5, 4.5);
448  TH1F* h_eclreg_D = new TH1F(h_eclreg_D_name.Data(), "h_eclreg_D", 5, -0.5, 4.5);
449 
450  TString h_phi_N_name = TString::Format("h_phi_N_%i", sigHypoPdgId);
451  TString h_phi_D_name = TString::Format("h_phi_D_%i", sigHypoPdgId);
452  TH1F* h_phi_N = new TH1F(h_phi_N_name.Data(), "h_phi_N", 5, -3.14159, 3.14159);
453  TH1F* h_phi_D = new TH1F(h_phi_D_name.Data(), "h_phi_D", 5, -3.14159, 3.14159);
454 
455  // Fill the histograms from the sample's TTree.
456 
457  sampleTree->Project(h_p_N_name.Data(), "p", pidSigCut.Data());
458  sampleTree->Project(h_p_D_name.Data(), "p");
459 
460  sampleTree->Project(h_th_N_name.Data(), "clusterTheta", pidSigCut.Data());
461  sampleTree->Project(h_th_D_name.Data(), "clusterTheta");
462 
463  sampleTree->Project(h_eclreg_N_name.Data(), "clusterReg", pidSigCut.Data());
464  sampleTree->Project(h_eclreg_D_name.Data(), "clusterReg");
465  paintUnderOverflow(h_eclreg_N);
466  paintUnderOverflow(h_eclreg_D);
467 
468  sampleTree->Project(h_phi_N_name.Data(), "clusterPhi", pidSigCut.Data());
469  sampleTree->Project(h_phi_D_name.Data(), "clusterPhi");
470 
471  // Compute the efficiency/fake rate.
472 
473  TString pid_glob_ratio_p_name = TString::Format("pid_glob_%i_%s__VS_p", sigHypoPdgId, ratioType.c_str());
474  TString pid_glob_ratio_th_name = TString::Format("pid_glob_%i_%s__VS_th", sigHypoPdgId, ratioType.c_str());
475  TString pid_glob_ratio_eclreg_name = TString::Format("pid_glob_%i_%s__VS_eclreg", sigHypoPdgId, ratioType.c_str());
476  TString pid_glob_ratio_phi_name = TString::Format("pid_glob_%i_%s__VS_phi", sigHypoPdgId, ratioType.c_str());
477 
478  // MetaOptions string.
479  std::string metaopts("pvalue-warn=0.01,pvalue-error=0.001,nostats");
480  std::string shifteropt("");
481  // Electron plots should be visible to the shifter by default.
482  if (sampleHypo == Const::electron || sigHypo == Const::electron) {
483  shifteropt = "shifter,";
484  }
485 
486  auto pdgIdDesc = (!mergeSampleCharge) ? std::to_string(sampleHypoPdgId) : std::to_string(std::abs(
487  sampleHypoPdgId)) + " and -" + std::to_string(std::abs(sampleHypoPdgId));
488 
489  if (TEfficiency::CheckConsistency(*h_p_N, *h_p_D)) {
490 
491  TEfficiency* t_pid_glob_ratio_p = new TEfficiency(*h_p_N, *h_p_D);
492  t_pid_glob_ratio_p->SetName(pid_glob_ratio_p_name.Data());
493  t_pid_glob_ratio_p->SetTitle(TString::Format("%s;p [GeV/c];#varepsilon/f", pid_glob_ratio_p_name.Data()).Data());
494 
495  t_pid_glob_ratio_p->SetConfidenceLevel(0.683);
496  t_pid_glob_ratio_p->SetStatisticOption(TEfficiency::kBUniform);
497  t_pid_glob_ratio_p->SetPosteriorMode();
498 
499  t_pid_glob_ratio_p->GetListOfFunctions()->Add(new TNamed("Description",
500  TString::Format("Sample PDG = %s ; %s of ECL global PID(%i) > %.2f as a function of $p_{trk}$.",
501  pdgIdDesc.c_str(),
502  ratioType.c_str(),
503  sigHypoPdgId,
504  c_PID).Data()));
505  t_pid_glob_ratio_p->GetListOfFunctions()->Add(new TNamed("Check",
506  "Shape should be consistent. Obviously, check for decreasing efficiency / increasing fake rate."));
507  t_pid_glob_ratio_p->GetListOfFunctions()->Add(new TNamed("Contact", "Marco Milesi. marco.milesi@desy.de"));
508  t_pid_glob_ratio_p->GetListOfFunctions()->Add(new TNamed("MetaOptions", (shifteropt + metaopts).c_str()));
509 
510  t_pid_glob_ratio_p->Write();
511 
512  delete t_pid_glob_ratio_p;
513 
514  }
515  if (TEfficiency::CheckConsistency(*h_th_N, *h_th_D)) {
516 
517  TEfficiency* t_pid_glob_ratio_th = new TEfficiency(*h_th_N, *h_th_D);
518  t_pid_glob_ratio_th->SetName(pid_glob_ratio_th_name.Data());
519  t_pid_glob_ratio_th->SetTitle(TString::Format("%s;#theta_{cluster} [rad];#varepsilon/f", pid_glob_ratio_th_name.Data()).Data());
520 
521  t_pid_glob_ratio_th->SetConfidenceLevel(0.683);
522  t_pid_glob_ratio_th->SetStatisticOption(TEfficiency::kBUniform);
523  t_pid_glob_ratio_th->SetPosteriorMode();
524 
525  t_pid_glob_ratio_th->GetListOfFunctions()->Add(new TNamed("Description",
526  TString::Format("Sample PDG = %s ; %s of ECL global PID(%i) > %.2f as a function of $\\theta_{cluster}$.",
527  pdgIdDesc.c_str(),
528  ratioType.c_str(),
529  sigHypoPdgId,
530  c_PID).Data()));
531  t_pid_glob_ratio_th->GetListOfFunctions()->Add(new TNamed("Check",
532  "Shape should be consistent. Obviously, check for decreasing efficiency / increasing fake rate."));
533  t_pid_glob_ratio_th->GetListOfFunctions()->Add(new TNamed("Contact", "Marco Milesi. marco.milesi@desy.de"));
534  t_pid_glob_ratio_th->GetListOfFunctions()->Add(new TNamed("MetaOptions", (shifteropt + metaopts).c_str()));
535 
536  t_pid_glob_ratio_th->Write();
537 
538  delete t_pid_glob_ratio_th;
539  }
540  if (TEfficiency::CheckConsistency(*h_eclreg_N, *h_eclreg_D)) {
541 
542  TEfficiency* t_pid_glob_ratio_eclreg = new TEfficiency(*h_eclreg_N, *h_eclreg_D);
543  t_pid_glob_ratio_eclreg->SetName(pid_glob_ratio_eclreg_name.Data());
544  t_pid_glob_ratio_eclreg->SetTitle(TString::Format("%s;ECL Region;#varepsilon/f", pid_glob_ratio_eclreg_name.Data()).Data());
545 
546  t_pid_glob_ratio_eclreg->SetConfidenceLevel(0.683);
547  t_pid_glob_ratio_eclreg->SetStatisticOption(TEfficiency::kBUniform);
548  t_pid_glob_ratio_eclreg->SetPosteriorMode();
549 
550  t_pid_glob_ratio_eclreg->GetListOfFunctions()->Add(new TNamed("Description",
551  TString::Format("Sample PDG = %s ; %s of ECL global PID(%i) > %.2f as a function of ECL cluster region ($\\theta_{cluster}$). Regions are labelled: 0 (outside ECL acceptance), 1 (ECL FWD), 2 (ECL Barrel), 3 (ECL BWD), 4 (ECL FWD/BWD gaps).",
552  pdgIdDesc.c_str(),
553  ratioType.c_str(),
554  sigHypoPdgId,
555  c_PID).Data()));
556  t_pid_glob_ratio_eclreg->GetListOfFunctions()->Add(new TNamed("Check",
557  "Shape should be consistent. Obviously, check for decreasing efficiency / increasing fake rate."));
558  t_pid_glob_ratio_eclreg->GetListOfFunctions()->Add(new TNamed("Contact", "Marco Milesi. marco.milesi@desy.de"));
559  t_pid_glob_ratio_eclreg->GetListOfFunctions()->Add(new TNamed("MetaOptions", metaopts.c_str()));
560 
561  t_pid_glob_ratio_eclreg->Write();
562 
563  delete t_pid_glob_ratio_eclreg;
564 
565  }
566  if (TEfficiency::CheckConsistency(*h_phi_N, *h_phi_D)) {
567 
568  TEfficiency* t_pid_glob_ratio_phi = new TEfficiency(*h_phi_N, *h_phi_D);
569  t_pid_glob_ratio_phi->SetName(pid_glob_ratio_phi_name.Data());
570  t_pid_glob_ratio_phi->SetTitle(TString::Format("%s;#phi_{cluster} [rad];#varepsilon/f", pid_glob_ratio_phi_name.Data()).Data());
571 
572  t_pid_glob_ratio_phi->SetConfidenceLevel(0.683);
573  t_pid_glob_ratio_phi->SetStatisticOption(TEfficiency::kBUniform);
574  t_pid_glob_ratio_phi->SetPosteriorMode();
575 
576  t_pid_glob_ratio_phi->GetListOfFunctions()->Add(new TNamed("Description",
577  TString::Format("Sample PDG = %s ; %s of ECL global PID(%i) > %.2f as a function of $\\phi_{cluster}$.",
578  pdgIdDesc.c_str(),
579  ratioType.c_str(),
580  sigHypoPdgId,
581  c_PID).Data()));
582  t_pid_glob_ratio_phi->GetListOfFunctions()->Add(new TNamed("Check",
583  "Shape should be consistent. Obviously, check for decreasing efficiency / increasing fake rate."));
584  t_pid_glob_ratio_phi->GetListOfFunctions()->Add(new TNamed("Contact", "Marco Milesi. marco.milesi@desy.de"));
585  t_pid_glob_ratio_phi->GetListOfFunctions()->Add(new TNamed("MetaOptions", (shifteropt + metaopts).c_str()));
586 
587  t_pid_glob_ratio_phi->Write();
588 
589  delete t_pid_glob_ratio_phi;
590  }
591 
592  delete h_p_N;
593  delete h_p_D;
594  delete h_th_N;
595  delete h_th_D;
596  delete h_eclreg_N;
597  delete h_eclreg_D;
598  delete h_phi_N;
599  delete h_phi_D;
600 
601 }
602 
603 
604 void ECLChargedPIDDataAnalysisValidationModule::dumpTrkClusMatchingEfficiency(TTree* tree, const Const::ChargedStable& sampleHypo,
605  const int sampleCharge, bool mergeSampleCharge)
606 {
607 
608  // Get the (unsigned) pdgId of the input sample particle.
609  const auto sampleHypoPdgId = sampleHypo.getPDGCode();
610 
611  // Histograms of pt, clusterTheta, clusterPhi... for "pass" (N, numerator) and "all" (D, denominator) events.
612  TString h_pt_N_name = TString::Format("h_pt_N_%i", sampleHypoPdgId);
613  TString h_pt_D_name = TString::Format("h_pt_D_%i", sampleHypoPdgId);
614  TH1F* h_pt_N = new TH1F(h_pt_N_name.Data(), "h_pt_N", 10, 0.0, 5.0);
615  TH1F* h_pt_D = new TH1F(h_pt_D_name.Data(), "h_pt_D", 10, 0.0, 5.0);
616 
617  TString h_th_N_name = TString::Format("h_th_N_%i", sampleHypoPdgId);
618  TString h_th_D_name = TString::Format("h_th_D_%i", sampleHypoPdgId);
619  TH1F* h_th_N = new TH1F(h_th_N_name.Data(), "h_th_N", m_th_binedges.size() - 1, m_th_binedges.data());
620  TH1F* h_th_D = new TH1F(h_th_D_name.Data(), "h_th_D", m_th_binedges.size() - 1, m_th_binedges.data());
621 
622  TString h_phi_N_name = TString::Format("h_phi_N_%i", sampleHypoPdgId);
623  TString h_phi_D_name = TString::Format("h_phi_D_%i", sampleHypoPdgId);
624  TH1F* h_phi_N = new TH1F(h_phi_N_name.Data(), "h_phi_N", 5, -3.14159, 3.14159);
625  TH1F* h_phi_D = new TH1F(h_phi_D_name.Data(), "h_phi_D", 5, -3.14159, 3.14159);
626 
627  TString match_cut_N("trackClusterMatch == 1");
628  TString match_cut_D("trackClusterMatch >= 0");
629 
630  // Fill the histograms from the sample's TTree.
631 
632  tree->Project(h_pt_N_name.Data(), "pt", match_cut_N.Data());
633  tree->Project(h_pt_D_name.Data(), "pt", match_cut_D.Data());
634 
635  tree->Project(h_th_N_name.Data(), "trkTheta", match_cut_N.Data());
636  tree->Project(h_th_D_name.Data(), "trkTheta", match_cut_D.Data());
637 
638  tree->Project(h_phi_N_name.Data(), "trkPhi", match_cut_N.Data());
639  tree->Project(h_phi_D_name.Data(), "trkPhi", match_cut_D.Data());
640 
641  // Compute the efficiency.
642 
643  TString match_eff_pt_name = TString::Format("trkclusmatch_%i_Efficiency__VS_pt", sampleHypoPdgId);
644  TString match_eff_th_name = TString::Format("trkclusmatch_%i_Efficiency__VS_th", sampleHypoPdgId);
645  TString match_eff_phi_name = TString::Format("trkclusmatch_%i_Efficiency__VS_phi", sampleHypoPdgId);
646 
647  // MetaOptions string.
648  std::string metaopts("pvalue-warn=0.01,pvalue-error=0.001,nostats");
649  std::string shifteropt("");
650  // Electron plots should be visible to the shifter by default.
651  if (sampleHypo == Const::electron) {
652  shifteropt = "shifter,";
653  }
654 
655  auto pdgIdDesc = (!mergeSampleCharge) ? std::to_string(sampleHypoPdgId * sampleCharge) : std::to_string(
656  sampleHypoPdgId) + " and -" + std::to_string(sampleHypoPdgId);
657 
658  if (TEfficiency::CheckConsistency(*h_pt_N, *h_pt_D)) {
659 
660  TEfficiency* t_match_eff_pt = new TEfficiency(*h_pt_N, *h_pt_D);
661  t_match_eff_pt->SetName(match_eff_pt_name.Data());
662  t_match_eff_pt->SetTitle(TString::Format("%s;p_{T}^{trk} [GeV/c];#varepsilon", match_eff_pt_name.Data()).Data());
663  t_match_eff_pt->SetTitle(match_eff_pt_name.Data());
664 
665  t_match_eff_pt->SetConfidenceLevel(0.683);
666  t_match_eff_pt->SetStatisticOption(TEfficiency::kBUniform);
667  t_match_eff_pt->SetPosteriorMode();
668 
669  t_match_eff_pt->GetListOfFunctions()->Add(new TNamed("Description",
670  TString::Format("Sample PDG = %s ; Efficiency of track-ECL-cluster matching as a function of $p_{T}^{trk}$.",
671  pdgIdDesc.c_str()).Data()));
672  t_match_eff_pt->GetListOfFunctions()->Add(new TNamed("Check",
673  "Shape should be consistent. Obviously, check for decreasing efficiency."));
674  t_match_eff_pt->GetListOfFunctions()->Add(new TNamed("Contact", "Frank Meier. frank.meier@desy.de"));
675  t_match_eff_pt->GetListOfFunctions()->Add(new TNamed("MetaOptions", (shifteropt + metaopts).c_str()));
676 
677  t_match_eff_pt->Write();
678 
679  delete t_match_eff_pt;
680 
681  }
682  if (TEfficiency::CheckConsistency(*h_th_N, *h_th_D)) {
683 
684  TEfficiency* t_match_eff_th = new TEfficiency(*h_th_N, *h_th_D);
685  t_match_eff_th->SetName(match_eff_th_name.Data());
686  t_match_eff_th->SetTitle(TString::Format("%s;#theta_{trk} [rad];#varepsilon", match_eff_th_name.Data()).Data());
687  t_match_eff_th->SetTitle(match_eff_th_name.Data());
688 
689  t_match_eff_th->SetConfidenceLevel(0.683);
690  t_match_eff_th->SetStatisticOption(TEfficiency::kBUniform);
691  t_match_eff_th->SetPosteriorMode();
692 
693  t_match_eff_th->GetListOfFunctions()->Add(new TNamed("Description",
694  TString::Format("Sample PDG = %s ; Efficiency of track-ECL-cluster matching as a function of $\\theta_{trk}$.",
695  pdgIdDesc.c_str()).Data()));
696  t_match_eff_th->GetListOfFunctions()->Add(new TNamed("Check",
697  "Shape should be consistent. Obviously, check for decreasing efficiency."));
698  t_match_eff_th->GetListOfFunctions()->Add(new TNamed("Contact", "Frank Meier. frank.meier@desy.de"));
699  t_match_eff_th->GetListOfFunctions()->Add(new TNamed("MetaOptions", (shifteropt + metaopts).c_str()));
700 
701  t_match_eff_th->Write();
702 
703  delete t_match_eff_th;
704 
705  }
706  if (TEfficiency::CheckConsistency(*h_phi_N, *h_phi_D)) {
707 
708  TEfficiency* t_match_eff_phi = new TEfficiency(*h_phi_N, *h_phi_D);
709  t_match_eff_phi->SetName(match_eff_phi_name.Data());
710  t_match_eff_phi->SetTitle(TString::Format("%s;#phi_{trk} [rad];#varepsilon", match_eff_phi_name.Data()).Data());
711 
712  t_match_eff_phi->SetConfidenceLevel(0.683);
713  t_match_eff_phi->SetStatisticOption(TEfficiency::kBUniform);
714  t_match_eff_phi->SetPosteriorMode();
715 
716  t_match_eff_phi->GetListOfFunctions()->Add(new TNamed("Description",
717  TString::Format("Sample PDG = %s ; Efficiency of track-ECL-cluster matching as a function of $\\phi_{trk}$.",
718  pdgIdDesc.c_str()).Data()));
719  t_match_eff_phi->GetListOfFunctions()->Add(new TNamed("Check",
720  "Shape should be consistent. Obviously, check for decreasing efficiency."));
721  t_match_eff_phi->GetListOfFunctions()->Add(new TNamed("Contact", "Frank Meier. frank.meier@desy.de"));
722  t_match_eff_phi->GetListOfFunctions()->Add(new TNamed("MetaOptions", (shifteropt + metaopts).c_str()));
723 
724  t_match_eff_phi->Write();
725 
726  delete t_match_eff_phi;
727  }
728 
729  delete h_pt_N;
730  delete h_pt_D;
731  delete h_th_N;
732  delete h_th_D;
733  delete h_phi_N;
734  delete h_phi_D;
735 
736 }
737 
738 
739 void ECLChargedPIDDataAnalysisValidationModule::paintUnderOverflow(TH1F* h)
740 {
741 
742  auto nentries = h->GetEntries();
743  auto nbins_vis = h->GetNbinsX();
744 
745  // Get the content and error of first/last visible bin.
746  float bin_vis_first = h->GetBinContent(1);
747  float bin_vis_last = h->GetBinContent(nbins_vis);
748  float bin_vis_first_err = h->GetBinError(1);
749  float bin_vis_last_err = h->GetBinError(nbins_vis);
750 
751  // Get the content and error of u/oflow bins.
752  float bin_uflow = h->GetBinContent(0);
753  float bin_oflow = h->GetBinContent(nbins_vis + 1);
754  float bin_uflow_err = h->GetBinError(0);
755  float bin_oflow_err = h->GetBinError(nbins_vis + 1);
756 
757  // Reset first/last visible bins to include u/oflow.
758  h->SetBinContent(1, bin_vis_first + bin_uflow);
759  h->SetBinError(1, sqrt(bin_vis_first_err * bin_vis_first_err + bin_uflow_err * bin_uflow_err));
760  h->SetBinContent(nbins_vis, bin_vis_last + bin_oflow);
761  h->SetBinError(nbins_vis, sqrt(bin_vis_last_err * bin_vis_last_err + bin_oflow_err * bin_oflow_err));
762 
763  // Reset total entries to the original value.
764  h->SetEntries(nentries);
765 
766 }
Belle2::Const::ChargedStable::c_SetSize
static const unsigned int c_SetSize
Number of elements (for use in array bounds etc.)
Definition: Const.h:491
Belle2::ECLChargedPIDDataAnalysisValidationModule::terminate
virtual void terminate() override
Termination action.
Definition: ECLChargedPIDDataAnalysisValidationModule.cc:261
Belle2::ECLCluster
ECL cluster data.
Definition: ECLCluster.h:39
Belle2::ECLChargedPIDDataAnalysisValidationModule::m_mergeChargeFlagByHypo
std::map< Const::ChargedStable, bool > m_mergeChargeFlagByHypo
A map to tell for each charged stable particle hypothesis whether particle and antiparticle should be...
Definition: ECLChargedPIDDataAnalysisValidationModule.h:110
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::MCParticle::c_IsVirtual
@ c_IsVirtual
bit 4: Particle is virtual and not going to Geant4.
Definition: MCParticle.h:66
Belle2::ECLChargedPIDDataAnalysisValidationModule::m_outputFile
std::vector< TFile * > m_outputFile
Output ROOT::TFile that contains the info to plot.
Definition: ECLChargedPIDDataAnalysisValidationModule.h:123
Belle2::Const::chargedStableSet
static const ParticleSet chargedStableSet
set of charged stable particles
Definition: Const.h:494
Belle2::ECLCluster::EHypothesisBit::c_nPhotons
@ c_nPhotons
CR is split into n photons (N1)
Belle2::ECLChargedPIDDataAnalysisValidationModule::m_inputPdgIdList
std::vector< int > m_inputPdgIdList
The pdgId list of the charged stable particles of interest.
Definition: ECLChargedPIDDataAnalysisValidationModule.h:99
Belle2::Const::ParticleType::getPDGCode
int getPDGCode() const
PDG code.
Definition: Const.h:349
Belle2::ECLChargedPIDDataAnalysisValidationModule::event
virtual void event() override
Called once for each event.
Definition: ECLChargedPIDDataAnalysisValidationModule.cc:124
Belle2::ECLChargedPIDDataAnalysisValidationModule::m_logl_bkg
std::vector< float > m_logl_bkg
Log-likelihood for the "background" particle hypothesis.
Definition: ECLChargedPIDDataAnalysisValidationModule.h:235
Belle2::ECLChargedPIDDataAnalysisValidationModule::m_pt
std::vector< float > m_pt
Track transverse momentum in [GeV/c].
Definition: ECLChargedPIDDataAnalysisValidationModule.h:160
Belle2::Const::kaon
static const ChargedStable kaon
charged kaon particle
Definition: Const.h:536
Belle2::ECLChargedPIDDataAnalysisValidationModule::m_trkPhi
std::vector< float > m_trkPhi
Track azimuthal angle in [rad].
Definition: ECLChargedPIDDataAnalysisValidationModule.h:178
Belle2::Module
Base class for Modules.
Definition: Module.h:74
Belle2::ECLChargedPIDDataAnalysisValidationModule::m_clusterReg
std::vector< float > m_clusterReg
Cluster ECL region.
Definition: ECLChargedPIDDataAnalysisValidationModule.h:198
Belle2::ECLChargedPIDDataAnalysisValidationModule::m_deltalogl_sig_bkg
std::vector< float > m_deltalogl_sig_bkg
Delta Log-likelihood "signal" vs.
Definition: ECLChargedPIDDataAnalysisValidationModule.h:251
Belle2::merge
std::vector< std::vector< double > > merge(std::vector< std::vector< std::vector< double >>> toMerge)
merge { vector<double> a, vector<double> b} into {a, b}
Definition: tools.h:44
Belle2::Const::pion
static const ChargedStable pion
charged pion particle
Definition: Const.h:535
Belle2::ECLChargedPIDDataAnalysisValidationModule::m_p
std::vector< float > m_p
Track momentum in [GeV/c].
Definition: ECLChargedPIDDataAnalysisValidationModule.h:151
Belle2::ECLChargedPIDDataAnalysisValidationModule::m_trackClusterMatch
std::vector< float > m_trackClusterMatch
Flag for track-cluster matching condition.
Definition: ECLChargedPIDDataAnalysisValidationModule.h:215
Belle2::Const::ParticleSet::find
const ParticleType & find(int pdg) const
Returns particle in set with given PDG code, or invalidParticle if not found.
Definition: Const.h:447
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::ECLChargedPIDDataAnalysisValidationModule::m_logl_sig
std::vector< float > m_logl_sig
Log-likelihood for the "signal" particle hypothesis.
Definition: ECLChargedPIDDataAnalysisValidationModule.h:224
Belle2::ECLChargedPIDDataAnalysisValidationModule::~ECLChargedPIDDataAnalysisValidationModule
virtual ~ECLChargedPIDDataAnalysisValidationModule()
Destructor of the module.
Definition: ECLChargedPIDDataAnalysisValidationModule.cc:49
Belle2::MCParticle::c_Initial
@ c_Initial
bit 5: Particle is initial such as e+ or e- and not going to Geant4
Definition: MCParticle.h:68
Belle2::ECLPidLikelihood
Container for likelihoods with ECL PID (ECLChargedPIDModule)
Definition: ECLPidLikelihood.h:34
Belle2::ECLChargedPIDDataAnalysisValidationModule::m_mergeChargeOfPdgIds
std::vector< unsigned int > m_mergeChargeOfPdgIds
The (unsigned) pdgId list of the charged stable particles for which particle and antiparticle should ...
Definition: ECLChargedPIDDataAnalysisValidationModule.h:105
Belle2::ECLChargedPIDDataAnalysisValidationModule::initialize
virtual void initialize() override
Initializes the module.
Definition: ECLChargedPIDDataAnalysisValidationModule.cc:54
Belle2::ECLChargedPIDDataAnalysisValidationModule::beginRun
virtual void beginRun() override
Called once before a new run begins.
Definition: ECLChargedPIDDataAnalysisValidationModule.cc:120
Belle2::ECLChargedPIDDataAnalysisValidationModule::m_trkTheta
std::vector< float > m_trkTheta
Track polar angle in [rad].
Definition: ECLChargedPIDDataAnalysisValidationModule.h:169
Belle2::Const::ChargedStable
Provides a type-safe way to pass members of the chargedStableSet set.
Definition: Const.h:465
Belle2::ECLChargedPIDDataAnalysisValidationModule::m_saveValidationTree
bool m_saveValidationTree
Save the TTree in the output file alongside the histograms.
Definition: ECLChargedPIDDataAnalysisValidationModule.h:135
Belle2::Const::muon
static const ChargedStable muon
muon particle
Definition: Const.h:534
Belle2::ECLChargedPIDDataAnalysisValidationModule::m_tree
std::vector< TTree * > m_tree
A ROOT::TTree filled with the info to make control plots.
Definition: ECLChargedPIDDataAnalysisValidationModule.h:142
Belle2::Track
Class that bundles various TrackFitResults.
Definition: Track.h:35
Belle2::ECLChargedPIDDataAnalysisValidationModule::m_outputFileName
std::string m_outputFileName
Base name of the output ROOT::TFile.
Definition: ECLChargedPIDDataAnalysisValidationModule.h:129
Belle2::StoreArray< MCParticle >
Belle2::ECLChargedPIDDataAnalysisValidationModule::endRun
virtual void endRun() override
Called once when a run ends.
Definition: ECLChargedPIDDataAnalysisValidationModule.cc:257
Belle2::ECLChargedPIDDataAnalysisValidationModule
This module dumps a tree and a set of histograms of ECL PID-related info used for validation,...
Definition: ECLChargedPIDDataAnalysisValidationModule.h:42
Belle2::ECLChargedPIDDataAnalysisValidationModule::c_PID
static constexpr float c_PID
Defintion of the PID cut threshold to compute the efficiency.
Definition: ECLChargedPIDDataAnalysisValidationModule.h:87
Belle2::ECLChargedPIDDataAnalysisValidationModule::m_clusterPhi
std::vector< float > m_clusterPhi
Cluster azimuthal angle in [rad].
Definition: ECLChargedPIDDataAnalysisValidationModule.h:208
Belle2::ECLChargedPIDDataAnalysisValidationModule::m_inputPdgIdSet
std::set< int > m_inputPdgIdSet
The pdgId set of the charged stable particles of interest.
Definition: ECLChargedPIDDataAnalysisValidationModule.h:116
Belle2::Const::ParticleType::getIndex
int getIndex() const
This particle's index in the associated set.
Definition: Const.h:337
Belle2::MCParticle::c_PrimaryParticle
@ c_PrimaryParticle
bit 0: Particle is primary particle.
Definition: MCParticle.h:58
Belle2::ECLChargedPIDDataAnalysisValidationModule::m_clusterTheta
std::vector< float > m_clusterTheta
Cluster polar angle in [rad].
Definition: ECLChargedPIDDataAnalysisValidationModule.h:188