Logo ROOT   6.18/05
Reference Guide
vo007_PhysicsHelpers.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_vecops
3/// \notebook -nodraw
4/// In this tutorial we demonstrate RVec helpers for physics computations such
5/// as angle differences \f$\Delta \phi\f$, the distance in the \f$\eta\f$-\f$\phi\f$
6/// plane \f$\Delta R\f$ and the invariant mass.
7///
8/// \macro_code
9/// \macro_output
10///
11/// \date March 2019
12/// \author Stefan Wunsch
13
14using namespace ROOT::VecOps;
15
16void vo007_PhysicsHelpers()
17{
18 // The DeltaPhi helper computes the closest angle between angles.
19 // This means that the resulting value is in the range [-pi, pi].
20
21 // Note that the helper also supports to compute the angle difference of an
22 // RVec and a scalar and two scalars. In addition, the computation of the
23 // difference and the behaviour at the boundary can be adjusted to radian and
24 // degrees.
25 RVec<float> phis = {0.0, 1.0, -0.5, M_PI + 1.0};
26 auto idx = Combinations(phis, 2);
27
28 auto phi1 = Take(phis, idx[0]);
29 auto phi2 = Take(phis, idx[1]);
30 auto dphi = DeltaPhi(phi1, phi2);
31
32 std::cout << "DeltaPhi(phi1 = " << phi1 << ",\n"
33 << " phi2 = " << phi2 << ")\n"
34 << " = " << dphi << "\n";
35
36 // The DeltaR helper is similar to the DeltaPhi helper and computes the distance
37 // in the \f$\eta\f$-\f$\phi\f$ plane.
38 RVec<float> etas = {2.4, -1.5, 1.0, 0.0};
39
40 auto eta1 = Take(etas, idx[0]);
41 auto eta2 = Take(etas, idx[1]);
42 auto dr = DeltaR(eta1, eta2, phi1, phi2);
43
44 std::cout << "\nDeltaR(eta1 = " << eta1 << ",\n"
45 << " eta2 = " << eta2 << ",\n"
46 << " phi1 = " << phi1 << ",\n"
47 << " phi2 = " << phi2 << ")\n"
48 << " = " << dr << "\n";
49
50 // The InvariantMasses helper computes the invariant mass of a two particle system
51 // given the properties transverse momentum (pt), rapidity (eta), azimuth (phi)
52 // and mass.
53 RVec<float> pt3 = {40, 20, 30};
54 RVec<float> eta3 = {2.5, 0.5, -1.0};
55 RVec<float> phi3 = {-0.5, 0.0, 1.0};
56 RVec<float> mass3 = {10, 10, 10};
57
58 RVec<float> pt4 = {20, 10, 40};
59 RVec<float> eta4 = {0.5, -0.5, 1.0};
60 RVec<float> phi4 = {0.0, 1.0, -1.0};
61 RVec<float> mass4 = {2, 2, 2};
62
63 auto invMass = InvariantMasses(pt3, eta3, phi3, mass3, pt4, eta4, phi4, mass4);
64
65 std::cout << "\nInvariantMass(pt1 = " << pt3 << ",\n"
66 << " eta1 = " << eta3 << ",\n"
67 << " phi1 = " << phi3 << ",\n"
68 << " mass1 = " << mass3 << ",\n"
69 << " pt2 = " << pt4 << ",\n"
70 << " eta2 = " << eta4 << ",\n"
71 << " phi2 = " << phi4 << ",\n"
72 << " mass2 = " << mass4 << ")\n"
73 << " = " << invMass << "\n";
74
75 // The InvariantMass helper also accepts a single set of (pt, eta, phi, mass) vectors. Then,
76 // the invariant mass of all particles in the collection is computed.
77
78 auto invMass2 = InvariantMass(pt3, eta3, phi3, mass3);
79
80 std::cout << "\nInvariantMass(pt = " << pt3 << ",\n"
81 << " eta = " << eta3 << ",\n"
82 << " phi = " << phi3 << ",\n"
83 << " mass = " << mass3 << ")\n"
84 << " = " << invMass2 << "\n";
85}
#define M_PI
Definition: Rotated.cxx:105
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:95
Vector1::Scalar DeltaPhi(const Vector1 &v1, const Vector2 &v2)
Find aximutal Angle difference between two generic vectors ( v2.Phi() - v1.Phi() ) The only requireme...
Definition: VectorUtil.h:59
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:215
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:1554
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:1023
RVec< RVec< std::size_t > > Combinations(const std::size_t size1, const std::size_t size2)
Return the indices that represent all combinations of the elements of two RVecs.
Definition: RVec.hxx:1146