Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
vo007_PhysicsHelpers.C File Reference

Detailed Description

View in nbviewer Open in SWAN In this tutorial we demonstrate RVec helpers for physics computations such as angle differences \(\Delta \phi\), the distance in the \(\eta\)- \(\phi\) plane \(\Delta R\) and the invariant mass.

using namespace ROOT::VecOps;
void vo007_PhysicsHelpers()
{
// The DeltaPhi helper computes the closest angle between angles.
// This means that the resulting value is in the range [-pi, pi].
// Note that the helper also supports to compute the angle difference of an
// RVec and a scalar and two scalars. In addition, the computation of the
// difference and the behaviour at the boundary can be adjusted to radian and
// degrees.
RVec<float> phis = {0.0, 1.0, -0.5, M_PI + 1.0};
auto idx = Combinations(phis, 2);
auto phi1 = Take(phis, idx[0]);
auto phi2 = Take(phis, idx[1]);
auto dphi = DeltaPhi(phi1, phi2);
std::cout << "DeltaPhi(phi1 = " << phi1 << ",\n"
<< " phi2 = " << phi2 << ")\n"
<< " = " << dphi << "\n";
// The DeltaR helper is similar to the DeltaPhi helper and computes the distance
// in the \f$\eta\f$-\f$\phi\f$ plane.
RVec<float> etas = {2.4, -1.5, 1.0, 0.0};
auto eta1 = Take(etas, idx[0]);
auto eta2 = Take(etas, idx[1]);
auto dr = DeltaR(eta1, eta2, phi1, phi2);
std::cout << "\nDeltaR(eta1 = " << eta1 << ",\n"
<< " eta2 = " << eta2 << ",\n"
<< " phi1 = " << phi1 << ",\n"
<< " phi2 = " << phi2 << ")\n"
<< " = " << dr << "\n";
// The InvariantMasses helper computes the invariant mass of a two particle system
// given the properties transverse momentum (pt), rapidity (eta), azimuth (phi)
// and mass.
RVec<float> pt3 = {40, 20, 30};
RVec<float> eta3 = {2.5, 0.5, -1.0};
RVec<float> phi3 = {-0.5, 0.0, 1.0};
RVec<float> mass3 = {10, 10, 10};
RVec<float> pt4 = {20, 10, 40};
RVec<float> eta4 = {0.5, -0.5, 1.0};
RVec<float> phi4 = {0.0, 1.0, -1.0};
RVec<float> mass4 = {2, 2, 2};
auto invMass = InvariantMasses(pt3, eta3, phi3, mass3, pt4, eta4, phi4, mass4);
std::cout << "\nInvariantMass(pt1 = " << pt3 << ",\n"
<< " eta1 = " << eta3 << ",\n"
<< " phi1 = " << phi3 << ",\n"
<< " mass1 = " << mass3 << ",\n"
<< " pt2 = " << pt4 << ",\n"
<< " eta2 = " << eta4 << ",\n"
<< " phi2 = " << phi4 << ",\n"
<< " mass2 = " << mass4 << ")\n"
<< " = " << invMass << "\n";
// The InvariantMass helper also accepts a single set of (pt, eta, phi, mass) vectors. Then,
// the invariant mass of all particles in the collection is computed.
auto invMass2 = InvariantMass(pt3, eta3, phi3, mass3);
std::cout << "\nInvariantMass(pt = " << pt3 << ",\n"
<< " eta = " << eta3 << ",\n"
<< " phi = " << phi3 << ",\n"
<< " mass = " << mass3 << ")\n"
<< " = " << invMass2 << "\n";
}
#define M_PI
Definition Rotated.cxx:105
A "std::vector"-like collection of values implementing handy operation to analyse them.
Definition RVec.hxx:296
Vector1::Scalar DeltaR(const Vector1 &v1, const Vector2 &v2)
Find difference in pseudorapidity (Eta) and Phi betwen two generic vectors The only requirements on t...
Definition VectorUtil.h:110
Vector1::Scalar InvariantMass(const Vector1 &v1, const Vector2 &v2)
return the invariant mass of two LorentzVector The only requirement on the LorentzVector is that they...
Definition VectorUtil.h:246
RVec< T > InvariantMasses(const RVec< T > &pt1, const RVec< T > &eta1, const RVec< T > &phi1, const RVec< T > &mass1, const RVec< T > &pt2, const RVec< T > &eta2, const RVec< T > &phi2, const RVec< T > &mass2)
Return the invariant mass of two particles given the collections of the quantities transverse momentu...
Definition RVec.hxx:1586
RVec< T > Take(const RVec< T > &v, const RVec< typename RVec< T >::size_type > &i)
Return elements of a vector at given indices.
Definition RVec.hxx:1058
DeltaPhi(phi1 = { 0, 0, 0, 1, 1, -0.5 },
phi2 = { 1, -0.5, 4.14159, -0.5, 4.14159, 4.14159 })
= { 1, -0.5, -2.14159, -1.5, 3.14159, -1.64159 }
DeltaR(eta1 = { 2.4, 2.4, 2.4, -1.5, -1.5, 1 },
eta2 = { -1.5, 1, 0, 1, 0, 0 },
phi1 = { 0, 0, 0, 1, 1, -0.5 },
phi2 = { 1, -0.5, 4.14159, -0.5, 4.14159, 4.14159 })
= { 4.02616, 1.48661, 3.21659, 2.91548, 3.48132, 1.92219 }
InvariantMass(pt1 = { 40, 20, 30 },
eta1 = { 2.5, 0.5, -1 },
phi1 = { -0.5, 0, 1 },
mass1 = { 10, 10, 10 },
pt2 = { 20, 10, 40 },
eta2 = { 0.5, -0.5, 1 },
phi2 = { 0, 1, -1 },
mass2 = { 2, 2, 2 })
= { 69.0799, 23.6971, 101.326 }
InvariantMass(pt = { 40, 20, 30 },
eta = { 2.5, 0.5, -1 },
phi = { -0.5, 0, 1 },
mass = { 10, 10, 10 })
= 220.308
Date
March 2019
Author
Stefan Wunsch

Definition in file vo007_PhysicsHelpers.C.