Belle II Software development
ECLChargedPIDDataAnalysisValidationModule.h
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9#pragma once
10
11// Belle 2
12#include <mdst/dataobjects/MCParticle.h>
13
14// ROOT
15#include <TFile.h>
16#include <TH1F.h>
17#include <TTree.h>
18
19// FRAMEWORK
20#include <framework/core/Module.h>
21#include <framework/datastore/StoreArray.h>
22#include <framework/gearbox/Const.h>
23
24#include <set>
25
26namespace Belle2 {
37
38 public:
39
44
49
53 virtual void initialize() override;
54
58 virtual void beginRun() override;
59
63 virtual void event() override;
64
68 virtual void endRun() override;
69
73 virtual void terminate() override;
74
75 private:
76
81 static constexpr float c_PID = 0.5;
82
87 static constexpr unsigned int c_chargedStableHypos = 2 * Const::ChargedStable::c_SetSize;
88
93 std::vector<int> m_inputPdgIdList;
94
99 std::vector<unsigned int> m_mergeChargeOfPdgIds;
100
104 std::map<Const::ChargedStable, bool> m_mergeChargeFlagByHypo;
105
110 std::set<int> m_inputPdgIdSet;
111
117 std::vector<TFile*> m_outputFile = std::vector<TFile*>(c_chargedStableHypos);
118
123 std::string m_outputFileName;
124
130
136 std::vector<TTree*> m_tree = std::vector<TTree*>(c_chargedStableHypos);
137
145 std::vector<float> m_p = std::vector<float>(c_chargedStableHypos);
146
154 std::vector<float> m_pt = std::vector<float>(c_chargedStableHypos);
155
163 std::vector<float> m_trkTheta = std::vector<float>(c_chargedStableHypos);
164
172 std::vector<float> m_trkPhi = std::vector<float>(c_chargedStableHypos);
173
182 std::vector<float> m_clusterTheta = std::vector<float>(c_chargedStableHypos);
183
192 std::vector<float> m_clusterReg = std::vector<float>(c_chargedStableHypos);
193
202 std::vector<float> m_clusterPhi = std::vector<float>(c_chargedStableHypos);
203
209 std::vector<float> m_trackClusterMatch = std::vector<float>(c_chargedStableHypos);
210
218 std::vector<float> m_logl_sig = std::vector<float>(c_chargedStableHypos);
219
229 std::vector<float> m_logl_bkg = std::vector<float>(c_chargedStableHypos);
230
245 std::vector<float> m_deltalogl_sig_bkg = std::vector<float>(c_chargedStableHypos);
246
258 std::vector<std::vector<float>> m_pids_glob = std::vector<std::vector<float>>(c_chargedStableHypos,
259 std::vector<float>(Const::ChargedStable::c_SetSize));
260
265 std::vector<float> m_p_binedges = {0.0, 0.5, 0.75, 1.0, 3.0, 5.0};
266
271 std::vector<float> m_th_binedges = {0.0, 0.2164208, 0.385, 0.561996, 1.13, 1.57, 1.88, 2.2462387, 2.47, 2.7070057, 3.1415926};
272
277
281 void dumpPIDVars(TTree* sampleTree, const Const::ChargedStable& sigHypo, const int sigCharge, const Const::ChargedStable& bkgHypo,
282 bool mergeSigCharge = false);
283
295 void dumpPIDEfficiencyFakeRate(TTree* sampleTree, const Const::ChargedStable& sampleHypo, const int sampleCharge,
296 const Const::ChargedStable& sigHypo, bool mergeSampleCharge = false);
297
307 void dumpTrkClusMatchingEfficiency(TTree* sampleTree, const Const::ChargedStable& sampleHypo, const int sampleCharge,
308 bool mergeSampleCharge = false);
309
313 inline bool isValidChargedPdg(const int pdg) const
314 {
316 }
317
321 void paintUnderOverflow(TH1F* h);
322
323 };
325}
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
static const ParticleSet chargedStableSet
set of charged stable particles
Definition: Const.h:618
static const ParticleType invalidParticle
Invalid particle, used internally.
Definition: Const.h:681
This module dumps a tree and a set of histograms of ECL PID-related info used for validation,...
static constexpr unsigned int c_chargedStableHypos
The maximal number of charged stable particle hypotheses.
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.
std::vector< float > m_p_binedges
Binning w/ variable bin size for track momentum (in [GeV/c]).
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.
Base class for Modules.
Definition: Module.h:72
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
Abstract base class for different kinds of events.