ATLAS Offline Software
Loading...
Searching...
No Matches
MvaTESVariableDecorator.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5// local include(s)
8
12
14 : TauRecToolBase(name) {
15}
16
17
19
20 ATH_CHECK(m_aveIntPerXKey.initialize());
23
24 return StatusCode::SUCCESS;
25}
26
27
28
30 // Tell clang to optimize assuming that FP operations may trap.
32
33 int mu = 0;
35 if (!eventInfoDecorHandle.isPresent()) {
36 ATH_MSG_WARNING ( "EventInfo decoration not available! Will set mu=0." );
37 }
38 else {
39 // convert from float to int to ignore peculiar values used in MC
40 mu = static_cast<int>(eventInfoDecorHandle(0));
41 }
42 static const SG::Accessor<float> acc_mu("mu");
43 acc_mu(xTau) = mu;
44
45 if (!m_vertexContainerKey.empty()) {
46 int nVtxPU = 0;
48 if (!vertexInHandle.isValid()) {
49 ATH_MSG_WARNING ("Could not retrieve HiveDataObj with key " << vertexInHandle.key() << ", will set nVtxPU=0.");
50 }
51 else {
52 const xAOD::VertexContainer* vertexContainer = vertexInHandle.cptr();
53 for (const auto *xVertex : *vertexContainer){
54 if (xVertex->vertexType() == xAOD::VxType::PileUp)
55 ++nVtxPU;
56 }
57 }
58 static const SG::Accessor<int> acc_nVtxPU("nVtxPU");
59 acc_nVtxPU(xTau) = nVtxPU;
60 }
61
62 if (!m_eventShapeKey.empty()) {
63 double rho = 0.;
65 if (!eventShape.isValid()) {
66 ATH_MSG_WARNING ("Could not retrieve EventShape with key " << m_eventShapeKey );
67 }
68 else if (!eventShape->getDensity(xAOD::EventShape::Density, rho)) {
69 ATH_MSG_WARNING ("Could not retrieve rho.");
70 }
71 static const SG::Accessor<float> acc_rho("rho");
72 acc_rho(xTau) = static_cast<float>(rho);
73 }
74
75 double center_lambda=0. , first_eng_dens=0. , em_probability=0. , second_lambda=0. ;
76 double mean_center_lambda=0. , mean_first_eng_dens=0. , mean_em_probability=0. , mean_second_lambda=0. ;
77 double mean_presampler_frac=0., lead_cluster_frac=0. , second_cluster_frac=0. , third_cluster_frac=0. ;
78 double clE=0., Etot=0.;
79
80 // approximate Upsilon based on clusters, not PFOs (for online TES)
81 TLorentzVector clusters_EM_P4;
82 clusters_EM_P4.SetPtEtaPhiM(0,0,0,0);
83 TLorentzVector clusters_had_P4;
84 clusters_had_P4.SetPtEtaPhiM(0,0,0,0);
85 TLorentzVector tauIntermediateAxisEM;
86 tauIntermediateAxisEM.SetPtEtaPhiM(0,0,0,0);
87
88 TLorentzVector tauAxis = tauRecTools::getTauAxis(xTau, m_doVertexCorrection);
89
90 // in the trigger, we must ignore the tau vertex, otherwise ptFinalCalib computed at calo-only step and precision step will differ
91 const xAOD::Vertex* vertex = nullptr;
92 if (m_doVertexCorrection) vertex = xTau.vertex();
93
94 // loop over tau clusters
95 std::vector<xAOD::CaloVertexedTopoCluster> vertexedClusters = xTau.vertexedClusters();
96
97 for (const xAOD::CaloVertexedTopoCluster& vertexedCluster : vertexedClusters) {
98 const xAOD::CaloCluster& cluster = vertexedCluster.clust();
99
100 TLorentzVector clusterP4 = m_doVertexCorrection? vertexedCluster.p4() : cluster.p4();
101 if (clusterP4.DeltaR(tauAxis) > 0.2) continue;
102
103 // FIXME: should we use calE for EMTopo clusters ?
104 // what's the energy scale when calculating the cluster momentum
105 clE = cluster.calE();
106 Etot += clE;
107
108 if (clE > lead_cluster_frac) {
109 third_cluster_frac = second_cluster_frac;
110 second_cluster_frac = lead_cluster_frac;
111 lead_cluster_frac = clE;
112 }
113 else if (clE > second_cluster_frac) {
114 third_cluster_frac = second_cluster_frac;
115 second_cluster_frac = clE;
116 }
117 else if (clE > third_cluster_frac) {
118 third_cluster_frac = clE;
119 }
120
122 mean_center_lambda += clE*center_lambda;
123 else ATH_MSG_WARNING("Failed to retrieve moment: CENTER_LAMBDA");
124
126 mean_first_eng_dens += clE*first_eng_dens;
127 else ATH_MSG_WARNING("Failed to retrieve moment: FIRST_ENG_DENS");
128
130 mean_em_probability += clE*em_probability;
131
132 // FIXME: should we use calE for EMTopo clusters ?
133 // what's the energy scale when calculating the cluster momentum
134 if (em_probability>0.5) clusters_EM_P4 += cluster.p4(xAOD::CaloCluster::State::CALIBRATED);
135 else clusters_had_P4 += cluster.p4(xAOD::CaloCluster::State::CALIBRATED);
136 }
137 else ATH_MSG_WARNING("Failed to retrieve moment: EM_PROBABILITY");
138
140 mean_second_lambda += clE*second_lambda;
141 else ATH_MSG_WARNING("Failed to retrieve moment: SECOND_LAMBDA");
142
143 mean_presampler_frac += (cluster.eSample(CaloSampling::PreSamplerB) + cluster.eSample(CaloSampling::PreSamplerE));
144
145 // EM-scale equivalent of IntermediateAxis p4
146 if (vertex) {
147 xAOD::CaloVertexedTopoCluster vertexedClusterEM(cluster, xAOD::CaloCluster::State::UNCALIBRATED, vertex->position());
148 tauIntermediateAxisEM += vertexedClusterEM.p4();
149 }
150 else {
151 tauIntermediateAxisEM += cluster.p4(xAOD::CaloCluster::State::UNCALIBRATED);
152 }
153 }
154
155 // calculate mean values
156 if (Etot>0.) {
157 mean_center_lambda /= Etot;
158 mean_first_eng_dens /= Etot;
159 mean_em_probability /= Etot;
160 mean_second_lambda /= Etot;
161 mean_presampler_frac /= Etot;
162 lead_cluster_frac /= Etot;
163 second_cluster_frac /= Etot;
164 third_cluster_frac /= Etot;
165
166 if (mean_first_eng_dens>0.) mean_first_eng_dens = TMath::Log10(mean_first_eng_dens/Etot);
167 }
168
169 // cluster-based upsilon, ranges from -1 to 1
170 double upsilon_cluster = -2.;
171 if (clusters_had_P4.E()+clusters_EM_P4.E()!=0.)
172 upsilon_cluster = (clusters_had_P4.E()-clusters_EM_P4.E())/(clusters_had_P4.E()+clusters_EM_P4.E());
173
174 xTau.setDetail(xAOD::TauJetParameters::ClustersMeanCenterLambda, static_cast<float>(mean_center_lambda));
175 xTau.setDetail(xAOD::TauJetParameters::ClustersMeanFirstEngDens, static_cast<float>(mean_first_eng_dens));
176 xTau.setDetail(xAOD::TauJetParameters::ClustersMeanEMProbability, static_cast<float>(mean_em_probability));
177 xTau.setDetail(xAOD::TauJetParameters::ClustersMeanSecondLambda, static_cast<float>(mean_second_lambda));
178 xTau.setDetail(xAOD::TauJetParameters::ClustersMeanPresamplerFrac, static_cast<float>(mean_presampler_frac));
179
180 static const SG::Accessor<float> acc_ClusterTotalEnergy("ClusterTotalEnergy");
181 acc_ClusterTotalEnergy(xTau) = static_cast<float>(Etot);
182
183 static const SG::Accessor<float> acc_ptIntermediateAxisEM("ptIntermediateAxisEM");
184 acc_ptIntermediateAxisEM(xTau) = static_cast<float>(tauIntermediateAxisEM.Pt());
185
186 // online-specific, not defined in TauDefs enum
187 static const SG::Accessor<float> acc_LeadClusterFrac("LeadClusterFrac");
188 static const SG::Accessor<float> acc_UpsilonCluster("UpsilonCluster");
189 acc_LeadClusterFrac(xTau) = static_cast<float>(lead_cluster_frac);
190 acc_UpsilonCluster(xTau) = static_cast<float>(upsilon_cluster);
191
192 if (inTrigger()) {
193 // for now only used by trigger, but could be used by offline 0p in the 2022 reprocessing
194 static const SG::Accessor<float> acc_SecondClusterFrac("SecondClusterFrac");
195 static const SG::Accessor<float> acc_ThirdClusterFrac("ThirdClusterFrac");
196 acc_SecondClusterFrac(xTau) = static_cast<float>(second_cluster_frac);
197 acc_ThirdClusterFrac(xTau) = static_cast<float>(third_cluster_frac);
198
199 return StatusCode::SUCCESS;
200 }
201
202 // summing corrected Pi0 PFO energies
203 TLorentzVector Pi0_totalP4;
204 Pi0_totalP4.SetPtEtaPhiM(0,0,0,0);
205
206 for (size_t i=0; i<xTau.nPi0PFOs(); i++){
207 Pi0_totalP4 += xTau.pi0PFO(i)->p4();
208 }
209
210 double Pi0_totalE = Pi0_totalP4.E();
211
212 // summing tau track momenta
213 TLorentzVector charged_totalP4;
214 charged_totalP4.SetPtEtaPhiM(0,0,0,0);
215
216 for (const xAOD::TauTrack* track : xTau.tracks()) {
217 charged_totalP4 += track->p4();
218 }
219
220 double charged_totalE = charged_totalP4.E();
221
222 // calculate relative difference and decorate onto tau
223 double relDiff=0.;
224 if (Pi0_totalE+charged_totalE != 0.){
225 relDiff = (charged_totalE - Pi0_totalE) / (charged_totalE + Pi0_totalE) ;
226 }
227 xTau.setDetail(xAOD::TauJetParameters::PFOEngRelDiff, static_cast<float>(relDiff));
228
229 return StatusCode::SUCCESS;
230}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_WARNING(x)
Handle class for reading a decoration on an object.
Handle class for reading from StoreGate.
virtual StatusCode execute(xAOD::TauJet &xTau) const override
Execute - called for each tau candidate.
SG::ReadHandleKey< xAOD::VertexContainer > m_vertexContainerKey
SG::ReadDecorHandleKey< xAOD::EventInfo > m_aveIntPerXKey
virtual StatusCode initialize() override
Tool initializer.
SG::ReadHandleKey< xAOD::EventShape > m_eventShapeKey
MvaTESVariableDecorator(const std::string &name="MvaTESVariableDecorator")
Gaudi::Property< bool > m_doVertexCorrection
Helper class to provide type-safe access to aux data.
Handle class for reading a decoration on an object.
bool isPresent() const
Is the referenced container present in SG?
virtual bool isValid() override final
Can the handle be successfully dereferenced?
const_pointer_type cptr()
Dereference the pointer.
virtual const std::string & key() const override final
Return the StoreGate ID for the referenced object.
TauRecToolBase(const std::string &name)
bool inTrigger() const
bool retrieveMoment(MomentType type, double &value) const
Retrieve individual moment.
float eSample(const CaloSample sampling) const
@ SECOND_LAMBDA
Second Moment in .
@ EM_PROBABILITY
Classification probability to be em-like.
@ FIRST_ENG_DENS
First Moment in E/V.
@ CENTER_LAMBDA
Shower depth at Cluster Centroid.
flt_t calE() const
Geet Energy in signal state CALIBRATED.
virtual FourMom_t p4() const
The full 4-momentum of the particle.
virtual FourMom_t p4() const final
The full 4-momentum of the particle.
Evaluate cluster kinematics with a different vertex / signal state.
virtual FourMom_t p4() const
The full 4-momentum of the particle.
Definition PFO_v1.cxx:95
std::vector< xAOD::CaloVertexedTopoCluster > vertexedClusters() const
const PFO * pi0PFO(size_t i) const
Get the pointer to a given pi0 PFO associated with this tau.
size_t nPi0PFOs() const
Get the number of pi0 PFO particles associated with this tau.
void setDetail(TauJetParameters::Detail detail, int value)
const Vertex * vertex() const
std::vector< const TauTrack * > tracks(TauJetParameters::TauTrackFlag flag=TauJetParameters::TauTrackFlag::classifiedCharged) const
Get the v<const pointer> to a given tauTrack collection associated with this tau.
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 ...
@ PileUp
Pile-up vertex.
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".
Vertex_v1 Vertex
Define the latest version of the vertex class.
TauTrack_v1 TauTrack
Definition of the current version.
Definition TauTrack.h:16
TauJet_v3 TauJet
Definition of the current "tau version".
Tell the compiler to optimize assuming that FP may trap.
#define CXXUTILS_TRAPPING_FP
Definition trapping_fp.h:24