Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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
14void vo007_PhysicsHelpers()
15{
16 // The DeltaPhi helper computes the closest angle between angles.
17 // This means that the resulting value is in the range [-pi, pi].
18
19 // Note that the helper also supports to compute the angle difference of an
20 // RVec and a scalar and two scalars. In addition, the computation of the
21 // difference and the behaviour at the boundary can be adjusted to radian and
22 // degrees.
23 ROOT::RVecF phis = {0.0, 1.0, -0.5, M_PI + 1.0};
24 auto idx = Combinations(phis, 2);
25
26 auto phi1 = Take(phis, idx[0]);
27 auto phi2 = Take(phis, idx[1]);
28 auto dphi = DeltaPhi(phi1, phi2);
29
30 std::cout << "DeltaPhi(phi1 = " << phi1 << ",\n"
31 << " phi2 = " << phi2 << ")\n"
32 << " = " << dphi << "\n";
33
34 // The DeltaR helper is similar to the DeltaPhi helper and computes the distance
35 // in the \f$\eta\f$-\f$\phi\f$ plane.
36 ROOT::RVecF etas = {2.4, -1.5, 1.0, 0.0};
37
38 auto eta1 = Take(etas, idx[0]);
39 auto eta2 = Take(etas, idx[1]);
40 auto dr = DeltaR(eta1, eta2, phi1, phi2);
41
42 std::cout << "\nDeltaR(eta1 = " << eta1 << ",\n"
43 << " eta2 = " << eta2 << ",\n"
44 << " phi1 = " << phi1 << ",\n"
45 << " phi2 = " << phi2 << ")\n"
46 << " = " << dr << "\n";
47
48 // The InvariantMasses helper computes the invariant mass of a two particle system
49 // given the properties transverse momentum (pt), rapidity (eta), azimuth (phi)
50 // and mass.
51 ROOT::RVecF pt3 = {40, 20, 30};
52 ROOT::RVecF eta3 = {2.5, 0.5, -1.0};
53 ROOT::RVecF phi3 = {-0.5, 0.0, 1.0};
54 ROOT::RVecF mass3 = {10, 10, 10};
55
56 ROOT::RVecF pt4 = {20, 10, 40};
57 ROOT::RVecF eta4 = {0.5, -0.5, 1.0};
58 ROOT::RVecF phi4 = {0.0, 1.0, -1.0};
59 ROOT::RVecF mass4 = {2, 2, 2};
60
61 auto invMass = InvariantMasses(pt3, eta3, phi3, mass3, pt4, eta4, phi4, mass4);
62
63 std::cout << "\nInvariantMass(pt1 = " << pt3 << ",\n"
64 << " eta1 = " << eta3 << ",\n"
65 << " phi1 = " << phi3 << ",\n"
66 << " mass1 = " << mass3 << ",\n"
67 << " pt2 = " << pt4 << ",\n"
68 << " eta2 = " << eta4 << ",\n"
69 << " phi2 = " << phi4 << ",\n"
70 << " mass2 = " << mass4 << ")\n"
71 << " = " << invMass << "\n";
72
73 // The InvariantMass helper also accepts a single set of (pt, eta, phi, mass) vectors. Then,
74 // the invariant mass of all particles in the collection is computed.
75
76 auto invMass2 = InvariantMass(pt3, eta3, phi3, mass3);
77
78 std::cout << "\nInvariantMass(pt = " << pt3 << ",\n"
79 << " eta = " << eta3 << ",\n"
80 << " phi = " << phi3 << ",\n"
81 << " mass = " << mass3 << ")\n"
82 << " = " << invMass2 << "\n";
83}
#define M_PI
Definition Rotated.cxx:105
RVec< Common_t > InvariantMasses(const RVec< T0 > &pt1, const RVec< T1 > &eta1, const RVec< T2 > &phi1, const RVec< T3 > &mass1, const RVec< T4 > &pt2, const RVec< T5 > &eta2, const RVec< T6 > &phi2, const RVec< T7 > &mass2)
Return the invariant mass of two particles given the collections of the quantities transverse momentu...
Definition RVec.hxx:2992
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:2302
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:2569
Vector1::Scalar DeltaR(const Vector1 &v1, const Vector2 &v2)
Find difference in pseudorapidity (Eta) and Phi between two generic vectors The only requirements on ...
Definition VectorUtil.h:112
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:61
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:248