|
ATLAS Offline Software
|
Go to the documentation of this file.
11 #ifndef FOURMOMUTILS_P4HELPERS_H
12 #define FOURMOMUTILS_P4HELPERS_H
34 inline double deltaPhi(
double phiA,
double phiB) {
39 template <
class Iterator_t,
class Predicate_t>
40 inline void sort(Iterator_t itr, Iterator_t itrEnd, Predicate_t
p) {
44 return sort(itr, itrEnd,
p);
48 template <
class Container_t,
class Predicate_t>
49 inline void sort(Container_t& container, Predicate_t
p) {
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 )
75 * ( mom2 + pz2 ) / ( mom2 - pz2 ) );
79 const double mom1 =
p1.p();
const double pz1 =
p1.pz();
81 dEta = -0.5 *
std::log( ( mom1 - pz1 ) / ( mom1 + pz1 ) ) -
p2.eta();
85 const double mom2 =
p2.p();
const double pz2 =
p2.pz();
87 dEta =
p1.eta() - 0.5 *
std::log( ( mom2 + pz2 ) / ( mom2 - pz2 ) );
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 )
107 dEtaSq =
std::log( ( mom1 - pz1 ) / ( mom1 + pz1 )
108 * ( mom2 + pz2 ) / ( mom2 - pz2 ) );
109 dEtaSq = 0.25 * dEtaSq * dEtaSq;
114 const double mom1 =
p1.p();
const double pz1 =
p1.pz();
115 if ( mom1 + pz1 > 0 )
117 dEtaSq = - 0.5 *
std::log( ( mom1 - pz1 ) / ( mom1 + pz1 ) ) -
p2.eta();
118 dEtaSq = dEtaSq * dEtaSq;
123 const double mom2 =
p2.p();
const double pz2 =
p2.pz();
124 if ( mom2 - pz2 > 0 )
126 dEtaSq =
p1.eta() - 0.5 *
std::log( ( mom2 + pz2 ) / ( mom2 - pz2 ) );
127 dEtaSq = dEtaSq * dEtaSq;
132 dEtaSq =
p1.eta() -
p2.eta();
133 dEtaSq = dEtaSq * dEtaSq;
159 double xp = pA.
px() * pB.
px() ;
163 xyp1 = 1.0 + pA.
py() * pB.
py() / xp ;
169 if ( pA.
py() / pA.
px() > 0 )
190 double cosPhi = ( pA.
px() * pB.
px() + pA.
py() * pB.
py() ) / pA.
pt() / pB.
pt();
191 double phi = acos(cosPhi);
197 dPhiSq = dPhiSq * dPhiSq;
212 const double dEta = p4.
eta() - eta;
224 return sqrt( dEtaSq + dPhiSq );
242 if (
dPhi > dR )
return false;
244 if (
dEta > dR )
return false;
246 if (
deltaR2 > dR*dR )
return false;
253 {
return (pA.
hlv()+pB.
hlv()).m(); }
258 {
return invMass( *pA, *pB ); }
270 {
return invMass( *pA, *pB, *pC, *pD ); }
278 #endif // FOURMOMUTILS_P4HELPERS_H
virtual double py() const =0
y component of momentum
double deltaEtaSq(const I4Momentum &p1, const I4Momentum &p2)
Computes efficiently .
double invMass(const I4Momentum &pA, const I4Momentum &pB)
invariant mass from two I4momentum references
P4Helpers provides static helper functions for kinematic calculation on objects deriving from I4Momen...
double deltaR2(double rapidity1, double phi1, double rapidity2, double phi2)
from bare rapidity,phi
bool isInDeltaR(const I4Momentum &p1, const I4Momentum &p2, double dR)
Check if 2 I4Momentum are in a cone.
virtual double pt() const =0
transverse momentum
virtual CLHEP::HepLorentzVector hlv() const =0
CLHEP HepLorentzVector.
double deltaPhi(double phiA, double phiB)
delta Phi in range [-pi,pi[
double deltaEta(const I4Momentum &p1, const I4Momentum &p2)
Computes efficiently .
bool dPhi(const xAOD::TauJet &tau, const xAOD::TauTrack &track, double &out)
Workaround x86 precision issues for FP inequality comparisons.
virtual double eta() const =0
pseudo rapidity
double deltaR(const I4Momentum &p4, double eta, double phi)
from 1 I4Momentum
virtual double phi() const =0
phi in [-pi,pi[
virtual Kind kind() const =0
tells what kind of P4XYZT this is
double deltaPhiSq(const I4Momentum &pA, const I4Momentum &pB)
delta Phi squared in range ([-pi,pi[)^2 from two I4momentum references
std::vector< std::string > remainder(const std::vector< std::string > &v1, const std::vector< std::string > &v2)
list of entries in a vector that are not in another
void sort(Iterator_t itr, Iterator_t itrEnd, Predicate_t p)
sort a container according to the given predicate
void sort(Container_t &container, Predicate_t p)
sort a container according to the given predicate
bool dEta(const xAOD::TauJet &tau, const xAOD::TauTrack &track, double &out)
virtual double px() const =0
x component of momentum