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