ATLAS Offline Software
Loading...
Searching...
No Matches
P4Helpers Namespace Reference

P4Helpers provides static helper functions for kinematic calculation on objects deriving from I4Momentum. More...

Functions

double deltaPhi (double phiA, double phiB)
 delta Phi in range [-pi,pi[
template<class Iterator_t, class Predicate_t>
void sort (Iterator_t itr, Iterator_t itrEnd, Predicate_t p)
 sort a container according to the given predicate
template<class Container_t, class Predicate_t>
void sort (Container_t &container, Predicate_t p)
 sort a container according to the given predicate
double deltaEta (const I4Momentum &p1, const I4Momentum &p2)
 Computes efficiently \( \Delta{\eta} \).
double deltaEtaSq (const I4Momentum &p1, const I4Momentum &p2)
 Computes efficiently \( \Delta{\eta}^2 \).
double deltaEta (const I4Momentum *const p1, const I4Momentum *const p2)
 Computes efficiently \( \Delta{\eta} \), pointer args.
double deltaPhi (const I4Momentum &p4, const double phi)
 delta Phi in range [-pi,pi[ from one I4momentum reference
double deltaPhi (const I4Momentum &pA, const I4Momentum &pB)
 delta Phi in range [-pi,pi[ from two I4momentum references
double deltaPhiSq (const I4Momentum &pA, const I4Momentum &pB)
 delta Phi squared in range ([-pi,pi[)^2 from two I4momentum references
double deltaPhi (const I4Momentum *const pA, const I4Momentum *const pB)
 delta Phi in range [-pi,pi[ from two I4momentum pointers
double deltaR (const I4Momentum &p4, double eta, double phi)
 \( \Delta{R} \) from 1 I4Momentum
double deltaR (const I4Momentum &pA, const I4Momentum &pB)
 delta R from two I4momentum reference
double deltaR (const I4Momentum *const pA, const I4Momentum *const pB)
 delta R from two I4momentum pointers
bool isInDeltaR (const I4Momentum &p1, const I4Momentum &p2, double dR)
 Check if 2 I4Momentum are in a \( \Delta{R} \) cone.
double invMass (const I4Momentum &pA, const I4Momentum &pB)
 invariant mass from two I4momentum references
double invMass (const I4Momentum *const pA, const I4Momentum *const pB)
 invariant mass from two I4momentum pointers
double invMass (const I4Momentum &pA, const I4Momentum &pB, const I4Momentum &pC, const I4Momentum &pD)
 invariant mass from four I4momentum references
double invMass (const I4Momentum *const pA, const I4Momentum *const pB, const I4Momentum *const pC, const I4Momentum *const pD)
 invariant mass from four I4momentum pointers

Detailed Description

P4Helpers provides static helper functions for kinematic calculation on objects deriving from I4Momentum.

Author
David Rousseau rouss.nosp@m.eau@.nosp@m.lal.i.nosp@m.n2p3.nosp@m..fr
Sebastien Binet binet.nosp@m.@cer.nosp@m.n.ch
Tadashi Maeno

Function Documentation

◆ deltaEta() [1/2]

double P4Helpers::deltaEta ( const I4Momentum & p1,
const I4Momentum & p2 )
inline

Computes efficiently \( \Delta{\eta} \).

Definition at line 66 of file P4Helpers.h.

67 {
68 double dEta(999);
69 if ( p1.kind() == I4Momentum::P4PXPYPZE && p2.kind() == I4Momentum::P4PXPYPZE )
70 {
71 const double mom1 = p1.p(); const double pz1 = p1.pz();
72 const double mom2 = p2.p(); const double pz2 = p2.pz();
73 if ( ( mom1 + pz1 ) * ( mom2 - pz2 ) > 0 )
74 dEta = -0.5 * std::log( ( mom1 - pz1 ) / ( mom1 + pz1 )
75 * ( mom2 + pz2 ) / ( mom2 - pz2 ) );
76 }
77 else if ( p1.kind() == I4Momentum::P4PXPYPZE )
78 {
79 const double mom1 = p1.p(); const double pz1 = p1.pz();
80 if ( mom1 + pz1 > 0 )
81 dEta = -0.5 * std::log( ( mom1 - pz1 ) / ( mom1 + pz1 ) ) - p2.eta();
82 }
83 else if ( p2.kind() == I4Momentum::P4PXPYPZE )
84 {
85 const double mom2 = p2.p(); const double pz2 = p2.pz();
86 if ( mom2 - pz2 > 0 )
87 dEta = p1.eta() - 0.5 * std::log( ( mom2 + pz2 ) / ( mom2 - pz2 ) );
88 }
89 else
90 {
91 dEta = p1.eta() - p2.eta();
92 }
93 return dEta;
94 }
bool dEta(const xAOD::TauJet &tau, const xAOD::CaloVertexedTopoCluster &cluster, float &out)

◆ deltaEta() [2/2]

double P4Helpers::deltaEta ( const I4Momentum *const p1,
const I4Momentum *const p2 )
inline

Computes efficiently \( \Delta{\eta} \), pointer args.

Definition at line 140 of file P4Helpers.h.

141 {
142 return deltaEta(*p1, *p2);
143 }
double deltaEta(const I4Momentum &p1, const I4Momentum &p2)
Computes efficiently .
Definition P4Helpers.h:66

◆ deltaEtaSq()

double P4Helpers::deltaEtaSq ( const I4Momentum & p1,
const I4Momentum & p2 )
inline

Computes efficiently \( \Delta{\eta}^2 \).

Definition at line 98 of file P4Helpers.h.

99 {
100 double dEtaSq(999);
101 if ( p1.kind() == I4Momentum::P4PXPYPZE && p2.kind() == I4Momentum::P4PXPYPZE )
102 {
103 const double mom1 = p1.p(); const double pz1 = p1.pz();
104 const double mom2 = p2.p(); const double pz2 = p2.pz();
105 if ( ( mom1 + pz1 ) * ( mom2 - pz2 ) > 0 )
106 {
107 dEtaSq = std::log( ( mom1 - pz1 ) / ( mom1 + pz1 )
108 * ( mom2 + pz2 ) / ( mom2 - pz2 ) );
109 dEtaSq = 0.25 * dEtaSq * dEtaSq;
110 }
111 }
112 else if ( p1.kind() == I4Momentum::P4PXPYPZE )
113 {
114 const double mom1 = p1.p(); const double pz1 = p1.pz();
115 if ( mom1 + pz1 > 0 )
116 {
117 dEtaSq = - 0.5 * std::log( ( mom1 - pz1 ) / ( mom1 + pz1 ) ) - p2.eta();
118 dEtaSq = dEtaSq * dEtaSq;
119 }
120 }
121 else if ( p2.kind() == I4Momentum::P4PXPYPZE )
122 {
123 const double mom2 = p2.p(); const double pz2 = p2.pz();
124 if ( mom2 - pz2 > 0 )
125 {
126 dEtaSq = p1.eta() - 0.5 * std::log( ( mom2 + pz2 ) / ( mom2 - pz2 ) );
127 dEtaSq = dEtaSq * dEtaSq;
128 }
129 }
130 else
131 {
132 dEtaSq = p1.eta() - p2.eta();
133 dEtaSq = dEtaSq * dEtaSq;
134 }
135 return dEtaSq;
136 }

◆ deltaPhi() [1/4]

double P4Helpers::deltaPhi ( const I4Momentum & p4,
const double phi )
inline

delta Phi in range [-pi,pi[ from one I4momentum reference

Definition at line 147 of file P4Helpers.h.

148 {
149 return P4Helpers::deltaPhi( p4.phi(), phi );
150 }
Scalar phi() const
phi method
virtual double phi() const =0
phi in [-pi,pi[
double deltaPhi(double phiA, double phiB)
delta Phi in range [-pi,pi[
Definition P4Helpers.h:34

◆ deltaPhi() [2/4]

double P4Helpers::deltaPhi ( const I4Momentum & pA,
const I4Momentum & pB )
inline

delta Phi in range [-pi,pi[ from two I4momentum references

Definition at line 154 of file P4Helpers.h.

155 {
156 double dPhi(999);
157 if ( pA.kind() == I4Momentum::P4PXPYPZE && pB.kind() == I4Momentum::P4PXPYPZE )
158 {
159 double xp = pA.px() * pB.px() ;
160 double xyp1 = 0;
161 if ( xp != 0 )
162 {
163 xyp1 = 1.0 + pA.py() * pB.py() / xp ;
164 dPhi = 0;
165 if ( xyp1 != 0 )
166 dPhi = atan( ( pA.py() / pA.px() - pB.py() / pB.px() ) / xyp1 );
167 if ( xyp1 < 0 )
168 {
169 if ( pA.py() / pA.px() > 0 )
170 dPhi += M_PI;
171 else
172 dPhi -= M_PI;
173 }
174 }
175 }
176 else
177 {
178 dPhi = remainder( -pA.phi() + pB.phi(), 2*M_PI );
179 }
180 return dPhi;
181 }
#define M_PI
virtual double px() const =0
x component of momentum
virtual Kind kind() const =0
tells what kind of P4XYZT this is
virtual double py() const =0
y component of momentum
std::vector< std::string > remainder(const std::vector< std::string > &v1, const std::vector< std::string > &v2)
bool dPhi(const xAOD::TauJet &tau, const xAOD::CaloVertexedTopoCluster &cluster, float &out)

◆ deltaPhi() [3/4]

double P4Helpers::deltaPhi ( const I4Momentum *const pA,
const I4Momentum *const pB )
inline

delta Phi in range [-pi,pi[ from two I4momentum pointers

Definition at line 204 of file P4Helpers.h.

205 { return deltaPhi( pA->phi(), pB->phi() ); }

◆ deltaPhi() [4/4]

double P4Helpers::deltaPhi ( double phiA,
double phiB )
inline

delta Phi in range [-pi,pi[

Definition at line 34 of file P4Helpers.h.

34 {
35 return -remainder(-phiA + phiB, 2 * M_PI);
36}

◆ deltaPhiSq()

double P4Helpers::deltaPhiSq ( const I4Momentum & pA,
const I4Momentum & pB )
inline

delta Phi squared in range ([-pi,pi[)^2 from two I4momentum references

Definition at line 185 of file P4Helpers.h.

186 {
187 double dPhiSq(999);
188 if ( pA.kind() == I4Momentum::P4PXPYPZE && pB.kind() == I4Momentum::P4PXPYPZE )
189 {
190 double cosPhi = ( pA.px() * pB.px() + pA.py() * pB.py() ) / pA.pt() / pB.pt();
191 double phi = acos(cosPhi);
192 dPhiSq = phi*phi;
193 }
194 else
195 {
196 dPhiSq = remainder( -pA.phi() + pB.phi(), 2*M_PI );
197 dPhiSq = dPhiSq * dPhiSq;
198 }
199 return dPhiSq;
200 }
virtual double pt() const =0
transverse momentum

◆ deltaR() [1/3]

double P4Helpers::deltaR ( const I4Momentum & p4,
double eta,
double phi )
inline

\( \Delta{R} \) from 1 I4Momentum

Definition at line 209 of file P4Helpers.h.

210 {
211 using std::sqrt;
212 const double dEta = p4.eta() - eta;
213 const double dPhi = P4Helpers::deltaPhi( p4, phi );
214 return sqrt( dEta*dEta + dPhi*dPhi );
215 }
Scalar eta() const
pseudorapidity method
virtual double eta() const =0
pseudo rapidity

◆ deltaR() [2/3]

double P4Helpers::deltaR ( const I4Momentum & pA,
const I4Momentum & pB )
inline

delta R from two I4momentum reference

Definition at line 219 of file P4Helpers.h.

220 {
221 using std::sqrt;
222 const double dEtaSq = P4Helpers::deltaEtaSq( pA, pB );
223 const double dPhiSq = P4Helpers::deltaPhiSq( pA, pB );
224 return sqrt( dEtaSq + dPhiSq );
225 }
double deltaEtaSq(const I4Momentum &p1, const I4Momentum &p2)
Computes efficiently .
Definition P4Helpers.h:98
double deltaPhiSq(const I4Momentum &pA, const I4Momentum &pB)
delta Phi squared in range ([-pi,pi[)^2 from two I4momentum references
Definition P4Helpers.h:185

◆ deltaR() [3/3]

double P4Helpers::deltaR ( const I4Momentum *const pA,
const I4Momentum *const pB )
inline

delta R from two I4momentum pointers

Definition at line 229 of file P4Helpers.h.

230 { return P4Helpers::deltaR( *pA, *pB ); }
double deltaR(const I4Momentum &p4, double eta, double phi)
from 1 I4Momentum
Definition P4Helpers.h:209

◆ invMass() [1/4]

double P4Helpers::invMass ( const I4Momentum & pA,
const I4Momentum & pB )
inline

invariant mass from two I4momentum references

Definition at line 252 of file P4Helpers.h.

253 { return (pA.hlv()+pB.hlv()).m(); }
virtual CLHEP::HepLorentzVector hlv() const =0
CLHEP HepLorentzVector.

◆ invMass() [2/4]

double P4Helpers::invMass ( const I4Momentum & pA,
const I4Momentum & pB,
const I4Momentum & pC,
const I4Momentum & pD )
inline

invariant mass from four I4momentum references

Definition at line 262 of file P4Helpers.h.

264 { return (pA.hlv()+pB.hlv()+pC.hlv()+pD.hlv()).m(); }

◆ invMass() [3/4]

double P4Helpers::invMass ( const I4Momentum *const pA,
const I4Momentum *const pB )
inline

invariant mass from two I4momentum pointers

Definition at line 257 of file P4Helpers.h.

258 { return invMass( *pA, *pB ); }
double invMass(const I4Momentum &pA, const I4Momentum &pB)
invariant mass from two I4momentum references
Definition P4Helpers.h:252

◆ invMass() [4/4]

double P4Helpers::invMass ( const I4Momentum *const pA,
const I4Momentum *const pB,
const I4Momentum *const pC,
const I4Momentum *const pD )
inline

invariant mass from four I4momentum pointers

Definition at line 268 of file P4Helpers.h.

270 { return invMass( *pA, *pB, *pC, *pD ); }

◆ isInDeltaR()

bool P4Helpers::isInDeltaR ( const I4Momentum & p1,
const I4Momentum & p2,
double dR )
inline

Check if 2 I4Momentum are in a \( \Delta{R} \) cone.

Parameters
dR[in] \( \Delta{R} \)
Returns
true if they are

Definition at line 236 of file P4Helpers.h.

238 {
239 using std::abs;
240 using std::sqrt;
241 const double dPhi = abs( P4Helpers::deltaPhi(p1,p2) ); // in [-pi,pi)
242 if ( dPhi > dR ) return false; // <==
243 const double dEta = abs( P4Helpers::deltaEta(p1,p2) );
244 if ( dEta > dR ) return false; // <==
245 const double deltaR2 = dEta*dEta + dPhi*dPhi;
246 if ( deltaR2 > dR*dR ) return false; // <==
247 return true;
248 }

◆ sort() [1/2]

template<class Container_t, class Predicate_t>
void P4Helpers::sort ( Container_t & container,
Predicate_t p )
inline

sort a container according to the given predicate

Definition at line 49 of file P4Helpers.h.

49 {
50 return P4Helpers::sort(container.begin(), container.end(), p);
51}
void sort(Iterator_t itr, Iterator_t itrEnd, Predicate_t p)
sort a container according to the given predicate
Definition P4Helpers.h:40

◆ sort() [2/2]

template<class Iterator_t, class Predicate_t>
void P4Helpers::sort ( Iterator_t itr,
Iterator_t itrEnd,
Predicate_t p )
inline

sort a container according to the given predicate

Definition at line 40 of file P4Helpers.h.

40 {
41 // Koenig's look-up at our rescue: handle correctly DataVector's
42 // special sort method. => We inject the std::sort into our namespace
43 using std::sort;
44 return sort(itr, itrEnd, p);
45}
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.