Belle II Software development
PIDDetectorWeights Class Reference

Class for handling the PID weights per detector, used to calculate the track helix isolation score per particle. More...

#include <PIDDetectorWeights.h>

Inheritance diagram for PIDDetectorWeights:

Classes

class  WeightsTable
 Nested class acting as a container the per-detector weights. More...
 

Public Member Functions

 PIDDetectorWeights ()
 Default constructor, necessary for ROOT to stream the object.
 
 ~PIDDetectorWeights ()
 Destructor.
 
 PIDDetectorWeights (const std::string &weightsCSVFileName)
 Constructor from CSV file of weights.
 
 PIDDetectorWeights (const std::string &treeName, const std::string &weightsROOTFileName)
 Constructor from ROOT file w/ TTree of weights.
 
ROOT::RDataFrame getWeightsRDF () const
 Get the RDataFrame of detector weights.
 
double getWeight (Const::ChargedStable hypo, Const::EDetector det, double p, double theta) const
 Lookup the weight from the internal map structures.
 

Private Member Functions

void fillWeightsTablePerHypoFromRDF ()
 Fill the internal weights container class per particle hypo, based on the content of the RDataFrame.
 
 ClassDef (PIDDetectorWeights, 4)
 ClassDef as this is a TObject.
 

Private Attributes

std::map< int, WeightsTablem_weightsTablePerHypo
 Map containing a WeightsTable object per particle hypo.
 
ROOT::RDataFrame m_weightsRDataFrame = ROOT::RDataFrame(1)
 The RDataFrame containing the detector weights per particle hypo, per phase space bin.
 
std::map< std::string, std::string > m_weightNames
 The names of the per-detector weight columns in the RDataFrame.
 

Detailed Description

Class for handling the PID weights per detector, used to calculate the track helix isolation score per particle.

Definition at line 31 of file PIDDetectorWeights.h.

Constructor & Destructor Documentation

◆ PIDDetectorWeights() [1/3]

PIDDetectorWeights ( )
inline

Default constructor, necessary for ROOT to stream the object.

Definition at line 75 of file PIDDetectorWeights.h.

75{};

◆ ~PIDDetectorWeights()

~PIDDetectorWeights ( )
inline

Destructor.

Definition at line 80 of file PIDDetectorWeights.h.

80{};

◆ PIDDetectorWeights() [2/3]

PIDDetectorWeights ( const std::string &  weightsCSVFileName)
inline

Constructor from CSV file of weights.

NB: please ensure all numeric types are stored as double in the CSV!

Parameters
weightsCSVFileNamethe path to the CSV file containing the detector weights per std charged particle hypothesis, in bins of p and theta.

Definition at line 88 of file PIDDetectorWeights.h.

89 {
90 m_weightsRDataFrame = ROOT::RDF::MakeCsvDataFrame(weightsCSVFileName);
92 };
void fillWeightsTablePerHypoFromRDF()
Fill the internal weights container class per particle hypo, based on the content of the RDataFrame.
ROOT::RDataFrame m_weightsRDataFrame
The RDataFrame containing the detector weights per particle hypo, per phase space bin.

◆ PIDDetectorWeights() [3/3]

PIDDetectorWeights ( const std::string &  treeName,
const std::string &  weightsROOTFileName 
)
inline

Constructor from ROOT file w/ TTree of weights.

Parameters
treeNamethe name of the TTree with detector weights per std charged particle hypothesis, in bins of p and theta.
weightsROOTFileNamethe path to the ROOT file containing the TTree.

Definition at line 99 of file PIDDetectorWeights.h.

99 :
100 m_weightsRDataFrame(treeName, weightsROOTFileName)
101 {
103 };

Member Function Documentation

◆ fillWeightsTablePerHypoFromRDF()

void fillWeightsTablePerHypoFromRDF ( )
private

Fill the internal weights container class per particle hypo, based on the content of the RDataFrame.

Definition at line 14 of file PIDDetectorWeights.cc.

15{
16
17 auto dummyBinEdges = std::set<double>();
18
19 // Flatten the RDataFrame content into the tabular structure for each std charged hypothesis.
20 for (const auto& hypo : Const::chargedStableSet) {
21
22 WeightsTable weightsTable;
23
24 auto cut_pdgId = [&](double pdg) { return (pdg == hypo.getPDGCode()); };
25 auto filteredRDF = m_weightsRDataFrame.Filter(cut_pdgId, {"pdgId"});
26
27 if (*filteredRDF.Count()) {
28
29 // Get the lower and upper bin edges from the RDF as C++ vectors
30 auto pMinEdges = filteredRDF.Take<double>("p_min").GetValue();
31 auto pMaxEdges = filteredRDF.Take<double>("p_max").GetValue();
32 // Get the bin indexes.
33 auto pBinIdxs = filteredRDF.Take<double>("p_bin_idx").GetValue();
34 // Convert to a set and take the union. This way duplicates are removed.
35 std::set<double> pMinEdges_set(std::begin(pMinEdges), std::end(pMinEdges));
36 std::set<double> pMaxEdges_set(std::begin(pMaxEdges), std::end(pMaxEdges));
37 std::set<double> pBinEdges;
38 std::set_union(std::begin(pMinEdges_set), std::end(pMinEdges_set),
39 std::begin(pMaxEdges_set), std::end(pMaxEdges_set),
40 std::inserter(pBinEdges, std::begin(pBinEdges)));
41 // Store the set of bin edges into the internal tabular structure for this particle hypothesis.
42 weightsTable.m_pBinEdges = pBinEdges;
43 // Store the nr. of p bins. Used for bin idx linearisation.
44 weightsTable.m_nPBins = pBinEdges.size() - 1;
45 // Check thqt the min, max bin edges have no gaps.
46 // Revert to a std::vector for easy index-based iteration.
47 std::vector<double> pMinEdgesVec(std::begin(pMinEdges_set), std::end(pMinEdges_set));
48 std::vector<double> pMaxEdgesVec(std::begin(pMaxEdges_set), std::end(pMaxEdges_set));
49 for (unsigned int iEdgeIdx(0); iEdgeIdx < pMaxEdgesVec.size() - 1; ++iEdgeIdx) {
50 if (pMaxEdgesVec[iEdgeIdx] != pMinEdgesVec[iEdgeIdx + 1]) {
51 B2FATAL("The p bins in the weights table are not contiguous. Please check the input CSV file.");
52 }
53 }
54
55 auto thetaMinEdges = filteredRDF.Take<double>("theta_min").GetValue();
56 auto thetaMaxEdges = filteredRDF.Take<double>("theta_max").GetValue();
57 auto thetaBinIdxs = filteredRDF.Take<double>("theta_bin_idx").GetValue();
58 std::set<double> thetaMinEdges_set(std::begin(thetaMinEdges), std::end(thetaMinEdges));
59 std::set<double> thetaMaxEdges_set(std::begin(thetaMaxEdges), std::end(thetaMaxEdges));
60 std::set<double> thetaBinEdges;
61 std::set_union(std::begin(thetaMinEdges_set), std::end(thetaMinEdges_set),
62 std::begin(thetaMaxEdges_set), std::end(thetaMaxEdges_set),
63 std::inserter(thetaBinEdges, std::begin(thetaBinEdges)));
64 weightsTable.m_thetaBinEdges = thetaBinEdges;
65 weightsTable.m_nThetaBins = thetaBinEdges.size() - 1;
66 std::vector<double> thetaMinEdgesVec(std::begin(thetaMinEdges_set), std::end(thetaMinEdges_set));
67 std::vector<double> thetaMaxEdgesVec(std::begin(thetaMaxEdges_set), std::end(thetaMaxEdges_set));
68 for (unsigned int iEdge(0); iEdge < thetaMaxEdgesVec.size() - 1; ++iEdge) {
69 if (thetaMaxEdgesVec[iEdge] != thetaMinEdgesVec[iEdge + 1]) {
70 B2FATAL("The theta bins in the weights table are not contiguous. Please check the input CSV file.");
71 }
72 }
73
74 // Get the p, theta bin index columns and zip them together.
75 std::vector<std::tuple<double, double>> pAndThetaIdxs;
76 std::transform(std::begin(pBinIdxs), std::end(pBinIdxs),
77 std::begin(thetaBinIdxs),
78 std::back_inserter(pAndThetaIdxs), [](const auto & pIdx, const auto & thetaIdx) { return std::make_tuple(pIdx, thetaIdx); });
79
80 // Fill a std::unordered_map with the linearised (p, theta) bin index as key and the vector index
81 // (i.e. the filtered RDataFrame row idx) as value for fast lookup.
82 unsigned int iRow(0);
83 for (const auto& tup : pAndThetaIdxs) {
84 double linBinIdx = (std::get<0>(tup) - 1.0) + (std::get<1>(tup) - 1.0) * weightsTable.m_nPBins;
85 weightsTable.m_linBinIdxsToRowIdxs.insert(std::make_pair(linBinIdx, iRow));
86 ++iRow;
87 }
88
89 // Store the vector (column) of weights per detector.
90 for (const Const::EDetector& det : Const::PIDDetectorSet::set()) {
91 auto detName = Const::parseDetectors(det);
92 auto colName = "ablat_s_" + detName;
93 auto weights = filteredRDF.Take<double>(colName.c_str()).GetValue();
94 weightsTable.m_weightsPerDet.insert(std::make_pair(detName, weights));
95 }
96
97 } else {
98
99 B2WARNING("Couldn't find detector weights in input ROOT::RDataFrame for std charged particle hypothesis: " << hypo.getPDGCode());
100 weightsTable.m_isEmpty = true;
101
102 }
103
104 m_weightsTablePerHypo.insert(std::make_pair(hypo.getPDGCode(), weightsTable));
105
106 };
107
108}
static DetectorSet set()
Accessor for the set of valid detector IDs.
Definition: Const.h:333
static const ParticleSet chargedStableSet
set of charged stable particles
Definition: Const.h:618
EDetector
Enum for identifying the detector components (detector and subdetector).
Definition: Const.h:42
static std::string parseDetectors(EDetector det)
Converts Const::EDetector object to string.
Definition: UnitConst.cc:161
std::map< int, WeightsTable > m_weightsTablePerHypo
Map containing a WeightsTable object per particle hypo.

◆ getWeight()

double getWeight ( Const::ChargedStable  hypo,
Const::EDetector  det,
double  p,
double  theta 
) const

Lookup the weight from the internal map structures.

Parameters
hypothe input std charged particle.
detthe input PID detector.
pthe particle momentum in [GeV/c].
thetathe particle polar angle in [rad].

Definition at line 111 of file PIDDetectorWeights.cc.

112{
113
114 const auto weightsTable = &m_weightsTablePerHypo.at(hypo.getPDGCode());
115
116 if (weightsTable->m_isEmpty) {
117 // Weights not available for this pdgId.
118 return std::numeric_limits<float>::quiet_NaN();
119 }
120
121 // Get the bin index of the input p value.
122 auto pBinUpperEdge_it = std::lower_bound(std::begin(weightsTable->m_pBinEdges),
123 std::end(weightsTable->m_pBinEdges),
124 p);
125 auto pBinIdx = std::distance(std::begin(weightsTable->m_pBinEdges), pBinUpperEdge_it);
126 // Get the bin index of the input theta value.
127 auto thetaBinUpperEdge_it = std::lower_bound(std::begin(weightsTable->m_thetaBinEdges),
128 std::end(weightsTable->m_thetaBinEdges),
129 theta);
130 auto thetaBinIdx = std::distance(std::begin(weightsTable->m_thetaBinEdges), thetaBinUpperEdge_it);
131
132 // Get the linearised bin index to look up.
133 double linBinIdx = (pBinIdx - 1.0) + (thetaBinIdx - 1.0) * weightsTable->m_nPBins;
134
135 if (!weightsTable->m_linBinIdxsToRowIdxs.count(linBinIdx)) {
136 // Out-of-bin range p or theta
137 B2DEBUG(11, "\n"
138 << "p = " << p << " [GeV/c], theta = " << theta << " [rad].\n"
139 << "Bin indexes: (" << pBinIdx << ", " << thetaBinIdx << ").\n"
140 << "Either input value is outside of bin range for PID detector weights:\n"
141 << "p : [" << *weightsTable->m_pBinEdges.begin() << ", " << *weightsTable->m_pBinEdges.rbegin() << "],\n"
142 << "theta: [" << *weightsTable->m_thetaBinEdges.begin() << ", " << *weightsTable->m_thetaBinEdges.rbegin() << "]");
143 return std::numeric_limits<float>::quiet_NaN();
144 }
145
146 auto rowIdx = weightsTable->m_linBinIdxsToRowIdxs.at(linBinIdx);
147
148 return weightsTable->m_weightsPerDet.at(Const::parseDetectors(det)).at(rowIdx);
149
150};
int getPDGCode() const
PDG code.
Definition: Const.h:473

◆ getWeightsRDF()

ROOT::RDataFrame getWeightsRDF ( ) const
inline

Get the RDataFrame of detector weights.

To be used for testing/debugging only when creating the payload.

Definition at line 109 of file PIDDetectorWeights.h.

110 {
111 return m_weightsRDataFrame;
112 };

Member Data Documentation

◆ m_weightNames

std::map<std::string, std::string> m_weightNames
private
Initial value:
= {
{Const::parseDetectors(Const::CDC), "ablat_s_CDC"},
{Const::parseDetectors(Const::TOP), "ablat_s_TOP"},
{Const::parseDetectors(Const::ARICH), "ablat_s_ARICH"},
{Const::parseDetectors(Const::ECL), "ablat_s_ECL"},
{Const::parseDetectors(Const::KLM), "ablat_s_KLM"}
}

The names of the per-detector weight columns in the RDataFrame.

Definition at line 141 of file PIDDetectorWeights.h.

◆ m_weightsRDataFrame

ROOT::RDataFrame m_weightsRDataFrame = ROOT::RDataFrame(1)
private

The RDataFrame containing the detector weights per particle hypo, per phase space bin.

Using an RDataFrame is convenient when creating the payload from an input CSV file. However, an RDataFrame cannot be streamed in the persistent ROOT file, thus we must prevent this from happening.

Definition at line 136 of file PIDDetectorWeights.h.

◆ m_weightsTablePerHypo

std::map<int, WeightsTable> m_weightsTablePerHypo
private

Map containing a WeightsTable object per particle hypo.

Use to lookup the index corresponding to a (p, theta) pair.

Definition at line 129 of file PIDDetectorWeights.h.


The documentation for this class was generated from the following files: