ATLAS Offline Software
Loading...
Searching...
No Matches
xAODP4Helpers.h
Go to the documentation of this file.
1
2
3/*
4 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
5*/
6
7// P4Helpers.h
8// Header file for P4 helper functions
9// Author: karsten.koeneke@cern.ch
11#ifndef FOURMOMUTILS_XAODP4HELPERS_H
12#define FOURMOMUTILS_XAODP4HELPERS_H
13
20
21
22// STL includes
23#include <cmath>
24
25// CxxUtils includes
26#include "CxxUtils/fpcompare.h" // for fpcompare::less and greater
27
28// EDM includes
29#include "xAODBase/IParticle.h"
31
32
33namespace xAOD
34{
35 namespace P4Helpers
36 {
38 inline
39 double deltaRapidity( const xAOD::IParticle& p1, const xAOD::IParticle& p2 )
40 {
41 return p1.rapidity() - p2.rapidity();
42 }
43
45 inline
46 double deltaRapidity( const xAOD::IParticle * const p1, const xAOD::IParticle * const p2 )
47 {
48 return xAOD::P4Helpers::deltaRapidity(*p1, *p2);
49 }
50
51
53 inline
54 double deltaEta( const xAOD::IParticle& p1, const xAOD::IParticle& p2 )
55 {
56 return p1.eta() - p2.eta();
57 }
58
60 inline
61 double deltaEta( const xAOD::IParticle * const p1, const xAOD::IParticle * const p2 )
62 {
63 return xAOD::P4Helpers::deltaEta(*p1, *p2);
64 }
65
66
68 inline
69 double deltaPhi( double phiA, double phiB )
70 {
71 return -remainder( -phiA + phiB, 2*M_PI );
72 }
73
75 inline
76 double deltaPhi( const xAOD::IParticle& p4, const double phi )
77 {
78 return xAOD::P4Helpers::deltaPhi( p4.phi(), phi );
79 }
80
82 inline
83 double deltaPhi( const xAOD::IParticle & pA, const xAOD::IParticle & pB )
84 {
85 return xAOD::P4Helpers::deltaPhi( pA.phi(), pB.phi() );
86 }
87
89 inline
90 double deltaPhi( const xAOD::IParticle * const pA, const xAOD::IParticle * const pB )
91 { return xAOD::P4Helpers::deltaPhi( pA->phi(), pB->phi() ); }
92
93
94
96 inline
97 double deltaPhi( const xAOD::IParticle & pA, const xAOD::MissingET & met )
98 {
99 return xAOD::P4Helpers::deltaPhi( pA.phi(), met.phi() );
100 }
101
103 inline
104 double deltaPhi( const xAOD::IParticle * const pA, const xAOD::MissingET * const met )
105 { return xAOD::P4Helpers::deltaPhi( pA->phi(), met->phi() ); }
106
107
108
110 inline
111 double deltaR2( double rapidity1, double phi1, double rapidity2, double phi2 )
112 {
113 const double dPhi = xAOD::P4Helpers::deltaPhi( phi1, phi2 );
114 const double dRapidity = rapidity1-rapidity2;
115 return dRapidity*dRapidity + dPhi*dPhi;
116 }
117
119 inline
120 double deltaR2( const xAOD::IParticle& p4, double rapidity, double phi, bool useRapidity=true )
121 {
122 if (useRapidity) {
123 return xAOD::P4Helpers::deltaR2(p4.rapidity(),p4.phi(),rapidity,phi);
124 }
125 else {
126 return xAOD::P4Helpers::deltaR2(p4.eta(),p4.phi(),rapidity,phi);
127 }
128 }
129
131 inline
132 double deltaR2( const xAOD::IParticle& pA, const xAOD::IParticle& pB, bool useRapidity=true )
133 {
134 if (useRapidity) {
135 return xAOD::P4Helpers::deltaR2(pA.rapidity(),pA.phi(),pB.rapidity(),pB.phi());
136 }
137 else {
138 return xAOD::P4Helpers::deltaR2(pA.eta(),pA.phi(),pB.eta(),pB.phi());
139 }
140 }
141
143 inline
144 double deltaR2( const xAOD::IParticle * const pA, const xAOD::IParticle * const pB, bool useRapidity=true )
145 { return xAOD::P4Helpers::deltaR2( *pA, *pB, useRapidity ); }
146
147
149 inline
150 double deltaR( double rapidity1, double phi1, double rapidity2, double phi2 )
151 { return std::sqrt( xAOD::P4Helpers::deltaR2( rapidity1, phi1, rapidity2, phi2 ) ); }
152
154 inline
155 double deltaR( const xAOD::IParticle& p4, double rapidity, double phi, bool useRapidity=true )
156 { return std::sqrt( xAOD::P4Helpers::deltaR2( p4, rapidity, phi, useRapidity ) ); }
157
159 inline
160 double deltaR( const xAOD::IParticle& pA, const xAOD::IParticle& pB, bool useRapidity=true )
161 { return std::sqrt( xAOD::P4Helpers::deltaR2( pA, pB, useRapidity ) ); }
162
164 inline
165 double deltaR( const xAOD::IParticle * const pA, const xAOD::IParticle * const pB, bool useRapidity=true )
166 { return std::sqrt( xAOD::P4Helpers::deltaR2( *pA, *pB, useRapidity ) ); }
167
168
169
173 inline
174 bool isInDeltaR( const xAOD::IParticle& p1, const xAOD::IParticle& p2,
175 double dR, bool useRapidity=true )
176 {
177 const double dPhi = std::abs( xAOD::P4Helpers::deltaPhi(p1,p2) ); // in [-pi,pi)
178 if ( CxxUtils::fpcompare::greater(dPhi,dR) ) return false; // <==
179 if (useRapidity) {
180 const double dRapidity = std::abs( xAOD::P4Helpers::deltaRapidity(p1,p2) );
181 if ( CxxUtils::fpcompare::greater(dRapidity,dR) ) return false; // <==
182 const double deltaR2 = dRapidity*dRapidity + dPhi*dPhi;
183 if ( CxxUtils::fpcompare::greater(deltaR2,dR*dR) ) return false; // <==
184 return true;
185 }
186 else {
187 const double dEta = std::abs( xAOD::P4Helpers::deltaEta(p1,p2) );
188 if ( CxxUtils::fpcompare::greater(dEta,dR) ) return false; // <==
189 const double deltaR2 = dEta*dEta + dPhi*dPhi;
190 if ( CxxUtils::fpcompare::greater(deltaR2,dR*dR) ) return false; // <==
191 return true;
192 }
193 }
194
195
196
198 inline
199 double deltaAbsP( const xAOD::IParticle& pA, const xAOD::IParticle& pB )
200 {
201 const double absPA = std::abs( pA.p4().P() );
202 const double absPB = std::abs( pB.p4().P() );
203 return absPA - absPB;
204 }
205
207 inline
208 double deltaAbsP( const xAOD::IParticle * const pA, const xAOD::IParticle * const pB )
209 { return xAOD::P4Helpers::deltaAbsP( *pA, *pB ); }
210
211
212 } //> namespace P4Helpers
213
214} // End: namespace xAOD
215
216
217#endif // FOURMOMUTILS_XAODP4HELPERS_H
#define M_PI
Class providing the definition of the 4-vector interface.
virtual double eta() const =0
The pseudorapidity ( ) of the particle.
virtual FourMom_t p4() const =0
The full 4-momentum of the particle.
virtual double rapidity() const =0
The true rapidity (y) of the particle.
virtual double phi() const =0
The azimuthal angle ( ) of the particle.
std::vector< std::string > remainder(const std::vector< std::string > &v1, const std::vector< std::string > &v2)
Workaround x86 precision issues for FP inequality comparisons.
bool greater(double a, double b)
Compare two FP numbers, working around x87 precision issues.
Definition fpcompare.h:140
P4Helpers provides static helper functions for kinematic calculation on objects deriving from I4Momen...
Definition P4Helpers.h:32
double deltaPhi(double phiA, double phiB)
delta Phi in range [-pi,pi[
double deltaRapidity(const xAOD::IParticle &p1, const xAOD::IParticle &p2)
Computes efficiently .
double deltaAbsP(const xAOD::IParticle &pA, const xAOD::IParticle &pB)
Get the delta-| | between two particles.
bool isInDeltaR(const xAOD::IParticle &p1, const xAOD::IParticle &p2, double dR, bool useRapidity=true)
Check if 2 xAOD::IParticle are in a cone.
double deltaR(double rapidity1, double phi1, double rapidity2, double phi2)
from bare bare rapidity,phi
double deltaEta(const xAOD::IParticle &p1, const xAOD::IParticle &p2)
Computes efficiently .
double deltaR2(double rapidity1, double phi1, double rapidity2, double phi2)
from bare rapidity,phi
ICaloAffectedTool is abstract interface for tools checking if 4 mom is in calo affected region.
MissingET_v1 MissingET
Version control by type defintion.
setSAddress setEtaMS setDirPhiMS setDirZMS setBarrelRadius setEndcapAlpha setEndcapRadius setInterceptInner setEtaMap setEtaBin setIsTgcFailure setDeltaPt deltaPhi