Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TVector3.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#ifndef ROOT_TVector3
12#define ROOT_TVector3
13
14#include "TError.h"
15#include "TVector2.h"
16#include "TMatrix.h"
17#include "TMath.h"
18
19class TRotation;
20
21
22class TVector3 : public TObject {
23
24public:
25
26 typedef Double_t Scalar; // to be able to use it with the ROOT::Math::VectorUtil functions
27
28 TVector3();
29
31 // The constructor.
32
33 TVector3(const Double_t *);
34 TVector3(const Float_t *);
35 // Constructors from an array
36
37 TVector3(const TVector3 &) noexcept;
38 // The copy constructor.
39
40 ~TVector3() override = default;
41 // Destructor
42
43 /// The length is always 3. For compatibility with the standard library.
44 constexpr std::size_t size() const { return 3; }
45
46 Double_t operator () (int) const;
47 inline Double_t operator [] (int) const;
48 // Get components by index (Geant4).
49
51 inline Double_t & operator [] (int);
52 // Set components by index.
53
54 inline Double_t x() const;
55 inline Double_t y() const;
56 inline Double_t z() const;
57 inline Double_t X() const;
58 inline Double_t Y() const;
59 inline Double_t Z() const;
60 inline Double_t Px() const;
61 inline Double_t Py() const;
62 inline Double_t Pz() const;
63 // The components in cartesian coordinate system.
64
65 inline void SetX(Double_t);
66 inline void SetY(Double_t);
67 inline void SetZ(Double_t);
68 inline void SetXYZ(Double_t x, Double_t y, Double_t z);
69 void SetPtEtaPhi(Double_t pt, Double_t eta, Double_t phi);
70 void SetPtThetaPhi(Double_t pt, Double_t theta, Double_t phi);
71
72 inline void GetXYZ(Double_t *carray) const;
73 inline void GetXYZ(Float_t *carray) const;
74 // Get the components into an array
75 // not checked!
76
77 Double_t Phi() const;
78 // The azimuth angle. returns phi from -pi to pi
79
80 Double_t Theta() const;
81 // The polar angle.
82
83 inline Double_t CosTheta() const;
84 // Cosine of the polar angle.
85
86 inline Double_t Mag2() const;
87 // The magnitude squared (rho^2 in spherical coordinate system).
88
89 Double_t Mag() const { return TMath::Sqrt(Mag2()); }
90 // The magnitude (rho in spherical coordinate system).
91
92 void SetPhi(Double_t);
93 // Set phi keeping mag and theta constant (BaBar).
94
95 void SetTheta(Double_t);
96 // Set theta keeping mag and phi constant (BaBar).
97
98 inline void SetMag(Double_t);
99 // Set magnitude keeping theta and phi constant (BaBar).
100
101 inline Double_t Perp2() const;
102 // The transverse component squared (R^2 in cylindrical coordinate system).
103
104 inline Double_t Pt() const;
105 Double_t Perp() const;
106 // The transverse component (R in cylindrical coordinate system).
107
108 inline void SetPerp(Double_t);
109 // Set the transverse component keeping phi and z constant.
110
111 inline Double_t Perp2(const TVector3 &) const;
112 // The transverse component w.r.t. given axis squared.
113
114 inline Double_t Pt(const TVector3 &) const;
115 Double_t Perp(const TVector3 &) const;
116 // The transverse component w.r.t. given axis.
117
118 inline Double_t DeltaPhi(const TVector3 &) const;
119 Double_t DeltaR(const TVector3 &) const;
120 inline Double_t DrEtaPhi(const TVector3 &) const;
121 inline TVector2 EtaPhiVector() const;
122 void SetMagThetaPhi(Double_t mag, Double_t theta, Double_t phi);
123
124 inline TVector3 & operator = (const TVector3 &);
125 // Assignment.
126
127 inline Bool_t operator == (const TVector3 &) const;
128 inline Bool_t operator != (const TVector3 &) const;
129 // Comparisons (Geant4).
130
131 inline TVector3 & operator += (const TVector3 &);
132 // Addition.
133
134 inline TVector3 & operator -= (const TVector3 &);
135 // Subtraction.
136
137 inline TVector3 operator - () const;
138 // Unary minus.
139
141 // Scaling with real numbers.
142
143 TVector3 Unit() const;
144 // Unit vector parallel to this.
145
146 inline TVector3 Orthogonal() const;
147 // Vector orthogonal to this (Geant4).
148
149 inline Double_t Dot(const TVector3 &) const;
150 // Scalar product.
151
152 inline TVector3 Cross(const TVector3 &) const;
153 // Cross product.
154
155 Double_t Angle(const TVector3 &) const;
156 // The angle w.r.t. another 3-vector.
157
158 Double_t PseudoRapidity() const;
159 // Returns the pseudo-rapidity, i.e. -ln(tan(theta/2))
160
161 inline Double_t Eta() const;
162
163 void RotateX(Double_t);
164 // Rotates the Hep3Vector around the x-axis.
165
166 void RotateY(Double_t);
167 // Rotates the Hep3Vector around the y-axis.
168
169 void RotateZ(Double_t);
170 // Rotates the Hep3Vector around the z-axis.
171
172 void RotateUz(const TVector3&);
173 // Rotates reference frame from Uz to newUz (unit vector) (Geant4).
174
175 void Rotate(Double_t, const TVector3 &);
176 // Rotates around the axis specified by another Hep3Vector.
177
179 TVector3 & Transform(const TRotation &);
180 // Transformation with a Rotation matrix.
181
182 inline TVector2 XYvector() const;
183
184 void Print(Option_t* option="") const override;
185
186private:
187
189 // The components.
190
191 ClassDefOverride(TVector3,3) // A 3D physics vector
192
193 // make TLorentzVector a friend class
195};
196
198// Addition of 3-vectors.
199
201// Subtraction of 3-vectors.
202
204// Scalar product of 3-vectors.
205
208// Scaling of 3-vectors with a real number
209
211
212
213inline Double_t & TVector3::operator[] (int i) { return operator()(i); }
214inline Double_t TVector3::operator[] (int i) const { return operator()(i); }
215
216inline Double_t TVector3::x() const { return fX; }
217inline Double_t TVector3::y() const { return fY; }
218inline Double_t TVector3::z() const { return fZ; }
219inline Double_t TVector3::X() const { return fX; }
220inline Double_t TVector3::Y() const { return fY; }
221inline Double_t TVector3::Z() const { return fZ; }
222inline Double_t TVector3::Px() const { return fX; }
223inline Double_t TVector3::Py() const { return fY; }
224inline Double_t TVector3::Pz() const { return fZ; }
225
226inline void TVector3::SetX(Double_t xx) { fX = xx; }
227inline void TVector3::SetY(Double_t yy) { fY = yy; }
228inline void TVector3::SetZ(Double_t zz) { fZ = zz; }
229
231 fX = xx;
232 fY = yy;
233 fZ = zz;
234}
235
236inline void TVector3::GetXYZ(Double_t *carray) const {
237 carray[0] = fX;
238 carray[1] = fY;
239 carray[2] = fZ;
240}
241
242inline void TVector3::GetXYZ(Float_t *carray) const {
243 carray[0] = fX;
244 carray[1] = fY;
245 carray[2] = fZ;
246}
247
248////////////////////////////////////////////////////////////////////////////////
249/// Constructors
251: fX(0.0), fY(0.0), fZ(0.0) {}
252
254 fX(p.fX), fY(p.fY), fZ(p.fZ) {}
255
257: fX(xx), fY(yy), fZ(zz) {}
258
259inline TVector3::TVector3(const Double_t * x0)
260: fX(x0[0]), fY(x0[1]), fZ(x0[2]) {}
261
262inline TVector3::TVector3(const Float_t * x0)
263: fX(x0[0]), fY(x0[1]), fZ(x0[2]) {}
264
265
266inline Double_t TVector3::operator () (int i) const {
267 switch(i) {
268 case 0:
269 return fX;
270 case 1:
271 return fY;
272 case 2:
273 return fZ;
274 default:
275 Error("operator()(i)", "bad index (%d) returning 0",i);
276 }
277 return 0.;
278}
279
281 switch(i) {
282 case 0:
283 return fX;
284 case 1:
285 return fY;
286 case 2:
287 return fZ;
288 default:
289 Error("operator()(i)", "bad index (%d) returning &fX",i);
290 }
291 return fX;
292}
293
295 fX = p.fX;
296 fY = p.fY;
297 fZ = p.fZ;
298 return *this;
299}
300
302 return (v.fX==fX && v.fY==fY && v.fZ==fZ) ? kTRUE : kFALSE;
303}
304
306 return (v.fX!=fX || v.fY!=fY || v.fZ!=fZ) ? kTRUE : kFALSE;
307}
308
310 fX += p.fX;
311 fY += p.fY;
312 fZ += p.fZ;
313 return *this;
314}
315
317 fX -= p.fX;
318 fY -= p.fY;
319 fZ -= p.fZ;
320 return *this;
321}
322
324 return TVector3(-fX, -fY, -fZ);
325}
326
328 fX *= a;
329 fY *= a;
330 fZ *= a;
331 return *this;
332}
333
334inline Double_t TVector3::Dot(const TVector3 & p) const {
335 return fX*p.fX + fY*p.fY + fZ*p.fZ;
336}
337
338inline TVector3 TVector3::Cross(const TVector3 & p) const {
339 return TVector3(fY*p.fZ-p.fY*fZ, fZ*p.fX-p.fZ*fX, fX*p.fY-p.fX*fY);
340}
341
342inline Double_t TVector3::Mag2() const { return fX*fX + fY*fY + fZ*fZ; }
343
344
346 Double_t xx = fX < 0.0 ? -fX : fX;
347 Double_t yy = fY < 0.0 ? -fY : fY;
348 Double_t zz = fZ < 0.0 ? -fZ : fZ;
349 if (xx < yy) {
350 return xx < zz ? TVector3(0,fZ,-fY) : TVector3(fY,-fX,0);
351 } else {
352 return yy < zz ? TVector3(-fZ,0,fX) : TVector3(fY,-fX,0);
353 }
354}
355
356inline Double_t TVector3::Perp2() const { return fX*fX + fY*fY; }
357
358
359inline Double_t TVector3::Pt() const { return Perp(); }
360
361inline Double_t TVector3::Perp2(const TVector3 & p) const {
362 Double_t tot = p.Mag2();
363 Double_t ss = Dot(p);
364 Double_t per = Mag2();
365 if (tot > 0.0) per -= ss*ss/tot;
366 if (per < 0) per = 0;
367 return per;
368}
369
370inline Double_t TVector3::Pt(const TVector3 & p) const {
371 return Perp(p);
372}
373
375 Double_t ptot = Mag();
376 return ptot == 0.0 ? 1.0 : fZ/ptot;
377}
378
380 Double_t factor = Mag();
381 if (factor == 0) {
382 Warning("SetMag","zero vector can't be stretched");
383 } else {
384 factor = ma/factor;
385 SetX(fX*factor);
386 SetY(fY*factor);
387 SetZ(fZ*factor);
388 }
389}
390
392 Double_t p = Perp();
393 if (p != 0.0) {
394 fX *= r/p;
395 fY *= r/p;
396 }
397}
398
399inline Double_t TVector3::DeltaPhi(const TVector3 & v) const {
400 return TVector2::Phi_mpi_pi(Phi()-v.Phi());
401}
402
403inline Double_t TVector3::Eta() const {
404 return PseudoRapidity();
405}
406
407inline Double_t TVector3::DrEtaPhi(const TVector3 & v) const{
408 return DeltaR(v);
409}
410
411
413 return TVector2 (Eta(),Phi());
414}
415
417 return TVector2(fX,fY);
418}
419
420
421#endif
#define a(i)
Definition RSha256.hxx:99
bool Bool_t
Boolean (0=false, 1=true) (bool)
Definition RtypesCore.h:77
float Float_t
Float 4 bytes (float)
Definition RtypesCore.h:71
constexpr Bool_t kFALSE
Definition RtypesCore.h:108
double Double_t
Double 8 bytes.
Definition RtypesCore.h:73
constexpr Bool_t kTRUE
Definition RtypesCore.h:107
const char Option_t
Option string (const char)
Definition RtypesCore.h:80
#define ClassDefOverride(name, id)
Definition Rtypes.h:348
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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
TMatrixT.
Definition TMatrixT.h:40
Mother of all ROOT objects.
Definition TObject.h:41
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:1074
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1088
<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:103
Double_t Z() const
Definition TVector3.h:221
void SetMagThetaPhi(Double_t mag, Double_t theta, Double_t phi)
Setter with mag, theta, phi.
Definition TVector3.cxx:407
void SetY(Double_t)
Definition TVector3.h:227
void SetXYZ(Double_t x, Double_t y, Double_t z)
Definition TVector3.h:230
Double_t Eta() const
Definition TVector3.h:403
TVector3 operator-() const
Definition TVector3.h:323
void SetPhi(Double_t)
Set phi keeping mag and theta constant (BaBar).
Definition TVector3.cxx:387
Double_t Px() const
Definition TVector3.h:222
Double_t fZ
Definition TVector3.h:188
Double_t fX
Definition TVector3.h:188
void Print(Option_t *option="") const override
Print vector parameters.
Definition TVector3.cxx:485
Double_t Phi() const
Return the azimuth angle. Returns phi from -pi to pi.
Definition TVector3.cxx:235
Double_t Y() const
Definition TVector3.h:220
Double_t Py() const
Definition TVector3.h:223
Double_t Dot(const TVector3 &) const
Definition TVector3.h:334
void RotateZ(Double_t)
Rotate vector around Z.
Definition TVector3.cxx:284
void SetPerp(Double_t)
Definition TVector3.h:391
Double_t Mag2() const
Definition TVector3.h:342
Bool_t operator!=(const TVector3 &) const
Definition TVector3.h:305
Double_t X() const
Definition TVector3.h:219
Double_t operator[](int) const
Definition TVector3.h:214
TVector3 Unit() const
Return unit vector parallel to this.
Definition TVector3.cxx:251
void RotateX(Double_t)
Rotate vector around X.
Definition TVector3.cxx:262
TVector3 & operator+=(const TVector3 &)
Definition TVector3.h:309
constexpr std::size_t size() const
The length is always 3. For compatibility with the standard library.
Definition TVector3.h:44
Double_t Angle(const TVector3 &) const
Return the angle w.r.t. another 3-vector.
Definition TVector3.cxx:202
Double_t Pz() const
Definition TVector3.h:224
Double_t PseudoRapidity() const
Double_t m = Mag(); return 0.5*log( (m+fZ)/(m-fZ) ); guard against Pt=0.
Definition TVector3.cxx:345
Double_t operator()(int) const
Definition TVector3.h:266
void Rotate(Double_t, const TVector3 &)
Rotate vector.
Definition TVector3.cxx:295
void SetPtEtaPhi(Double_t pt, Double_t eta, Double_t phi)
Set Pt, Eta and Phi.
Definition TVector3.cxx:357
void GetXYZ(Double_t *carray) const
Definition TVector3.h:236
TVector3 & operator*=(Double_t)
Definition TVector3.h:327
TVector2 EtaPhiVector() const
Definition TVector3.h:412
Double_t Scalar
Definition TVector3.h:26
Double_t CosTheta() const
Definition TVector3.h:374
Double_t z() const
Definition TVector3.h:218
Double_t DrEtaPhi(const TVector3 &) const
Definition TVector3.h:407
TVector3 & operator=(const TVector3 &)
Definition TVector3.h:294
TVector3 Orthogonal() const
Definition TVector3.h:345
Double_t x() const
Definition TVector3.h:216
Double_t Mag() const
Definition TVector3.h:89
Double_t Theta() const
Return the polar angle.
Definition TVector3.cxx:243
Double_t fY
Definition TVector3.h:188
void SetMag(Double_t)
Definition TVector3.h:379
Bool_t operator==(const TVector3 &) const
Definition TVector3.h:301
Double_t Perp() const
Return the transverse component (R in cylindrical coordinate system)
Definition TVector3.cxx:218
~TVector3() override=default
TVector3 & operator-=(const TVector3 &)
Definition TVector3.h:316
TVector3()
Constructors.
Definition TVector3.h:250
Double_t y() const
Definition TVector3.h:217
TVector3 & Transform(const TRotation &)
Transform this vector with a TRotation.
Definition TVector3.cxx:195
Double_t DeltaPhi(const TVector3 &) const
Definition TVector3.h:399
void SetX(Double_t)
Definition TVector3.h:226
TVector2 XYvector() const
Definition TVector3.h:416
Double_t Perp2() const
Definition TVector3.h:356
TVector3 Cross(const TVector3 &) const
Definition TVector3.h:338
void RotateY(Double_t)
Rotate vector around Y.
Definition TVector3.cxx:273
void SetZ(Double_t)
Definition TVector3.h:228
Double_t Pt() const
Definition TVector3.h:359
void SetPtThetaPhi(Double_t pt, Double_t theta, Double_t phi)
Set Pt, Theta and Phi.
Definition TVector3.cxx:365
void RotateUz(const TVector3 &)
Express the coordinates of this vector in a new reference frame S' whose z' axis will lie in the dire...
Definition TVector3.cxx:324
Double_t DeltaR(const TVector3 &) const
Return deltaR with respect to v.
Definition TVector3.cxx:397
void SetTheta(Double_t)
Set theta keeping mag and phi constant (BaBar).
Definition TVector3.cxx:375
TPaveText * pt
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:673