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