ATLAS Offline Software
Loading...
Searching...
No Matches
ConstituentLoaderTauCluster.cxx
Go to the documentation of this file.
1/*
2Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
6
7namespace FlavorTagInference {
8 using FeatureFunc_t = std::function<float(const xAOD::CaloVertexedTopoCluster&, const xAOD::TauJet&)>;
9 using FeatureFuncAsReference_t = std::function<bool(const xAOD::TauJet&, const xAOD::CaloVertexedTopoCluster&, float&)>;
10 ConstituentLoaderTauCluster::ConstituentLoaderTauCluster(const ConstituentsInputConfig& cfg, double max_cluster_dr, bool doVertexCorrection) :
12 m_max_cluster_dr(max_cluster_dr),
13 m_doVertexCorrection(doVertexCorrection)
14 {
15 for (const InputVariableConfig& input_var : cfg.inputs) {
16 m_feature_extractors.push_back(getFeatureExtractor(input_var.name));
17 }
18 }
19
20 std::vector<xAOD::CaloVertexedTopoCluster> ConstituentLoaderTauCluster::getTauClusters(const xAOD::TauJet* tau) const {
21 std::vector<xAOD::CaloVertexedTopoCluster> clusters;
22 TLorentzVector tauAxis = tauRecTools::getTauAxis(*tau, m_doVertexCorrection);
23 for (const xAOD::CaloVertexedTopoCluster& vertexedCluster : tau->vertexedClusters()) {
24 TLorentzVector clusterP4 = vertexedCluster.p4();
25 if (clusterP4.DeltaR(tauAxis) > m_max_cluster_dr) continue;
26 clusters.push_back(vertexedCluster);
27 }
29 // Sort by descending et
31 return lhs.p4().Et() > rhs.p4().Et();
32 };
33 std::sort(clusters.begin(), clusters.end(), et_cmp);
34 } else {
35 // throw
36 throw std::runtime_error("Unsupported sorting order");
37 }
38 // Truncate clusters
39 if (clusters.size() > m_config.max_n_constituents) {
40 clusters.resize(m_config.max_n_constituents, clusters[0]);
41 }
42 return clusters;
43 }
44
45
46 Inputs ConstituentLoaderTauCluster::getFeatures(const xAOD::TauJet* tau, const std::vector<xAOD::CaloVertexedTopoCluster>& tau_clusters) const {
47 std::vector<int64_t> features_dim = {static_cast<int64_t>(tau_clusters.size()), static_cast<int64_t>(m_feature_extractors.size())};
48 std::vector<float> features;
49 features.reserve(tau_clusters.size() * m_feature_extractors.size());
50 for (const auto& cluster : tau_clusters) {
51 for (const auto& extractor : m_feature_extractors) {
52 features.push_back(extractor(cluster, *tau));
53 }
54 }
55 return Inputs{std::move(features), std::move(features_dim)};
56 }
57
58 std::tuple<Inputs, std::vector<const xAOD::IParticle*>> ConstituentLoaderTauCluster::getData(const xAOD::IParticle& i_tau) const {
59 auto tau = dynamic_cast<const xAOD::TauJet*>(&i_tau);
60 std::vector<xAOD::CaloVertexedTopoCluster> sorted_tau_cls = getTauClusters(tau);
61 return std::make_tuple(getFeatures(tau, sorted_tau_cls), std::vector<const xAOD::IParticle*>{} );
62 }
63
65 FeatureFuncAsReference_t func_as_ref = nullptr;
66 try {
67 func_as_ref = m_func_map.at(var_name);
68 } catch (const std::out_of_range &e) {
69 throw std::runtime_error("Variable '" + var_name + "' not defined");
70 }
71 return [func_as_ref](const xAOD::CaloVertexedTopoCluster& cls, const xAOD::TauJet& tau) {
72 float out;
73 bool success = func_as_ref(tau, cls, out);
74 if (!success) {
75 throw std::runtime_error("Error in cluster variable calculation");
76 }
77 return out;
78 };
79 }
80
81 const std::string& ConstituentLoaderTauCluster::getName() const {
82 return m_name;
83 }
85 return m_config.type;
86 }
90 const std::set<std::string>& ConstituentLoaderTauCluster::getUsedRemap() const {
91 return m_used_remap;
92 }
93}
94
95
96namespace TauClusterVars {
98
99bool et_log(const xAOD::TauJet& /*tau*/, const xAOD::CaloVertexedTopoCluster &cluster, float &out) {
100 out = std::log10(cluster.p4().Et());
101 return true;
102}
103
104bool pt_tau_log(const xAOD::TauJet &tau, const xAOD::CaloVertexedTopoCluster& /*cluster*/, float &out) {
105 out = std::log10(std::max(tau.pt(), 1e-6));
106 return true;
107}
108
109bool pt_jetseed_log(const xAOD::TauJet &tau, const xAOD::CaloVertexedTopoCluster& /*cluster*/, float &out) {
110 out = std::log10(tau.ptJetSeed());
111 return true;
112}
113
114bool dEta(const xAOD::TauJet &tau, const xAOD::CaloVertexedTopoCluster &cluster, float &out) {
115 out = cluster.eta() - tau.eta();
116 return true;
117}
118
119bool dPhi(const xAOD::TauJet &tau, const xAOD::CaloVertexedTopoCluster &cluster, float &out) {
120 out = cluster.p4().DeltaPhi(tau.p4());
121 return true;
122}
123
124bool SECOND_R(const xAOD::TauJet& /*tau*/, const xAOD::CaloVertexedTopoCluster &cluster, float &out) {
125 double double_out;
126 const auto success = cluster.clust().retrieveMoment(MomentType::SECOND_R, double_out);
127 out = static_cast<float>(double_out);
128 return success;
129}
130
131bool SECOND_LAMBDA(const xAOD::TauJet& /*tau*/, const xAOD::CaloVertexedTopoCluster &cluster, float &out) {
132 double double_out;
133 const auto success = cluster.clust().retrieveMoment(MomentType::SECOND_LAMBDA, double_out);
134 out = static_cast<float>(double_out);
135 return success;
136}
137
138bool CENTER_LAMBDA(const xAOD::TauJet& /*tau*/, const xAOD::CaloVertexedTopoCluster &cluster, float &out) {
139 double double_out;
140 const auto success = cluster.clust().retrieveMoment(MomentType::CENTER_LAMBDA, double_out);
141 out = static_cast<float>(double_out);
142 return success;
143}
144
146 static const SG::ConstAccessor<float> acc_ClustersMeanSecondLambda("ClustersMeanSecondLambda");
147 float ClustersMeanSecondLambda = acc_ClustersMeanSecondLambda(tau);
148 double secondLambda(0);
149 const auto success = cluster.clust().retrieveMoment(MomentType::SECOND_LAMBDA, secondLambda);
150 out = (ClustersMeanSecondLambda != 0.) ? secondLambda/ClustersMeanSecondLambda : 0.;
151 return success;
152}
153
155 static const SG::ConstAccessor<float> acc_ClustersMeanCenterLambda("ClustersMeanCenterLambda");
156 float ClustersMeanCenterLambda = acc_ClustersMeanCenterLambda(tau);
157 double centerLambda(0);
158 const auto success = cluster.clust().retrieveMoment(MomentType::CENTER_LAMBDA, centerLambda);
159 if (ClustersMeanCenterLambda == 0.){
160 out = 250.;
161 }else {
162 out = centerLambda/ClustersMeanCenterLambda;
163 }
164
165 out = std::min(out, 250.0f);
166
167 return success;
168}
169
170
172 // the ClustersMeanFirstEngDens is the log10 of the energy weighted average of the First_ENG_DENS
173 // divided by ETot to make it dimension-less,
174 // so we need to evaluate the difference of log10(clusterFirstEngDens/clusterTotalEnergy) and the ClustersMeanFirstEngDens
175 double clusterFirstEngDens = 0.0;
176 bool status = cluster.clust().retrieveMoment(MomentType::FIRST_ENG_DENS, clusterFirstEngDens);
177 if (clusterFirstEngDens < 1e-6) clusterFirstEngDens = 1e-6;
178
179 static const SG::ConstAccessor<float> acc_ClusterTotalEnergy("ClusterTotalEnergy");
180 float clusterTotalEnergy = acc_ClusterTotalEnergy(tau);
181 if (clusterTotalEnergy < 1e-6) clusterTotalEnergy = 1e-6;
182
183 static const SG::ConstAccessor<float> acc_ClustersMeanFirstEngDens("ClustersMeanFirstEngDens");
184 float clustersMeanFirstEngDens = acc_ClustersMeanFirstEngDens(tau);
185
186 out = std::log10(clusterFirstEngDens/clusterTotalEnergy) - clustersMeanFirstEngDens;
187
188 return status;
189}
190
191//Extension - Variables for GNTau
192bool e(const xAOD::TauJet& /*tau*/, const xAOD::CaloVertexedTopoCluster &cluster, float &out) {
193 out = cluster.p4().E();
194 return true;
195}
196
197bool et(const xAOD::TauJet& /*tau*/, const xAOD::CaloVertexedTopoCluster &cluster, float &out) {
198 out = cluster.p4().Et();
199 return true;
200}
201
202bool FIRST_ENG_DENS(const xAOD::TauJet& /*tau*/, const xAOD::CaloVertexedTopoCluster &cluster, float &out) {
203 double clusterFirstEngDens = 0.0;
204 bool status = cluster.clust().retrieveMoment(MomentType::FIRST_ENG_DENS, clusterFirstEngDens);
205 out = clusterFirstEngDens;
206 return status;
207}
208
209bool EM_PROBABILITY(const xAOD::TauJet& /*tau*/, const xAOD::CaloVertexedTopoCluster &cluster, float &out) {
210 double clusterEMprob = 0.0;
211 bool status = cluster.clust().retrieveMoment(MomentType::EM_PROBABILITY, clusterEMprob);
212 out = clusterEMprob;
213 return status;
214}
215
216bool CENTER_MAG(const xAOD::TauJet& /*tau*/, const xAOD::CaloVertexedTopoCluster &cluster, float &out) {
217 double clusterCenterMag = 0.0;
218 bool status = cluster.clust().retrieveMoment(MomentType::CENTER_MAG, clusterCenterMag);
219 out = clusterCenterMag;
220 return status;
221}
222
223} // namespace Cluster
std::function< float(const xAOD::CaloVertexedTopoCluster &, const xAOD::TauJet &)> FeatureFunc_t
ConstituentLoaderTauCluster(const ConstituentsInputConfig &cfg, double max_cluster_dr, bool doVertexCorrection)
static const std::unordered_map< std::string, FeatureFuncAsReference_t > m_func_map
const FTagDataDependencyNames & getDependencies() const override
Inputs getFeatures(const xAOD::TauJet *tau, const std::vector< xAOD::CaloVertexedTopoCluster > &tau_clusters) const
std::tuple< Inputs, std::vector< const xAOD::IParticle * > > getData(const xAOD::IParticle &p) const override
const std::set< std::string > & getUsedRemap() const override
FeatureFunc_t getFeatureExtractor(const std::string &var_name) const
std::vector< xAOD::CaloVertexedTopoCluster > getTauClusters(const xAOD::TauJet *tau) const
std::function< bool(const xAOD::TauJet &, const xAOD::CaloVertexedTopoCluster &, float &)> FeatureFuncAsReference_t
const ConstituentsType & getType() const override
Helper class to provide constant type-safe access to aux data.
bool retrieveMoment(MomentType type, double &value) const
Retrieve individual moment.
MomentType
Enums to identify different moments.
virtual FourMom_t p4() const final
The full 4-momentum of the particle.
const CaloCluster & clust() const
Return the cluster being proxied,.
virtual double eta() const final
The pseudorapidity ( ) of the particle.
Evaluate cluster kinematics with a different vertex / signal state.
Class providing the definition of the 4-vector interface.
std::vector< xAOD::CaloVertexedTopoCluster > vertexedClusters() const
virtual FourMom_t p4() const
The full 4-momentum of the particle.
Definition TauJet_v3.cxx:96
virtual double pt() const
The transverse momentum ( ) of the particle.
double ptJetSeed() const
virtual double eta() const
The pseudorapidity ( ) of the particle.
This file contains "getter" functions used for accessing tagger inputs from the EDM.
std::function< float(const xAOD::CaloVertexedTopoCluster &, const xAOD::TauJet &)> FeatureFunc_t
std::function< bool(const xAOD::TauJet &, const xAOD::CaloVertexedTopoCluster &, float &)> FeatureFuncAsReference_t
std::pair< std::vector< float >, std::vector< int64_t > > Inputs
bool et_log(const xAOD::TauJet &, const xAOD::CaloVertexedTopoCluster &cluster, float &out)
bool FIRST_ENG_DENS(const xAOD::TauJet &, const xAOD::CaloVertexedTopoCluster &cluster, float &out)
bool et(const xAOD::TauJet &, const xAOD::CaloVertexedTopoCluster &cluster, float &out)
bool SECOND_LAMBDAOverClustersMeanSecondLambda(const xAOD::TauJet &tau, const xAOD::CaloVertexedTopoCluster &cluster, float &out)
bool CENTER_MAG(const xAOD::TauJet &, const xAOD::CaloVertexedTopoCluster &cluster, float &out)
xAOD::CaloCluster::MomentType MomentType
bool CENTER_LAMBDA(const xAOD::TauJet &, const xAOD::CaloVertexedTopoCluster &cluster, float &out)
bool dPhi(const xAOD::TauJet &tau, const xAOD::CaloVertexedTopoCluster &cluster, float &out)
bool EM_PROBABILITY(const xAOD::TauJet &, const xAOD::CaloVertexedTopoCluster &cluster, float &out)
bool CENTER_LAMBDAOverClustersMeanCenterLambda(const xAOD::TauJet &tau, const xAOD::CaloVertexedTopoCluster &cluster, float &out)
bool FirstEngDensOverClustersMeanFirstEngDens(const xAOD::TauJet &tau, const xAOD::CaloVertexedTopoCluster &cluster, float &out)
bool SECOND_R(const xAOD::TauJet &, const xAOD::CaloVertexedTopoCluster &cluster, float &out)
bool SECOND_LAMBDA(const xAOD::TauJet &, const xAOD::CaloVertexedTopoCluster &cluster, float &out)
bool dEta(const xAOD::TauJet &tau, const xAOD::CaloVertexedTopoCluster &cluster, float &out)
bool pt_tau_log(const xAOD::TauJet &tau, const xAOD::CaloVertexedTopoCluster &, float &out)
bool e(const xAOD::TauJet &, const xAOD::CaloVertexedTopoCluster &cluster, float &out)
bool pt_jetseed_log(const xAOD::TauJet &tau, const xAOD::CaloVertexedTopoCluster &, float &out)
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
TLorentzVector getTauAxis(const xAOD::TauJet &tau, bool doVertexCorrection=true)
Return the four momentum of the tau axis The tau axis is widely used to select clusters and cells in ...
TauJet_v3 TauJet
Definition of the current "tau version".