Logo ROOT  
Reference Guide
VectorUtil.h
Go to the documentation of this file.
1// @(#)root/mathcore:$Id: 9ef2a4a7bd1b62c1293920c2af2f64791c75bdd8 $
2// Authors: W. Brown, M. Fischler, L. Moneta 2005
3
4
5/**********************************************************************
6 * *
7 * Copyright (c) 2005 , LCG ROOT MathLib Team *
8 * *
9 * *
10 **********************************************************************/
11
12// Header file for Vector Utility functions
13//
14// Created by: moneta at Tue May 31 21:10:29 2005
15//
16// Last update: Tue May 31 21:10:29 2005
17//
18#ifndef ROOT_Math_GenVector_VectorUtil
19#define ROOT_Math_GenVector_VectorUtil 1
20
21#include "Math/Math.h"
22
23
25
26namespace ROOT {
27
28 namespace Math {
29
30
31 // utility functions for vector classes
32
33
34
35 /**
36 Global Helper functions for generic Vector classes. Any Vector classes implementing some defined member functions,
37 like Phi() or Eta() or mag() can use these functions.
38 The functions returning a scalar value, returns always double precision number even if the vector are
39 based on another precision type
40
41 @ingroup GenVector
42 */
43
44
45 namespace VectorUtil {
46
47
48 // methods for 3D vectors
49
50 /**
51 Find aximutal Angle difference between two generic vectors ( v2.Phi() - v1.Phi() )
52 The only requirements on the Vector classes is that they implement the Phi() method
53 \param v1 Vector of any type implementing the Phi() operator
54 \param v2 Vector of any type implementing the Phi() operator
55 \return Phi difference
56 \f[ \Delta \phi = \phi_2 - \phi_1 \f]
57 */
58 template <class Vector1, class Vector2>
59 inline typename Vector1::Scalar DeltaPhi( const Vector1 & v1, const Vector2 & v2) {
60 typename Vector1::Scalar dphi = v2.Phi() - v1.Phi();
61 if ( dphi > M_PI ) {
62 dphi -= 2.0*M_PI;
63 } else if ( dphi <= -M_PI ) {
64 dphi += 2.0*M_PI;
65 }
66 return dphi;
67 }
68
69
70
71 /**
72 Find square of the difference in pseudorapidity (Eta) and Phi betwen two generic vectors
73 The only requirements on the Vector classes is that they implement the Phi() and Eta() method
74 \param v1 Vector 1
75 \param v2 Vector 2
76 \return Angle between the two vectors
77 \f[ \Delta R2 = ( \Delta \phi )^2 + ( \Delta \eta )^2 \f]
78 */
79 template <class Vector1, class Vector2>
80 inline typename Vector1::Scalar DeltaR2( const Vector1 & v1, const Vector2 & v2) {
81 typename Vector1::Scalar dphi = DeltaPhi(v1,v2);
82 typename Vector1::Scalar deta = v2.Eta() - v1.Eta();
83 return dphi*dphi + deta*deta;
84 }
85
86 /**
87 Find difference in pseudorapidity (Eta) and Phi betwen two generic vectors
88 The only requirements on the Vector classes is that they implement the Phi() and Eta() method
89 \param v1 Vector 1
90 \param v2 Vector 2
91 \return Angle between the two vectors
92 \f[ \Delta R = \sqrt{ ( \Delta \phi )^2 + ( \Delta \eta )^2 } \f]
93 */
94 template <class Vector1, class Vector2>
95 inline typename Vector1::Scalar DeltaR( const Vector1 & v1, const Vector2 & v2) {
96 return std::sqrt( DeltaR2(v1,v2) );
97 }
98
99
100
101 /**
102 Find CosTheta Angle between two generic 3D vectors
103 pre-requisite: vectors implement the X(), Y() and Z()
104 \param v1 Vector v1
105 \param v2 Vector v2
106 \return cosine of Angle between the two vectors
107 \f[ \cos \theta = \frac { \vec{v1} \cdot \vec{v2} }{ | \vec{v1} | | \vec{v2} | } \f]
108 */
109 // this cannot be made all generic since Mag2() for 2, 3 or 4 D is different
110 // need to have a specialization for polar Coordinates ??
111 template <class Vector1, class Vector2>
112 double CosTheta( const Vector1 & v1, const Vector2 & v2) {
113 double arg;
114 double v1_r2 = v1.X()*v1.X() + v1.Y()*v1.Y() + v1.Z()*v1.Z();
115 double v2_r2 = v2.X()*v2.X() + v2.Y()*v2.Y() + v2.Z()*v2.Z();
116 double ptot2 = v1_r2*v2_r2;
117 if(ptot2 <= 0) {
118 arg = 0.0;
119 }else{
120 double pdot = v1.X()*v2.X() + v1.Y()*v2.Y() + v1.Z()*v2.Z();
121 arg = pdot/std::sqrt(ptot2);
122 if(arg > 1.0) arg = 1.0;
123 if(arg < -1.0) arg = -1.0;
124 }
125 return arg;
126 }
127
128
129 /**
130 Find Angle between two vectors.
131 Use the CosTheta() function
132 \param v1 Vector v1
133 \param v2 Vector v2
134 \return Angle between the two vectors
135 \f[ \theta = \cos ^{-1} \frac { \vec{v1} \cdot \vec{v2} }{ | \vec{v1} | | \vec{v2} | } \f]
136 */
137 template <class Vector1, class Vector2>
138 inline double Angle( const Vector1 & v1, const Vector2 & v2) {
139 return std::acos( CosTheta(v1, v2) );
140 }
141
142 /**
143 Find the projection of v along the given direction u.
144 \param v Vector v for which the propjection is to be found
145 \param u Vector specifying the direction
146 \return Vector projection (same type of v)
147 \f[ \vec{proj} = \frac{ \vec{v} \cdot \vec{u} }{|\vec{u}|}\vec{u} \f]
148 Precondition is that Vector1 implements Dot function and Vector2 implements X(),Y() and Z()
149 */
150 template <class Vector1, class Vector2>
151 Vector1 ProjVector( const Vector1 & v, const Vector2 & u) {
152 double magU2 = u.X()*u.X() + u.Y()*u.Y() + u.Z()*u.Z();
153 if (magU2 == 0) return Vector1(0,0,0);
154 double d = v.Dot(u)/magU2;
155 return Vector1( u.X() * d, u.Y() * d, u.Z() * d);
156 }
157
158 /**
159 Find the vector component of v perpendicular to the given direction of u
160 \param v Vector v for which the perpendicular component is to be found
161 \param u Vector specifying the direction
162 \return Vector component of v which is perpendicular to u
163 \f[ \vec{perp} = \vec{v} - \frac{ \vec{v} \cdot \vec{u} }{|\vec{u}|}\vec{u} \f]
164 Precondition is that Vector1 implements Dot function and Vector2 implements X(),Y() and Z()
165 */
166 template <class Vector1, class Vector2>
167 inline Vector1 PerpVector( const Vector1 & v, const Vector2 & u) {
168 return v - ProjVector(v,u);
169 }
170
171 /**
172 Find the magnitude square of the vector component of v perpendicular to the given direction of u
173 \param v Vector v for which the perpendicular component is to be found
174 \param u Vector specifying the direction
175 \return square value of the component of v which is perpendicular to u
176 \f[ perp = | \vec{v} - \frac{ \vec{v} \cdot \vec{u} }{|\vec{u}|}\vec{u} |^2 \f]
177 Precondition is that Vector1 implements Dot function and Vector2 implements X(),Y() and Z()
178 */
179 template <class Vector1, class Vector2>
180 inline double Perp2( const Vector1 & v, const Vector2 & u) {
181 double magU2 = u.X()*u.X() + u.Y()*u.Y() + u.Z()*u.Z();
182 double prjvu = v.Dot(u);
183 double magV2 = v.Dot(v);
184 return magU2 > 0.0 ? magV2-prjvu*prjvu/magU2 : magV2;
185 }
186
187 /**
188 Find the magnitude of the vector component of v perpendicular to the given direction of u
189 \param v Vector v for which the perpendicular component is to be found
190 \param u Vector specifying the direction
191 \return value of the component of v which is perpendicular to u
192 \f[ perp = | \vec{v} - \frac{ \vec{v} \cdot \vec{u} }{|\vec{u}|}\vec{u} | \f]
193 Precondition is that Vector1 implements Dot function and Vector2 implements X(),Y() and Z()
194 */
195 template <class Vector1, class Vector2>
196 inline double Perp( const Vector1 & v, const Vector2 & u) {
197 return std::sqrt(Perp2(v,u) );
198 }
199
200
201
202 // Lorentz Vector functions
203
204
205 /**
206 return the invariant mass of two LorentzVector
207 The only requirement on the LorentzVector is that they need to implement the
208 X() , Y(), Z() and E() methods.
209 \param v1 LorenzVector 1
210 \param v2 LorenzVector 2
211 \return invariant mass M
212 \f[ M_{12} = \sqrt{ (\vec{v1} + \vec{v2} ) \cdot (\vec{v1} + \vec{v2} ) } \f]
213 */
214 template <class Vector1, class Vector2>
215 inline typename Vector1::Scalar InvariantMass( const Vector1 & v1, const Vector2 & v2) {
216 typedef typename Vector1::Scalar Scalar;
217 Scalar ee = (v1.E() + v2.E() );
218 Scalar xx = (v1.X() + v2.X() );
219 Scalar yy = (v1.Y() + v2.Y() );
220 Scalar zz = (v1.Z() + v2.Z() );
221 Scalar mm2 = ee*ee - xx*xx - yy*yy - zz*zz;
222 return mm2 < 0.0 ? -std::sqrt(-mm2) : std::sqrt(mm2);
223 // PxPyPzE4D<double> q(xx,yy,zz,ee);
224 // return q.M();
225 //return ( v1 + v2).mag();
226 }
227
228 template <class Vector1, class Vector2>
229 inline typename Vector1::Scalar InvariantMass2( const Vector1 & v1, const Vector2 & v2) {
230 typedef typename Vector1::Scalar Scalar;
231 Scalar ee = (v1.E() + v2.E() );
232 Scalar xx = (v1.X() + v2.X() );
233 Scalar yy = (v1.Y() + v2.Y() );
234 Scalar zz = (v1.Z() + v2.Z() );
235 Scalar mm2 = ee*ee - xx*xx - yy*yy - zz*zz;
236 return mm2 ; // < 0.0 ? -std::sqrt(-mm2) : std::sqrt(mm2);
237 // PxPyPzE4D<double> q(xx,yy,zz,ee);
238 // return q.M();
239 //return ( v1 + v2).mag();
240 }
241
242 // rotation and transformations
243
244
245#ifndef __CINT__
246 /**
247 rotation along X axis for a generic vector by an Angle alpha
248 returning a new vector.
249 The only pre requisite on the Vector is that it has to implement the X() , Y() and Z()
250 and SetXYZ methods.
251 */
252 template <class Vector>
253 Vector RotateX(const Vector & v, double alpha) {
254 double sina = sin(alpha);
255 double cosa = cos(alpha);
256 double y2 = v.Y() * cosa - v.Z()*sina;
257 double z2 = v.Z() * cosa + v.Y() * sina;
258 Vector vrot;
259 vrot.SetXYZ(v.X(), y2, z2);
260 return vrot;
261 }
262
263 /**
264 rotation along Y axis for a generic vector by an Angle alpha
265 returning a new vector.
266 The only pre requisite on the Vector is that it has to implement the X() , Y() and Z()
267 and SetXYZ methods.
268 */
269 template <class Vector>
270 Vector RotateY(const Vector & v, double alpha) {
271 double sina = sin(alpha);
272 double cosa = cos(alpha);
273 double x2 = v.X() * cosa + v.Z() * sina;
274 double z2 = v.Z() * cosa - v.X() * sina;
275 Vector vrot;
276 vrot.SetXYZ(x2, v.Y(), z2);
277 return vrot;
278 }
279
280 /**
281 rotation along Z axis for a generic vector by an Angle alpha
282 returning a new vector.
283 The only pre requisite on the Vector is that it has to implement the X() , Y() and Z()
284 and SetXYZ methods.
285 */
286 template <class Vector>
287 Vector RotateZ(const Vector & v, double alpha) {
288 double sina = sin(alpha);
289 double cosa = cos(alpha);
290 double x2 = v.X() * cosa - v.Y() * sina;
291 double y2 = v.Y() * cosa + v.X() * sina;
292 Vector vrot;
293 vrot.SetXYZ(x2, y2, v.Z());
294 return vrot;
295 }
296
297
298 /**
299 rotation on a generic vector using a generic rotation class.
300 The only requirement on the vector is that implements the
301 X(), Y(), Z() and SetXYZ methods.
302 The requirement on the rotation matrix is that need to implement the
303 (i,j) operator returning the matrix element with R(0,0) = xx element
304 */
305 template<class Vector, class RotationMatrix>
306 Vector Rotate(const Vector &v, const RotationMatrix & rot) {
307 double xX = v.X();
308 double yY = v.Y();
309 double zZ = v.Z();
310 double x2 = rot(0,0)*xX + rot(0,1)*yY + rot(0,2)*zZ;
311 double y2 = rot(1,0)*xX + rot(1,1)*yY + rot(1,2)*zZ;
312 double z2 = rot(2,0)*xX + rot(2,1)*yY + rot(2,2)*zZ;
313 Vector vrot;
314 vrot.SetXYZ(x2,y2,z2);
315 return vrot;
316 }
317
318 /**
319 Boost a generic Lorentz Vector class using a generic 3D Vector class describing the boost
320 The only requirement on the vector is that implements the
321 X(), Y(), Z(), T() and SetXYZT methods.
322 The requirement on the boost vector is that needs to implement the
323 X(), Y() , Z() retorning the vector elements describing the boost
324 The beta of the boost must be <= 1 or a nul Lorentz Vector will be returned
325 */
326 template <class LVector, class BoostVector>
327 LVector boost(const LVector & v, const BoostVector & b) {
328 double bx = b.X();
329 double by = b.Y();
330 double bz = b.Z();
331 double b2 = bx*bx + by*by + bz*bz;
332 if (b2 >= 1) {
333 GenVector::Throw ( "Beta Vector supplied to set Boost represents speed >= c");
334 return LVector();
335 }
336 double gamma = 1.0 / std::sqrt(1.0 - b2);
337 double bp = bx*v.X() + by*v.Y() + bz*v.Z();
338 double gamma2 = b2 > 0 ? (gamma - 1.0)/b2 : 0.0;
339 double x2 = v.X() + gamma2*bp*bx + gamma*bx*v.T();
340 double y2 = v.Y() + gamma2*bp*by + gamma*by*v.T();
341 double z2 = v.Z() + gamma2*bp*bz + gamma*bz*v.T();
342 double t2 = gamma*(v.T() + bp);
343 LVector lv;
344 lv.SetXYZT(x2,y2,z2,t2);
345 return lv;
346 }
347
348
349 /**
350 Boost a generic Lorentz Vector class along the X direction with a factor beta
351 The only requirement on the vector is that implements the
352 X(), Y(), Z(), T() and SetXYZT methods.
353 The beta of the boost must be <= 1 or a nul Lorentz Vector will be returned
354 */
355 template <class LVector, class T>
356 LVector boostX(const LVector & v, T beta) {
357 if (beta >= 1) {
358 GenVector::Throw ("Beta Vector supplied to set Boost represents speed >= c");
359 return LVector();
360 }
361 T gamma = 1.0/ std::sqrt(1.0 - beta*beta);
362 typename LVector::Scalar x2 = gamma * v.X() + gamma * beta * v.T();
363 typename LVector::Scalar t2 = gamma * beta * v.X() + gamma * v.T();
364
365 LVector lv;
366 lv.SetXYZT(x2,v.Y(),v.Z(),t2);
367 return lv;
368 }
369
370 /**
371 Boost a generic Lorentz Vector class along the Y direction with a factor beta
372 The only requirement on the vector is that implements the
373 X(), Y(), Z(), T() methods and be constructed from x,y,z,t values
374 The beta of the boost must be <= 1 or a nul Lorentz Vector will be returned
375 */
376 template <class LVector>
377 LVector boostY(const LVector & v, double beta) {
378 if (beta >= 1) {
379 GenVector::Throw ("Beta Vector supplied to set Boost represents speed >= c");
380 return LVector();
381 }
382 double gamma = 1.0/ std::sqrt(1.0 - beta*beta);
383 double y2 = gamma * v.Y() + gamma * beta * v.T();
384 double t2 = gamma * beta * v.Y() + gamma * v.T();
385 LVector lv;
386 lv.SetXYZT(v.X(),y2,v.Z(),t2);
387 return lv;
388 }
389
390 /**
391 Boost a generic Lorentz Vector class along the Z direction with a factor beta
392 The only requirement on the vector is that implements the
393 X(), Y(), Z(), T() methods and be constructed from x,y,z,t values
394 The beta of the boost must be <= 1 or a nul Lorentz Vector will be returned
395 */
396 template <class LVector>
397 LVector boostZ(const LVector & v, double beta) {
398 if (beta >= 1) {
399 GenVector::Throw ( "Beta Vector supplied to set Boost represents speed >= c");
400 return LVector();
401 }
402 double gamma = 1.0/ std::sqrt(1.0 - beta*beta);
403 double z2 = gamma * v.Z() + gamma * beta * v.T();
404 double t2 = gamma * beta * v.Z() + gamma * v.T();
405 LVector lv;
406 lv.SetXYZT(v.X(),v.Y(),z2,t2);
407 return lv;
408 }
409
410#endif
411
412
413
414
415 // MATRIX VECTOR MULTIPLICATION
416 // cannot define an operator * otherwise conflicts with rotations
417 // operations like Rotation3D * vector use Mult
418
419 /**
420 Multiplications of a generic matrices with a DisplacementVector3D of any coordinate system.
421 Assume that the matrix implements the operator( i,j) and that it has at least 3 columns and 3 rows. There is no check on the matrix size !!
422 */
423 template<class Matrix, class CoordSystem, class U>
424 inline
427 vret.SetXYZ( m(0,0) * v.x() + m(0,1) * v.y() + m(0,2) * v.z() ,
428 m(1,0) * v.x() + m(1,1) * v.y() + m(1,2) * v.z() ,
429 m(2,0) * v.x() + m(2,1) * v.y() + m(2,2) * v.z() );
430 return vret;
431 }
432
433
434 /**
435 Multiplications of a generic matrices with a generic PositionVector
436 Assume that the matrix implements the operator( i,j) and that it has at least 3 columns and 3 rows. There is no check on the matrix size !!
437 */
438 template<class Matrix, class CoordSystem, class U>
439 inline
442 pret.SetXYZ( m(0,0) * p.x() + m(0,1) * p.y() + m(0,2) * p.z() ,
443 m(1,0) * p.x() + m(1,1) * p.y() + m(1,2) * p.z() ,
444 m(2,0) * p.x() + m(2,1) * p.y() + m(2,2) * p.z() );
445 return pret;
446 }
447
448
449 /**
450 Multiplications of a generic matrices with a LorentzVector described
451 in any coordinate system.
452 Assume that the matrix implements the operator( i,j) and that it has at least 4 columns and 4 rows. There is no check on the matrix size !!
453 */
454 // this will not be ambigous with operator*(Scalar, LorentzVector) since that one // Scalar is passed by value
455 template<class CoordSystem, class Matrix>
456 inline
459 vret.SetXYZT( m(0,0)*v.x() + m(0,1)*v.y() + m(0,2)*v.z() + m(0,3)* v.t() ,
460 m(1,0)*v.x() + m(1,1)*v.y() + m(1,2)*v.z() + m(1,3)* v.t() ,
461 m(2,0)*v.x() + m(2,1)*v.y() + m(2,2)*v.z() + m(2,3)* v.t() ,
462 m(3,0)*v.x() + m(3,1)*v.y() + m(3,2)*v.z() + m(3,3)* v.t() );
463 return vret;
464 }
465
466
467
468 // non-template utility functions for all objects
469
470
471 /**
472 Return a phi angle in the interval (0,2*PI]
473 */
474 double Phi_0_2pi(double phi);
475 /**
476 Returns phi angle in the interval (-PI,PI]
477 */
478 double Phi_mpi_pi(double phi);
479
480
481
482 } // end namespace Vector Util
483
484
485
486 } // end namespace Math
487
488} // end namespace ROOT
489
490
491#endif /* ROOT_Math_GenVector_VectorUtil */
#define d(i)
Definition: RSha256.hxx:102
#define b(i)
Definition: RSha256.hxx:100
static const double x2[5]
#define M_PI
Definition: Rotated.cxx:105
double acos(double)
double cos(double)
double sqrt(double)
double sin(double)
Class describing a generic displacement vector in 3 dimensions.
DisplacementVector3D< CoordSystem, Tag > & SetXYZ(Scalar a, Scalar b, Scalar c)
set the values of the vector from the cartesian components (x,y,z) (if the vector is held in polar or...
LorentzVector< CoordSystem > & SetXYZT(Scalar xx, Scalar yy, Scalar zz, Scalar tt)
set the values of the vector from the cartesian components (x,y,z,t) (if the vector is held in anothe...
Class describing a generic position vector (point) in 3 dimensions.
double beta(double x, double y)
Calculates the beta function.
Namespace for new Math classes and functions.
double gamma(double x)
double T(double x)
Definition: ChebyshevPol.h:34
void Throw(const char *)
function throwing exception, by creating internally a GenVector_exception only when needed
DisplacementVector3D< CoordSystem, U > Mult(const Matrix &m, const DisplacementVector3D< CoordSystem, U > &v)
Multiplications of a generic matrices with a DisplacementVector3D of any coordinate system.
Definition: VectorUtil.h:425
Vector1::Scalar DeltaR2(const Vector1 &v1, const Vector2 &v2)
Find square of the difference in pseudorapidity (Eta) and Phi betwen two generic vectors The only req...
Definition: VectorUtil.h:80
double Perp(const Vector1 &v, const Vector2 &u)
Find the magnitude of the vector component of v perpendicular to the given direction of u.
Definition: VectorUtil.h:196
LVector boostY(const LVector &v, double beta)
Boost a generic Lorentz Vector class along the Y direction with a factor beta The only requirement on...
Definition: VectorUtil.h:377
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
LVector boost(const LVector &v, const BoostVector &b)
Boost a generic Lorentz Vector class using a generic 3D Vector class describing the boost The only re...
Definition: VectorUtil.h:327
double Angle(const Vector1 &v1, const Vector2 &v2)
Find Angle between two vectors.
Definition: VectorUtil.h:138
Vector RotateZ(const Vector &v, double alpha)
rotation along Z axis for a generic vector by an Angle alpha returning a new vector.
Definition: VectorUtil.h:287
Vector RotateY(const Vector &v, double alpha)
rotation along Y axis for a generic vector by an Angle alpha returning a new vector.
Definition: VectorUtil.h:270
Vector RotateX(const Vector &v, double alpha)
rotation along X axis for a generic vector by an Angle alpha returning a new vector.
Definition: VectorUtil.h:253
double CosTheta(const Vector1 &v1, const Vector2 &v2)
Find CosTheta Angle between two generic 3D vectors pre-requisite: vectors implement the X(),...
Definition: VectorUtil.h:112
Vector1 PerpVector(const Vector1 &v, const Vector2 &u)
Find the vector component of v perpendicular to the given direction of u.
Definition: VectorUtil.h:167
double Perp2(const Vector1 &v, const Vector2 &u)
Find the magnitude square of the vector component of v perpendicular to the given direction of u.
Definition: VectorUtil.h:180
double Phi_0_2pi(double phi)
Return a phi angle in the interval (0,2*PI].
Definition: VectorUtil.cxx:22
Vector Rotate(const Vector &v, const RotationMatrix &rot)
rotation on a generic vector using a generic rotation class.
Definition: VectorUtil.h:306
double Phi_mpi_pi(double phi)
Returns phi angle in the interval (-PI,PI].
Definition: VectorUtil.cxx:36
Vector1::Scalar InvariantMass2(const Vector1 &v1, const Vector2 &v2)
Definition: VectorUtil.h:229
Vector1 ProjVector(const Vector1 &v, const Vector2 &u)
Find the projection of v along the given direction u.
Definition: VectorUtil.h:151
LVector boostZ(const LVector &v, double beta)
Boost a generic Lorentz Vector class along the Z direction with a factor beta The only requirement on...
Definition: VectorUtil.h:397
LVector boostX(const LVector &v, T beta)
Boost a generic Lorentz Vector class along the X direction with a factor beta The only requirement on...
Definition: VectorUtil.h:356
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
Rotation3D::Scalar Scalar
VSD Structures.
Definition: StringConv.hxx:21
static constexpr double mm2
auto * lv
Definition: textalign.C:5
auto * m
Definition: textangle.C:8