Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TLorentzVector.h
Go to the documentation of this file.
1// @(#)root/physics:$Id$
2// Author: Pasha Murat , Peter Malzacher 12/02/99
3
4/*************************************************************************
5 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12#ifndef ROOT_TLorentzVector
13#define ROOT_TLorentzVector
14
15
16//////////////////////////////////////////////////////////////////////////
17// //
18// TLorentzVector //
19// //
20// Place holder for real lorentz vector class. //
21// //
22//////////////////////////////////////////////////////////////////////////
23
24#include "TMath.h"
25#include "TVector3.h"
26#include "TRotation.h"
27#include "Math/Vector4D.h"
28
30
31class TLorentzVector : public TObject {
32
33private:
34
35 TVector3 fP; // 3 vector component
36 Double_t fE; // time or energy of (x,y,z,t) or (px,py,pz,e)
37
38public:
39
40 typedef Double_t Scalar; // to be able to use it with the ROOT::Math::VectorUtil functions
41
42 enum { kX=0, kY=1, kZ=2, kT=3, kNUM_COORDINATES=4, kSIZE=kNUM_COORDINATES };
43 // Safe indexing of the coordinates when using with matrices, arrays, etc.
44
46
48 // Constructor giving the components x, y, z, t.
49
50 TLorentzVector(const Double_t * carray);
51 TLorentzVector(const Float_t * carray);
52 // Constructor from an array, not checked!
53
54 TLorentzVector(const TVector3 & vector3, Double_t t);
55 // Constructor giving a 3-Vector and a time component.
56
57 TLorentzVector(const TLorentzVector & lorentzvector);
58 // Copy constructor.
59
60 ~TLorentzVector() override{};
61 // Destructor
62
63 // inline operator TVector3 () const;
64 // inline operator TVector3 & ();
65 // Conversion (cast) to TVector3.
66
67 inline Double_t X() const;
68 inline Double_t Y() const;
69 inline Double_t Z() const;
70 inline Double_t T() const;
71 // Get position and time.
72
73 inline void SetX(Double_t a);
74 inline void SetY(Double_t a);
75 inline void SetZ(Double_t a);
76 inline void SetT(Double_t a);
77 // Set position and time.
78
79 inline Double_t Px() const;
80 inline Double_t Py() const;
81 inline Double_t Pz() const;
82 inline Double_t P() const;
83 inline Double_t E() const;
84 inline Double_t Energy() const;
85 // Get momentum and energy.
86
87 inline void SetPx(Double_t a);
88 inline void SetPy(Double_t a);
89 inline void SetPz(Double_t a);
90 inline void SetE(Double_t a);
91 // Set momentum and energy.
92
93 inline TVector3 Vect() const ;
94 // Get spatial component.
95
96 inline void SetVect(const TVector3 & vect3);
97 // Set spatial component.
98
99 inline Double_t Theta() const;
100 inline Double_t CosTheta() const;
101 inline Double_t Phi() const; //returns phi from -pi to pi
102 inline Double_t Rho() const;
103 // Get spatial vector components in spherical coordinate system.
104
105 inline void SetTheta(Double_t theta);
106 inline void SetPhi(Double_t phi);
107 inline void SetRho(Double_t rho);
108 // Set spatial vector components in spherical coordinate system.
109
110 inline void SetPxPyPzE(Double_t px, Double_t py, Double_t pz, Double_t e);
111 inline void SetXYZT(Double_t x, Double_t y, Double_t z, Double_t t);
112 inline void SetXYZM(Double_t x, Double_t y, Double_t z, Double_t m);
113 inline void SetPtEtaPhiM(Double_t pt, Double_t eta, Double_t phi, Double_t m);
114 inline void SetPtEtaPhiE(Double_t pt, Double_t eta, Double_t phi, Double_t e);
115 // Setters to provide the functionality (but a more meanigful name) of
116 // the previous version eg SetV4... PsetV4...
117
118 inline void GetXYZT(Double_t *carray) const;
119 inline void GetXYZT(Float_t *carray) const;
120 // Getters into an arry
121 // no checking!
122
123 Double_t operator () (int i) const;
124 inline Double_t operator [] (int i) const;
125 // Get components by index.
126
127 Double_t & operator () (int i);
128 inline Double_t & operator [] (int i);
129 // Set components by index.
130
131 inline TLorentzVector & operator = (const TLorentzVector &);
132 // Assignment.
133
134 inline TLorentzVector operator + (const TLorentzVector &) const;
136 // Additions.
137
138 inline TLorentzVector operator - (const TLorentzVector &) const;
140 // Subtractions.
141
142 inline TLorentzVector operator - () const;
143 // Unary minus.
144
145 inline TLorentzVector operator * (Double_t a) const;
147 // Scaling with real numbers.
148
149 inline Bool_t operator == (const TLorentzVector &) const;
150 inline Bool_t operator != (const TLorentzVector &) const;
151 // Comparisons.
152
153 inline Double_t Perp2() const;
154 // Transverse component of the spatial vector squared.
155
156 inline Double_t Pt() const;
157 inline Double_t Perp() const;
158 // Transverse component of the spatial vector (R in cylindrical system).
159
160 inline void SetPerp(Double_t);
161 // Set the transverse component of the spatial vector.
162
163 inline Double_t Perp2(const TVector3 & v) const;
164 // Transverse component of the spatial vector w.r.t. given axis squared.
165
166 inline Double_t Pt(const TVector3 & v) const;
167 inline Double_t Perp(const TVector3 & v) const;
168 // Transverse component of the spatial vector w.r.t. given axis.
169
170 inline Double_t Et2() const;
171 // Transverse energy squared.
172
173 inline Double_t Et() const;
174 // Transverse energy.
175
176 inline Double_t Et2(const TVector3 &) const;
177 // Transverse energy w.r.t. given axis squared.
178
179 inline Double_t Et(const TVector3 &) const;
180 // Transverse energy w.r.t. given axis.
181
182 inline Double_t DeltaPhi(const TLorentzVector &) const;
183 inline Double_t DeltaR(const TLorentzVector &, Bool_t useRapidity=kFALSE) const;
184 inline Double_t DrEtaPhi(const TLorentzVector &) const;
185 inline Double_t DrRapidityPhi(const TLorentzVector &) const;
186 inline TVector2 EtaPhiVector();
187
188 inline Double_t Angle(const TVector3 & v) const;
189 // Angle wrt. another vector.
190
191 inline Double_t Mag2() const;
192 inline Double_t M2() const;
193 // Invariant mass squared.
194
195 inline Double_t Mag() const;
196 inline Double_t M() const;
197 // Invariant mass. If mag2() is negative then -sqrt(-mag2()) is returned.
198
199 inline Double_t Mt2() const;
200 // Transverse mass squared.
201
202 inline Double_t Mt() const;
203 // Transverse mass.
204
205 inline Double_t Beta() const;
206 inline Double_t Gamma() const;
207
208 inline Double_t Dot(const TLorentzVector &) const;
209 inline Double_t operator * (const TLorentzVector &) const;
210 // Scalar product.
211
212 inline void SetVectMag(const TVector3 & spatial, Double_t magnitude);
213 inline void SetVectM(const TVector3 & spatial, Double_t mass);
214 // Copy spatial coordinates, and set energy = sqrt(mass^2 + spatial^2)
215
216 inline Double_t Plus() const;
217 inline Double_t Minus() const;
218 // Returns t +/- z.
219 // Related to the positive/negative light-cone component,
220 // which some define this way and others define as (t +/- z)/sqrt(2)
221
222 inline TVector3 BoostVector() const ;
223 // Returns the spatial components divided by the time component.
224
226 inline void Boost(const TVector3 &);
227 // Lorentz boost.
228
229 Double_t Rapidity() const;
230 // Returns the rapidity, i.e. 0.5*ln((E+pz)/(E-pz))
231
232 inline Double_t Eta() const;
233 inline Double_t PseudoRapidity() const;
234 // Returns the pseudo-rapidity, i.e. -ln(tan(theta/2))
235
236 inline void RotateX(Double_t angle);
237 // Rotate the spatial component around the x-axis.
238
239 inline void RotateY(Double_t angle);
240 // Rotate the spatial component around the y-axis.
241
242 inline void RotateZ(Double_t angle);
243 // Rotate the spatial component around the z-axis.
244
245 inline void RotateUz(const TVector3 & newUzVector);
246 // Rotates the reference frame from Uz to newUz (unit vector).
247
248 inline void Rotate(Double_t, const TVector3 &);
249 // Rotate the spatial component around specified axis.
250
251 inline TLorentzVector & operator *= (const TRotation &);
252 inline TLorentzVector & Transform(const TRotation &);
253 // Transformation with HepRotation.
254
257 // Transformation with HepLorenzRotation.
258
259 operator ROOT::Math::PxPyPzEVector() const {
260 return {Px(), Py(), Pz(), E()};
261 }
262
263 void Print(Option_t *option="") const override;
264
265 ClassDefOverride(TLorentzVector,4) // A four vector with (-,-,-,+) metric
266};
267
268
269//inline TLorentzVector operator * (const TLorentzVector &, Double_t a);
270// moved to TLorentzVector::operator * (Double_t a)
272// Scaling LorentzVector with a real number
273
274
275inline Double_t TLorentzVector::X() const { return fP.X(); }
276inline Double_t TLorentzVector::Y() const { return fP.Y(); }
277inline Double_t TLorentzVector::Z() const { return fP.Z(); }
278inline Double_t TLorentzVector::T() const { return fE; }
279
283inline void TLorentzVector::SetT(Double_t a) { fE = a; }
284
285inline Double_t TLorentzVector::Px() const { return X(); }
286inline Double_t TLorentzVector::Py() const { return Y(); }
287inline Double_t TLorentzVector::Pz() const { return Z(); }
288inline Double_t TLorentzVector::P() const { return fP.Mag(); }
289inline Double_t TLorentzVector::E() const { return T(); }
290inline Double_t TLorentzVector::Energy() const { return T(); }
291
296
297inline TVector3 TLorentzVector::Vect() const { return fP; }
298
299inline void TLorentzVector::SetVect(const TVector3 &p) { fP = p; }
300
302 return fP.Phi();
303}
304
306 return fP.Theta();
307}
308
310 return fP.CosTheta();
311}
312
313
315 return fP.Mag();
316}
317
319 fP.SetTheta(th);
320}
321
323 fP.SetPhi(phi);
324}
325
327 fP.SetMag(rho);
328}
329
331 fP.SetXYZ(x, y, z);
332 SetT(t);
333}
334
336 SetXYZT(px, py, pz, e);
337}
338
340 if ( m >= 0 )
341 SetXYZT( x, y, z, TMath::Sqrt(x*x+y*y+z*z+m*m) );
342 else
343 SetXYZT( x, y, z, TMath::Sqrt( TMath::Max((x*x+y*y+z*z-m*m), 0. ) ) );
344}
345
347 pt = TMath::Abs(pt);
348 SetXYZM(pt*TMath::Cos(phi), pt*TMath::Sin(phi), pt*sinh(eta) ,m);
349}
350
352 pt = TMath::Abs(pt);
353 SetXYZT(pt*TMath::Cos(phi), pt*TMath::Sin(phi), pt*sinh(eta) ,e);
354}
355
356inline void TLorentzVector::GetXYZT(Double_t *carray) const {
357 fP.GetXYZ(carray);
358 carray[3] = fE;
359}
360
361inline void TLorentzVector::GetXYZT(Float_t *carray) const{
362 fP.GetXYZ(carray);
363 carray[3] = fE;
364}
365
366inline Double_t & TLorentzVector::operator [] (int i) { return (*this)(i); }
367inline Double_t TLorentzVector::operator [] (int i) const { return (*this)(i); }
368
370 fP = q.Vect();
371 fE = q.T();
372 return *this;
373}
374
376 return TLorentzVector(fP+q.Vect(), fE+q.T());
377}
378
380 fP += q.Vect();
381 fE += q.T();
382 return *this;
383}
384
386 return TLorentzVector(fP-q.Vect(), fE-q.T());
387}
388
390 fP -= q.Vect();
391 fE -= q.T();
392 return *this;
393}
394
396 return TLorentzVector(-X(), -Y(), -Z(), -T());
397}
398
400 fP *= a;
401 fE *= a;
402 return *this;
403}
404
406 return TLorentzVector(a*X(), a*Y(), a*Z(), a*T());
407}
408
410 return (Vect() == q.Vect() && T() == q.T());
411}
412
414 return (Vect() != q.Vect() || T() != q.T());
415}
416
417inline Double_t TLorentzVector::Perp2() const { return fP.Perp2(); }
418
419inline Double_t TLorentzVector::Perp() const { return fP.Perp(); }
420
421inline Double_t TLorentzVector::Pt() const { return Perp(); }
422
424 fP.SetPerp(r);
425}
426
428 return fP.Perp2(v);
429}
430
432 return fP.Perp(v);
433}
434
435inline Double_t TLorentzVector::Pt(const TVector3 &v) const {
436 return Perp(v);
437}
438
440 Double_t pt2 = fP.Perp2();
441 return pt2 == 0 ? 0 : E()*E() * pt2/(pt2+Z()*Z());
442}
443
445 Double_t etet = Et2();
446 return E() < 0.0 ? -sqrt(etet) : sqrt(etet);
447}
448
449inline Double_t TLorentzVector::Et2(const TVector3 & v) const {
450 Double_t pt2 = fP.Perp2(v);
451 Double_t pv = fP.Dot(v.Unit());
452 return pt2 == 0 ? 0 : E()*E() * pt2/(pt2+pv*pv);
453}
454
455inline Double_t TLorentzVector::Et(const TVector3 & v) const {
456 Double_t etet = Et2(v);
457 return E() < 0.0 ? -sqrt(etet) : sqrt(etet);
458}
459
461 return TVector2::Phi_mpi_pi(Phi()-v.Phi());
462}
463
465 return PseudoRapidity();
466}
467
468inline Double_t TLorentzVector::DeltaR(const TLorentzVector & v, const Bool_t useRapidity) const {
469 if(useRapidity){
470 Double_t drap = Rapidity()-v.Rapidity();
471 Double_t dphi = TVector2::Phi_mpi_pi(Phi()-v.Phi());
472 return TMath::Sqrt( drap*drap+dphi*dphi );
473 } else {
474 Double_t deta = Eta()-v.Eta();
475 Double_t dphi = TVector2::Phi_mpi_pi(Phi()-v.Phi());
476 return TMath::Sqrt( deta*deta+dphi*dphi );
477 }
478}
479
481 return DeltaR(v);
482}
483
485 return DeltaR(v, kTRUE);
486}
487
489 return TVector2 (Eta(),Phi());
490}
491
492
494 return fP.Angle(v);
495}
496
498 return T()*T() - fP.Mag2();
499}
500
502 Double_t mm = Mag2();
503 return mm < 0.0 ? -TMath::Sqrt(-mm) : TMath::Sqrt(mm);
504}
505
506inline Double_t TLorentzVector::M2() const { return Mag2(); }
507inline Double_t TLorentzVector::M() const { return Mag(); }
508
510 return E()*E() - Z()*Z();
511}
512
514 Double_t mm = Mt2();
515 return mm < 0.0 ? -TMath::Sqrt(-mm) : TMath::Sqrt(mm);
516}
517
519 return fP.Mag() / fE;
520}
521
523 Double_t b = Beta();
524 return 1.0/TMath::Sqrt(1- b*b);
525}
526
527inline void TLorentzVector::SetVectMag(const TVector3 & spatial, Double_t magnitude) {
528 SetXYZM(spatial.X(), spatial.Y(), spatial.Z(), magnitude);
529}
530
531inline void TLorentzVector::SetVectM(const TVector3 & spatial, Double_t mass) {
532 SetVectMag(spatial, mass);
533}
534
536 return T()*q.T() - Z()*q.Z() - Y()*q.Y() - X()*q.X();
537}
538
540 return Dot(q);
541}
542
543//Member functions Plus() and Minus() return the positive and negative
544//light-cone components:
545//
546// Double_t pcone = v.Plus();
547// Double_t mcone = v.Minus();
548//
549//CAVEAT: The values returned are T{+,-}Z. It is known that some authors
550//find it easier to define these components as (T{+,-}Z)/sqrt(2). Thus
551//check what definition is used in the physics you're working in and adapt
552//your code accordingly.
553
555 return T() + Z();
556}
557
559 return T() - Z();
560}
561
563 return TVector3(X()/T(), Y()/T(), Z()/T());
564}
565
566inline void TLorentzVector::Boost(const TVector3 & b) {
567 Boost(b.X(), b.Y(), b.Z());
568}
569
571 return fP.PseudoRapidity();
572}
573
576}
577
580}
581
584}
585
586inline void TLorentzVector::RotateUz(const TVector3 &newUzVector) {
587 fP.RotateUz(newUzVector);
588}
589
591 fP.Rotate(a,v);
592}
593
595 fP *= m;
596 return *this;
597}
598
600 fP.Transform(m);
601 return *this;
602}
603
605 return TLorentzVector(a*p.X(), a*p.Y(), a*p.Z(), a*p.T());
606}
607
609 : fP(), fE(0.0) {}
610
612 : fP(x,y,z), fE(t) {}
613
615 : fP(x0), fE(x0[3]) {}
616
618 : fP(x0), fE(x0[3]) {}
619
621 : fP(p), fE(e) {}
622
624 , fP(p.Vect()), fE(p.T()) {}
625
626
627
629{
630 //dereferencing operator const
631 switch(i) {
632 case kX:
633 return fP.X();
634 case kY:
635 return fP.Y();
636 case kZ:
637 return fP.Z();
638 case kT:
639 return fE;
640 default:
641 Error("operator()()", "bad index (%d) returning 0",i);
642 }
643 return 0.;
644}
645
647{
648 //dereferencing operator
649 switch(i) {
650 case kX:
651 return fP.fX;
652 case kY:
653 return fP.fY;
654 case kZ:
655 return fP.fZ;
656 case kT:
657 return fE;
658 default:
659 Error("operator()()", "bad index (%d) returning &fE",i);
660 }
661 return fE;
662}
663
664#endif
#define b(i)
Definition RSha256.hxx:100
#define a(i)
Definition RSha256.hxx:99
#define e(i)
Definition RSha256.hxx:103
bool Bool_t
Definition RtypesCore.h:63
float Float_t
Definition RtypesCore.h:57
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
double Double_t
Definition RtypesCore.h:59
constexpr Bool_t kTRUE
Definition RtypesCore.h:93
const char Option_t
Definition RtypesCore.h:66
#define ClassDefOverride(name, id)
Definition Rtypes.h:346
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint angle
float * q
TLorentzVector operator*(Double_t a, const TLorentzVector &)
<div class="legacybox"><h2>Legacy Code</h2> TLorentzRotation is a legacy interface: there will be no ...
Double_t Energy() const
Double_t Minus() const
void RotateZ(Double_t angle)
void Rotate(Double_t, const TVector3 &)
Bool_t operator!=(const TLorentzVector &) const
Double_t Rapidity() const
TLorentzVector & operator+=(const TLorentzVector &)
TLorentzVector operator+(const TLorentzVector &) const
TLorentzVector & operator=(const TLorentzVector &)
void RotateX(Double_t angle)
void SetPy(Double_t a)
TLorentzVector & operator*=(Double_t a)
Bool_t operator==(const TLorentzVector &) const
void Print(Option_t *option="") const override
Print the TLorentz vector components as (x,y,z,t) and (P,eta,phi,E) representations.
Double_t Angle(const TVector3 &v) const
void GetXYZT(Double_t *carray) const
Double_t Rho() const
Double_t CosTheta() const
Double_t operator[](int i) const
void SetXYZT(Double_t x, Double_t y, Double_t z, Double_t t)
Double_t Plus() const
void SetY(Double_t a)
void RotateUz(const TVector3 &newUzVector)
void SetPerp(Double_t)
Double_t Beta() const
Double_t Y() const
TVector3 Vect() const
void SetT(Double_t a)
Double_t Px() const
Double_t Mag() const
Double_t Dot(const TLorentzVector &) const
Double_t M() const
Double_t M2() const
Double_t Pt() const
Double_t Perp2() const
TLorentzVector operator*(Double_t a) const
void SetPtEtaPhiE(Double_t pt, Double_t eta, Double_t phi, Double_t e)
void SetE(Double_t a)
Double_t Pz() const
Double_t X() const
Double_t DrRapidityPhi(const TLorentzVector &) const
Double_t Mt2() const
Double_t Py() const
void SetPtEtaPhiM(Double_t pt, Double_t eta, Double_t phi, Double_t m)
void SetPhi(Double_t phi)
void RotateY(Double_t angle)
void Boost(Double_t, Double_t, Double_t)
Double_t Eta() const
void SetPz(Double_t a)
TLorentzVector operator-() const
Double_t P() const
TVector2 EtaPhiVector()
void SetRho(Double_t rho)
TLorentzVector & Transform(const TRotation &)
void SetVectMag(const TVector3 &spatial, Double_t magnitude)
Double_t Et2() const
Double_t PseudoRapidity() const
Double_t Theta() const
void SetPx(Double_t a)
~TLorentzVector() override
void SetPxPyPzE(Double_t px, Double_t py, Double_t pz, Double_t e)
Double_t Mag2() const
Double_t Et() const
Double_t Gamma() const
void SetTheta(Double_t theta)
Double_t operator()(int i) const
Double_t Phi() const
Double_t DeltaPhi(const TLorentzVector &) const
Double_t Mt() const
Double_t Perp() const
void SetZ(Double_t a)
Double_t E() const
void SetVect(const TVector3 &vect3)
void SetVectM(const TVector3 &spatial, Double_t mass)
void SetXYZM(Double_t x, Double_t y, Double_t z, Double_t m)
Double_t DeltaR(const TLorentzVector &, Bool_t useRapidity=kFALSE) const
TVector3 BoostVector() const
Double_t DrEtaPhi(const TLorentzVector &) const
TLorentzVector & operator-=(const TLorentzVector &)
Double_t Z() const
Double_t T() const
void SetX(Double_t a)
Mother of all ROOT objects.
Definition TObject.h:41
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:993
<div class="legacybox"><h2>Legacy Code</h2> TRotation is a legacy interface: there will be no bug fix...
Definition TRotation.h:20
static Double_t Phi_mpi_pi(Double_t x)
Returns phi angle in the interval [-PI,PI)
Definition TVector2.cxx:109
Double_t Z() const
Definition TVector3.h:218
void SetY(Double_t)
Definition TVector3.h:224
void SetXYZ(Double_t x, Double_t y, Double_t z)
Definition TVector3.h:227
void SetPhi(Double_t)
Set phi keeping mag and theta constant (BaBar).
Definition TVector3.cxx:368
Double_t fZ
Definition TVector3.h:185
Double_t fX
Definition TVector3.h:185
Double_t Phi() const
Return the azimuth angle. Returns phi from -pi to pi.
Definition TVector3.cxx:236
Double_t Y() const
Definition TVector3.h:217
Double_t Dot(const TVector3 &) const
Definition TVector3.h:331
void RotateZ(Double_t)
Rotate vector around Z.
Definition TVector3.cxx:285
void SetPerp(Double_t)
Definition TVector3.h:388
Double_t Mag2() const
Definition TVector3.h:339
Double_t X() const
Definition TVector3.h:216
void RotateX(Double_t)
Rotate vector around X.
Definition TVector3.cxx:263
Double_t Angle(const TVector3 &) const
Return the angle w.r.t. another 3-vector.
Definition TVector3.cxx:203
Double_t PseudoRapidity() const
Double_t m = Mag(); return 0.5*log( (m+fZ)/(m-fZ) ); guard against Pt=0.
Definition TVector3.cxx:326
void Rotate(Double_t, const TVector3 &)
Rotate vector.
Definition TVector3.cxx:296
void GetXYZ(Double_t *carray) const
Definition TVector3.h:233
Double_t CosTheta() const
Definition TVector3.h:371
Double_t Mag() const
Definition TVector3.h:86
Double_t Theta() const
Return the polar angle.
Definition TVector3.cxx:244
Double_t fY
Definition TVector3.h:185
void SetMag(Double_t)
Definition TVector3.h:376
Double_t Perp() const
Return the transverse component (R in cylindrical coordinate system)
Definition TVector3.cxx:219
TVector3 & Transform(const TRotation &)
Transform this vector with a TRotation.
Definition TVector3.cxx:196
void SetX(Double_t)
Definition TVector3.h:223
Double_t Perp2() const
Definition TVector3.h:353
void RotateY(Double_t)
Rotate vector around Y.
Definition TVector3.cxx:274
void SetZ(Double_t)
Definition TVector3.h:225
void RotateUz(const TVector3 &)
NewUzVector must be normalized !
Definition TVector3.cxx:305
void SetTheta(Double_t)
Set theta keeping mag and phi constant (BaBar).
Definition TVector3.cxx:356
TPaveText * pt
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
LorentzVector< PxPyPzE4D< double > > PxPyPzEVector
Definition Vector4Dfwd.h:48
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:250
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:662
Double_t Cos(Double_t)
Returns the cosine of an angle of x radians.
Definition TMath.h:594
Double_t Sin(Double_t)
Returns the sine of an angle of x radians.
Definition TMath.h:588
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123
#define Dot(u, v)
Definition normal.c:49
TMarker m
Definition textangle.C:8