8 #include <ecl/modules/eclChargedPIDDataAnalysisExpert/ECLChargedPIDDataAnalysisValidationModule.h>
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>
16 #include <TEfficiency.h>
24 REG_MODULE(ECLChargedPIDDataAnalysisValidation)
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).");
36 std::vector<int> defaultChargedPdgIds;
38 defaultChargedPdgIds.push_back(hypo.getPDGCode());
39 defaultChargedPdgIds.push_back(-hypo.getPDGCode());
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.",
63 B2INFO(
"Initialising ROOT objects...");
74 B2WARNING(
"For (unsigned) hypothesis " << hypo.getPDGCode() <<
", validation plots will be merged for +/- charged particles.");
79 std::string chargedPdgIdStr;
85 if (!isValidChargedPdg(std::abs(chargedPdgId))) {
86 B2FATAL(
"PDG: " << chargedPdgId <<
" in m_inputPdgIdSet is not that of a valid particle in Const::chargedStableSet! Aborting...");
96 auto chargedIdx = chargedHypo.
getIndex();
98 if (chargedPdgId > 0) {
99 chargedPdgIdStr = std::to_string(chargedPdgId);
101 chargedPdgIdStr =
"anti" + std::to_string(std::abs(chargedPdgId));
108 m_outputFile[chargedIdx] =
new TFile(fname.c_str(),
"RECREATE");
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");
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");
122 m_tree[chargedIdx]->Branch(
"pids_glob", &m_pids_glob[chargedIdx]);
143 auto chargedIdx = chargedHypo.
getIndex();
145 if (chargedPdgId < 0) {
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();
159 m_logl_sig[chargedIdx] = std::numeric_limits<float>::quiet_NaN();
160 m_logl_bkg[chargedIdx] = std::numeric_limits<float>::quiet_NaN();
163 m_pids_glob[chargedIdx][chargedStable.getIndex()] = std::numeric_limits<float>::quiet_NaN();
166 for (
const auto& particle : m_MCParticles) {
176 if (std::abs(particle.getPDG()) != std::abs(chargedPdgId))
continue;
179 if (particle.getPDG() != chargedPdgId)
continue;
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();
195 if (itrack_max < 0)
continue;
197 const auto track = particle.getRelationsFrom<
Track>()[itrack_max];
198 const auto fitRes = track->getTrackFitResultWithClosestMass(
Const::pion);
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();
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];
211 if (!eclCluster->isTrack())
continue;
212 icluster_match = icluster;
216 if (icluster_match < 0) {
221 const auto eclCluster = eclClusters[icluster_match];
225 m_clusterReg[chargedIdx] = eclCluster->getDetectorRegion();
237 double lh_sig = eclLikelihood->getLikelihood(chargedStableSig);
238 double lh_bkg = eclLikelihood->getLikelihood(chargedStableBkg);
247 lh_all += eclLikelihood->getLikelihood(chargedStable);
250 m_pids_glob[chargedIdx][chargedStable.getIndex()] = eclLikelihood->getLikelihood(chargedStable) / lh_all;
255 m_tree[chargedIdx]->Fill();
275 if (mergeCharge && chargedPdgId <= 0)
continue;
278 const auto chargeSign =
static_cast<int>(chargedPdgId / std::abs(chargedPdgId));
281 const auto chargedStableSig = chargedStableSample;
288 auto chargedSampleIdx = (chargeSign > 0) ? chargedStableSig.
getIndex() : chargedStableSig.getIndex() +
293 auto pdgIdDesc = (mergeCharge) ? std::to_string(chargedPdgId) +
" and -" + std::to_string(chargedPdgId) : std::to_string(
297 TNamed(
"Description", TString::Format(
"ECL Charged PID control plots for charged stable particles/antiparticles ; Sample PDG = %s",
298 pdgIdDesc.c_str()).Data()).Write();
301 dumpPIDVars(
m_tree[chargedSampleIdx], chargedStableSig, chargeSign, chargedStableBkg, mergeCharge);
303 dumpPIDEfficiencyFakeRate(
m_tree[chargedSampleIdx], chargedStableSample, chargeSign, chargedStableSig, mergeCharge);
306 dumpPIDEfficiencyFakeRate(
m_tree[chargedSampleIdx], chargedStableSample, chargeSign,
Const::electron, mergeCharge);
307 dumpPIDEfficiencyFakeRate(
m_tree[chargedSampleIdx], chargedStableSample, chargeSign,
Const::muon, mergeCharge);
310 dumpTrkClusMatchingEfficiency(
m_tree[chargedSampleIdx], chargedStableSample, chargeSign, mergeCharge);
314 m_tree[chargedSampleIdx]->Write();
322 void ECLChargedPIDDataAnalysisValidationModule::dumpPIDVars(TTree* sampleTree,
const Const::ChargedStable& sigHypo,
328 const auto sigHypoIdx = sigHypo.
getIndex();
329 const auto sigHypoPdgId = sigHypo.
getPDGCode();
332 const auto bkgHypoPdgId = bkgHypo.
getPDGCode();
336 TString pidSigBranch = TString::Format(
"pids_glob[%i]", sigHypoIdx);
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());
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());
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());
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");
361 paintUnderOverflow(h_pid);
362 paintUnderOverflow(h_deltalogl);
363 paintUnderOverflow(h_trkclusmatch);
365 h_pid->SetOption(
"HIST");
366 h_deltalogl->SetOption(
"HIST");
367 h_trkclusmatch->SetOption(
"HIST");
370 std::string metaopts(
"pvalue-warn=0.1,pvalue-error=0.01");
371 std::string shifteropt(
"");
374 shifteropt =
"shifter,";
377 auto pdgIdDesc = (!mergeSigCharge) ? std::to_string(sigHypoPdgId * sigCharge) : std::to_string(
378 sigHypoPdgId) +
" and -" + std::to_string(sigHypoPdgId);
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.",
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()));
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.",
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()));
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()));
409 h_deltalogl->Write();
410 h_trkclusmatch->Write();
414 delete h_trkclusmatch;
419 void ECLChargedPIDDataAnalysisValidationModule::dumpPIDEfficiencyFakeRate(TTree* sampleTree,
const Const::ChargedStable& sampleHypo,
424 const std::string ratioType = (sampleHypo == sigHypo) ?
"Efficiency" :
"FakeRate";
427 const int sampleHypoPdgId = sampleHypo.
getPDGCode() * sampleCharge;
430 const auto sigHypoIdx = sigHypo.
getIndex();
431 const auto sigHypoPdgId = sigHypo.
getPDGCode();
435 TString pidSigCut = TString::Format(
"pids_glob[%i] > %f", sigHypoIdx,
c_PID);
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);
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());
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);
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);
460 sampleTree->Project(h_p_N_name.Data(),
"p", pidSigCut.Data());
461 sampleTree->Project(h_p_D_name.Data(),
"p");
463 sampleTree->Project(h_th_N_name.Data(),
"clusterTheta", pidSigCut.Data());
464 sampleTree->Project(h_th_D_name.Data(),
"clusterTheta");
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);
471 sampleTree->Project(h_phi_N_name.Data(),
"clusterPhi", pidSigCut.Data());
472 sampleTree->Project(h_phi_D_name.Data(),
"clusterPhi");
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());
482 std::string metaopts(
"pvalue-warn=0.01,pvalue-error=0.001,nostats");
483 std::string shifteropt(
"");
486 shifteropt =
"shifter,";
489 auto pdgIdDesc = (!mergeSampleCharge) ? std::to_string(sampleHypoPdgId) : std::to_string(std::abs(
490 sampleHypoPdgId)) +
" and -" + std::to_string(std::abs(sampleHypoPdgId));
492 if (TEfficiency::CheckConsistency(*h_p_N, *h_p_D)) {
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());
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();
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}$.",
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()));
513 t_pid_glob_ratio_p->Write();
515 delete t_pid_glob_ratio_p;
518 if (TEfficiency::CheckConsistency(*h_th_N, *h_th_D)) {
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());
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();
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}$.",
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()));
539 t_pid_glob_ratio_th->Write();
541 delete t_pid_glob_ratio_th;
543 if (TEfficiency::CheckConsistency(*h_eclreg_N, *h_eclreg_D)) {
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());
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();
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).",
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()));
564 t_pid_glob_ratio_eclreg->Write();
566 delete t_pid_glob_ratio_eclreg;
569 if (TEfficiency::CheckConsistency(*h_phi_N, *h_phi_D)) {
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());
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();
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}$.",
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()));
590 t_pid_glob_ratio_phi->Write();
592 delete t_pid_glob_ratio_phi;
607 void ECLChargedPIDDataAnalysisValidationModule::dumpTrkClusMatchingEfficiency(TTree* tree,
const Const::ChargedStable& sampleHypo,
608 const int sampleCharge,
bool mergeSampleCharge)
612 const auto sampleHypoPdgId = sampleHypo.
getPDGCode();
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);
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());
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);
630 TString match_cut_N(
"trackClusterMatch == 1");
631 TString match_cut_D(
"trackClusterMatch >= 0");
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());
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());
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());
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);
651 std::string metaopts(
"pvalue-warn=0.01,pvalue-error=0.001,nostats");
652 std::string shifteropt(
"");
655 shifteropt =
"shifter,";
658 auto pdgIdDesc = (!mergeSampleCharge) ? std::to_string(sampleHypoPdgId * sampleCharge) : std::to_string(
659 sampleHypoPdgId) +
" and -" + std::to_string(sampleHypoPdgId);
661 if (TEfficiency::CheckConsistency(*h_pt_N, *h_pt_D)) {
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());
668 t_match_eff_pt->SetConfidenceLevel(0.683);
669 t_match_eff_pt->SetStatisticOption(TEfficiency::kBUniform);
670 t_match_eff_pt->SetPosteriorMode();
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()));
680 t_match_eff_pt->Write();
682 delete t_match_eff_pt;
685 if (TEfficiency::CheckConsistency(*h_th_N, *h_th_D)) {
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());
692 t_match_eff_th->SetConfidenceLevel(0.683);
693 t_match_eff_th->SetStatisticOption(TEfficiency::kBUniform);
694 t_match_eff_th->SetPosteriorMode();
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()));
704 t_match_eff_th->Write();
706 delete t_match_eff_th;
709 if (TEfficiency::CheckConsistency(*h_phi_N, *h_phi_D)) {
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());
715 t_match_eff_phi->SetConfidenceLevel(0.683);
716 t_match_eff_phi->SetStatisticOption(TEfficiency::kBUniform);
717 t_match_eff_phi->SetPosteriorMode();
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()));
727 t_match_eff_phi->Write();
729 delete t_match_eff_phi;
742 void ECLChargedPIDDataAnalysisValidationModule::paintUnderOverflow(TH1F* h)
745 auto nentries = h->GetEntries();
746 auto nbins_vis = h->GetNbinsX();
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);
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);
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));
767 h->SetEntries(nentries);
Provides a type-safe way to pass members of the chargedStableSet set.
static const unsigned int c_SetSize
Number of elements (for use in array bounds etc.)
const ParticleType & find(int pdg) const
Returns particle in set with given PDG code, or invalidParticle if not found.
int getPDGCode() const
PDG code.
int getIndex() const
This particle's index in the associated set.
static const ChargedStable muon
muon particle
static const ParticleSet chargedStableSet
set of charged stable particles
static const ChargedStable pion
charged pion particle
static const ChargedStable kaon
charged kaon particle
static const ChargedStable electron
electron particle
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 initialize() override
Initializes the module.
virtual void event() override
Called once for each event.
std::vector< float > m_p
Track momentum in [GeV/c].
virtual void endRun() override
Called once when a run ends.
std::vector< float > m_deltalogl_sig_bkg
Delta Log-likelihood "signal" vs.
virtual void terminate() override
Termination action.
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 ...
virtual ~ECLChargedPIDDataAnalysisValidationModule()
Destructor of the module.
std::vector< float > m_clusterReg
Cluster ECL region.
std::vector< float > m_trkPhi
Track azimuthal angle in [rad].
std::vector< float > m_trkTheta
Track polar 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.
@ 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
@ c_PrimaryParticle
bit 0: Particle is primary particle.
@ c_IsVirtual
bit 4: Particle is virtual and not going to Geant4.
Class that bundles various TrackFitResults.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
std::vector< std::vector< double > > merge(std::vector< std::vector< std::vector< double >>> toMerge)
merge { vector<double> a, vector<double> b} into {a, b}
Abstract base class for different kinds of events.