ATLAS Offline Software
Loading...
Searching...
No Matches
TauAxisSetter.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#ifndef XAOD_ANALYSIS
6
7#include "TauAxisSetter.h"
9
10TauAxisSetter::TauAxisSetter(const std::string& name) :
11TauRecToolBase(name) {}
12
13StatusCode TauAxisSetter::execute(xAOD::TauJet& tau) const {
14
15 if (tau.jet() == nullptr) {
16 ATH_MSG_ERROR("Tau jet link is invalid.");
17 return StatusCode::FAILURE;
18 }
19
20 const xAOD::Jet* jetSeed = tau.jet();
21
22 // Barycenter is the sum of cluster p4 in the seed jet
23 TLorentzVector baryCenter;
24
25 xAOD::JetConstituentVector constituents = jetSeed->getConstituents();
26 for (const xAOD::JetConstituent* constituent : constituents) {
27 baryCenter += tauRecTools::GetConstituentP4(*constituent);
28 }
29
30 ATH_MSG_DEBUG("barycenter (eta, phi): " << baryCenter.Eta() << " " << baryCenter.Phi());
31
32 // Detector axis is the total p4 of clusters within m_clusterCone core of the barycenter
33 TLorentzVector tauDetectorAxis;
34
35 for (const xAOD::JetConstituent* constituent : constituents) {
36 TLorentzVector constituentP4 = tauRecTools::GetConstituentP4(*constituent);
37
38 if (baryCenter.DeltaR(constituentP4) > m_clusterCone) continue;
39
40 tauDetectorAxis += constituentP4;
41 }
42
43 if (tauDetectorAxis.Pt() == 0. && !m_doVertexCorrection) {
44 ATH_MSG_DEBUG("this tau candidate does not have any constituent clusters!");
45 return StatusCode::FAILURE;
46 }
47
48 ATH_MSG_DEBUG("detector axis:" << tauDetectorAxis.Pt()<< " " << tauDetectorAxis.Eta() << " " << tauDetectorAxis.Phi() << " " << tauDetectorAxis.E());
49 tau.setP4(tauDetectorAxis.Pt(), tauDetectorAxis.Eta(), tauDetectorAxis.Phi(), tau.m());
50 tau.setP4(xAOD::TauJetParameters::DetectorAxis, tauDetectorAxis.Pt(), tauDetectorAxis.Eta(), tauDetectorAxis.Phi(), tauDetectorAxis.M());
51
52
54 // Tau intermediate axis (corrected for tau vertex)
55 TLorentzVector tauInterAxis;
56
57 const xAOD::Vertex* jetVertex = tauRecTools::getJetVertex(*jetSeed);
58
59 // Redo the vertex correction when tau vertex is different from jet vertex
60 if (jetVertex != tau.vertex()) {
61
62 // If seed jet has a vertex, then tau must have one
63 if (tau.vertex() == nullptr) {
64 ATH_MSG_WARNING("The seed jet has a vertex, while the tau candidate does not. It should not happen.");
65 return StatusCode::FAILURE;
66 }
67
68 const xAOD::Vertex* tauVertex = tau.vertex();
69
70 // Relative position of the tau vertex and jet vertex
71 Amg::Vector3D position = tauVertex->position();
72 if (jetVertex != nullptr) {
73 position -= jetVertex->position();
74 }
75
76 // Barycenter at the tau vertex
77 TLorentzVector baryCenterTauVertex;
78
79 // Loop over the jet constituents, and calculate the barycenter using the four momentum
80 // corrected to point at tau vertex
81 for (const xAOD::JetConstituent* constituent : constituents) {
82 baryCenterTauVertex += getVertexCorrectedP4(*constituent, position);
83 }
84 ATH_MSG_DEBUG("barycenter (eta, phi) at tau vertex: " << baryCenterTauVertex.Eta() << " " << baryCenterTauVertex.Phi());
85
86 // Tau intermediate axis is the four momentum (corrected to point at tau vertex) of clusters
87 // within m_clusterCone of the barycenter
88 for (const xAOD::JetConstituent* constituent : constituents) {
89 TLorentzVector constituentP4 = getVertexCorrectedP4(*constituent, position);
90 if (baryCenterTauVertex.DeltaR(constituentP4) > m_clusterCone) continue;
91
92 tauInterAxis += constituentP4;
93 }
94 }
95 else {
96 tauInterAxis = tauDetectorAxis;
97 }
98
99 if (tauInterAxis.Pt() == 0.) {
100 ATH_MSG_DEBUG("this tau candidate does not have any constituent clusters!");
101 return StatusCode::FAILURE;
102 }
103
104 ATH_MSG_DEBUG("tau axis:" << tauInterAxis.Pt()<< " " << tauInterAxis.Eta() << " " << tauInterAxis.Phi() << " " << tauInterAxis.E() );
105 tau.setP4(tauInterAxis.Pt(), tauInterAxis.Eta(), tauInterAxis.Phi(), tau.m());
106 tau.setP4(xAOD::TauJetParameters::IntermediateAxis, tauInterAxis.Pt(), tauInterAxis.Eta(), tauInterAxis.Phi(), tauInterAxis.M());
107 } // End of m_doVertexCorrection
108
109 return StatusCode::SUCCESS;
110}
111
112
113
115 const Amg::Vector3D& position) const {
116 TLorentzVector vertexCorrectedP4;
117
118 if (constituent.type() == xAOD::Type::CaloCluster) {
119 const xAOD::CaloCluster* cluster = static_cast<const xAOD::CaloCluster*>(constituent.rawConstituent());
120 vertexCorrectedP4 = xAOD::CaloVertexedTopoCluster(*cluster, position).p4();;
121 }
122 else if ( constituent->type() == xAOD::Type::FlowElement ) {
123 const xAOD::FlowElement* fe = static_cast<const xAOD::FlowElement*>( constituent->rawConstituent() );
124 vertexCorrectedP4 = getVertexCorrectedP4(*fe, position);
125 }
126 else {
127 ATH_MSG_WARNING("Seed jet constituent type not supported, will not do vertex correction !");
128 vertexCorrectedP4 = tauRecTools::GetConstituentP4(constituent);
129 }
130
131 return vertexCorrectedP4;
132}
133
135 const Amg::Vector3D& position) const {
136
137 TLorentzVector vertexCorrectedP4;
138 // Only perfrom vertex corretion for neutral FlowElement
139 if (!fe.isCharged()) {
140 TVector3 pos(position.x(), position.y(), position.z());
141 vertexCorrectedP4 = FEHelpers::getVertexCorrectedFourVec(fe,pos);
142 }
143 else {
144 vertexCorrectedP4 = fe.p4();
145 }
146
147 ATH_MSG_DEBUG("Original fe four momentum, pt: " << fe.pt() <<
148 " eta: " << fe.eta() << " phi: " << fe.phi() << " e: " << fe.e());
149 ATH_MSG_DEBUG("Vertex corrected four momentum, pt: " << vertexCorrectedP4.Pt() <<
150 " eta: " << vertexCorrectedP4.Eta() << " phi: " << vertexCorrectedP4.Phi() << " e: " << vertexCorrectedP4.E());
151
152 return vertexCorrectedP4;
153
154}
155
156#endif
#define ATH_MSG_ERROR(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
TLorentzVector getVertexCorrectedP4(const xAOD::JetConstituent &constituent, const Amg::Vector3D &position) const
Get the vertex corrected four momentum.
Gaudi::Property< bool > m_doVertexCorrection
virtual StatusCode execute(xAOD::TauJet &tau) const override
Execution of this tool.
TauAxisSetter(const std::string &name)
Constructor.
Gaudi::Property< double > m_clusterCone
TauRecToolBase(const std::string &name)
virtual FourMom_t p4() const final
The full 4-momentum of the particle.
Evaluate cluster kinematics with a different vertex / signal state.
virtual double pt() const override
virtual double phi() const override
The azimuthal angle ( ) of the particle.
virtual double eta() const override
The pseudorapidity ( ) of the particle.
virtual double e() const override
The total energy of the particle.
virtual FourMom_t p4() const override
The full 4-momentum of the particle.
A vector of jet constituents at the scale used during jet finding.
4-vector of jet constituent at the scale used during jet finding.
Type::ObjectType type() const
The full 4-momentum of the particle.
const IParticle * rawConstituent() const
Access the real underlying IParticle.
JetConstituentVector getConstituents() const
Return a vector of consituents. The object behaves like vector<const IParticle*>. See JetConstituentV...
Definition Jet_v1.cxx:147
const Vertex * vertex() const
void setP4(double pt, double eta, double phi, double m)
Set methods for IParticle values.
virtual double m() const
The invariant mass of the particle.
const Jet * jet() const
const Amg::Vector3D & position() const
Returns the 3-pos.
Eigen::Matrix< double, 3, 1 > Vector3D
TLorentzVector getVertexCorrectedFourVec(const xAOD::FlowElement &fe, const xAOD::Vertex &vertexToCorrectTo)
Definition FEHelpers.cxx:13
const xAOD::Vertex * getJetVertex(const xAOD::Jet &jet)
Return the vertex of jet candidate.
TLorentzVector GetConstituentP4(const xAOD::JetConstituent &constituent)
@ FlowElement
The object is a track-calo-cluster.
Definition ObjectType.h:52
@ CaloCluster
The object is a calorimeter cluster.
Definition ObjectType.h:39
Jet_v1 Jet
Definition of the current "jet version".
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.
FlowElement_v1 FlowElement
Definition of the current "pfo version".
Definition FlowElement.h:16
Vertex_v1 Vertex
Define the latest version of the vertex class.
TauJet_v3 TauJet
Definition of the current "tau version".