ROOT  6.06/09
Reference Guide
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 #ifndef ROOT_TError
15 #include "TError.h"
16 #endif
17 #ifndef ROOT_TVector2
18 #include "TVector2.h"
19 #endif
20 #ifndef ROOT_TMatrix
21 #include "TMatrix.h"
22 #endif
23 
24 class TRotation;
25 
26 
27 class TVector3 : public TObject {
28 
29 public:
30 
31  typedef Double_t Scalar; // to be able to use it with the ROOT::Math::VectorUtil functions
32 
33  TVector3();
34 
36  // The constructor.
37 
38  TVector3(const Double_t *);
39  TVector3(const Float_t *);
40  // Constructors from an array
41 
42  TVector3(const TVector3 &);
43  // The copy constructor.
44 
45  virtual ~TVector3();
46  // Destructor
47 
48  Double_t operator () (int) const;
49  inline Double_t operator [] (int) const;
50  // Get components by index (Geant4).
51 
52  Double_t & operator () (int);
53  inline Double_t & operator [] (int);
54  // Set components by index.
55 
56  inline Double_t x() const;
57  inline Double_t y() const;
58  inline Double_t z() const;
59  inline Double_t X() const;
60  inline Double_t Y() const;
61  inline Double_t Z() const;
62  inline Double_t Px() const;
63  inline Double_t Py() const;
64  inline Double_t Pz() const;
65  // The components in cartesian coordinate system.
66 
67  inline void SetX(Double_t);
68  inline void SetY(Double_t);
69  inline void SetZ(Double_t);
70  inline void SetXYZ(Double_t x, Double_t y, Double_t z);
71  void SetPtEtaPhi(Double_t pt, Double_t eta, Double_t phi);
72  void SetPtThetaPhi(Double_t pt, Double_t theta, Double_t phi);
73 
74  inline void GetXYZ(Double_t *carray) const;
75  inline void GetXYZ(Float_t *carray) const;
76  // Get the components into an array
77  // not checked!
78 
79  Double_t Phi() const;
80  // The azimuth angle. returns phi from -pi to pi
81 
82  Double_t Theta() const;
83  // The polar angle.
84 
85  inline Double_t CosTheta() const;
86  // Cosine of the polar angle.
87 
88  inline Double_t Mag2() const;
89  // The magnitude squared (rho^2 in spherical coordinate system).
90 
91  Double_t Mag() const;
92  // The magnitude (rho in spherical coordinate system).
93 
94  void SetPhi(Double_t);
95  // Set phi keeping mag and theta constant (BaBar).
96 
97  void SetTheta(Double_t);
98  // Set theta keeping mag and phi constant (BaBar).
99 
100  inline void SetMag(Double_t);
101  // Set magnitude keeping theta and phi constant (BaBar).
102 
103  inline Double_t Perp2() const;
104  // The transverse component squared (R^2 in cylindrical coordinate system).
105 
106  inline Double_t Pt() const;
107  Double_t Perp() const;
108  // The transverse component (R in cylindrical coordinate system).
109 
110  inline void SetPerp(Double_t);
111  // Set the transverse component keeping phi and z constant.
112 
113  inline Double_t Perp2(const TVector3 &) const;
114  // The transverse component w.r.t. given axis squared.
115 
116  inline Double_t Pt(const TVector3 &) const;
117  Double_t Perp(const TVector3 &) const;
118  // The transverse component w.r.t. given axis.
119 
120  inline Double_t DeltaPhi(const TVector3 &) const;
121  Double_t DeltaR(const TVector3 &) const;
122  inline Double_t DrEtaPhi(const TVector3 &) const;
123  inline TVector2 EtaPhiVector() const;
124  void SetMagThetaPhi(Double_t mag, Double_t theta, Double_t phi);
125 
126  inline TVector3 & operator = (const TVector3 &);
127  // Assignment.
128 
129  inline Bool_t operator == (const TVector3 &) const;
130  inline Bool_t operator != (const TVector3 &) const;
131  // Comparisons (Geant4).
132 
133  inline TVector3 & operator += (const TVector3 &);
134  // Addition.
135 
136  inline TVector3 & operator -= (const TVector3 &);
137  // Subtraction.
138 
139  inline TVector3 operator - () const;
140  // Unary minus.
141 
142  inline TVector3 & operator *= (Double_t);
143  // Scaling with real numbers.
144 
145  TVector3 Unit() const;
146  // Unit vector parallel to this.
147 
148  inline TVector3 Orthogonal() const;
149  // Vector orthogonal to this (Geant4).
150 
151  inline Double_t Dot(const TVector3 &) const;
152  // Scalar product.
153 
154  inline TVector3 Cross(const TVector3 &) const;
155  // Cross product.
156 
157  Double_t Angle(const TVector3 &) const;
158  // The angle w.r.t. another 3-vector.
159 
160  Double_t PseudoRapidity() const;
161  // Returns the pseudo-rapidity, i.e. -ln(tan(theta/2))
162 
163  inline Double_t Eta() const;
164 
165  void RotateX(Double_t);
166  // Rotates the Hep3Vector around the x-axis.
167 
168  void RotateY(Double_t);
169  // Rotates the Hep3Vector around the y-axis.
170 
171  void RotateZ(Double_t);
172  // Rotates the Hep3Vector around the z-axis.
173 
174  void RotateUz(const TVector3&);
175  // Rotates reference frame from Uz to newUz (unit vector) (Geant4).
176 
177  void Rotate(Double_t, const TVector3 &);
178  // Rotates around the axis specified by another Hep3Vector.
179 
180  TVector3 & operator *= (const TRotation &);
181  TVector3 & Transform(const TRotation &);
182  // Transformation with a Rotation matrix.
183 
184  inline TVector2 XYvector() const;
185 
186  void Print(Option_t* option="") const;
187 
188 private:
189 
191  // The components.
192 
193  ClassDef(TVector3,3) // A 3D physics vector
194 
195 };
196 
197 
198 TVector3 operator + (const TVector3 &, const TVector3 &);
199 // Addition of 3-vectors.
200 
201 TVector3 operator - (const TVector3 &, const TVector3 &);
202 // Subtraction of 3-vectors.
203 
204 Double_t operator * (const TVector3 &, const TVector3 &);
205 // Scalar product of 3-vectors.
206 
209 // Scaling of 3-vectors with a real number
210 
211 TVector3 operator * (const TMatrix &, const TVector3 &);
212 
213 
215 Double_t TVector3::operator[] (int i) const { return operator()(i); }
216 
217 inline Double_t TVector3::x() const { return fX; }
218 inline Double_t TVector3::y() const { return fY; }
219 inline Double_t TVector3::z() const { return fZ; }
220 inline Double_t TVector3::X() const { return fX; }
221 inline Double_t TVector3::Y() const { return fY; }
222 inline Double_t TVector3::Z() const { return fZ; }
223 inline Double_t TVector3::Px() const { return fX; }
224 inline Double_t TVector3::Py() const { return fY; }
225 inline Double_t TVector3::Pz() const { return fZ; }
226 
227 inline void TVector3::SetX(Double_t xx) { fX = xx; }
228 inline void TVector3::SetY(Double_t yy) { fY = yy; }
229 inline void TVector3::SetZ(Double_t zz) { fZ = zz; }
230 
231 inline void TVector3::SetXYZ(Double_t xx, Double_t yy, Double_t zz) {
232  fX = xx;
233  fY = yy;
234  fZ = zz;
235 }
236 
237 inline void TVector3::GetXYZ(Double_t *carray) const {
238  carray[0] = fX;
239  carray[1] = fY;
240  carray[2] = fZ;
241 }
242 
243 inline void TVector3::GetXYZ(Float_t *carray) const {
244  carray[0] = fX;
245  carray[1] = fY;
246  carray[2] = fZ;
247 }
248 
249 
251  fX = p.fX;
252  fY = p.fY;
253  fZ = p.fZ;
254  return *this;
255 }
256 
257 inline Bool_t TVector3::operator == (const TVector3& v) const {
258  return (v.fX==fX && v.fY==fY && v.fZ==fZ) ? kTRUE : kFALSE;
259 }
260 
261 inline Bool_t TVector3::operator != (const TVector3& v) const {
262  return (v.fX!=fX || v.fY!=fY || v.fZ!=fZ) ? kTRUE : kFALSE;
263 }
264 
266  fX += p.fX;
267  fY += p.fY;
268  fZ += p.fZ;
269  return *this;
270 }
271 
273  fX -= p.fX;
274  fY -= p.fY;
275  fZ -= p.fZ;
276  return *this;
277 }
278 
280  return TVector3(-fX, -fY, -fZ);
281 }
282 
284  fX *= a;
285  fY *= a;
286  fZ *= a;
287  return *this;
288 }
289 
290 inline Double_t TVector3::Dot(const TVector3 & p) const {
291  return fX*p.fX + fY*p.fY + fZ*p.fZ;
292 }
293 
294 inline TVector3 TVector3::Cross(const TVector3 & p) const {
295  return TVector3(fY*p.fZ-p.fY*fZ, fZ*p.fX-p.fZ*fX, fX*p.fY-p.fX*fY);
296 }
297 
298 inline Double_t TVector3::Mag2() const { return fX*fX + fY*fY + fZ*fZ; }
299 
300 
302  Double_t xx = fX < 0.0 ? -fX : fX;
303  Double_t yy = fY < 0.0 ? -fY : fY;
304  Double_t zz = fZ < 0.0 ? -fZ : fZ;
305  if (xx < yy) {
306  return xx < zz ? TVector3(0,fZ,-fY) : TVector3(fY,-fX,0);
307  } else {
308  return yy < zz ? TVector3(-fZ,0,fX) : TVector3(fY,-fX,0);
309  }
310 }
311 
312 inline Double_t TVector3::Perp2() const { return fX*fX + fY*fY; }
313 
314 
315 inline Double_t TVector3::Pt() const { return Perp(); }
316 
317 inline Double_t TVector3::Perp2(const TVector3 & p) const {
318  Double_t tot = p.Mag2();
319  Double_t ss = Dot(p);
320  Double_t per = Mag2();
321  if (tot > 0.0) per -= ss*ss/tot;
322  if (per < 0) per = 0;
323  return per;
324 }
325 
326 inline Double_t TVector3::Pt(const TVector3 & p) const {
327  return Perp(p);
328 }
329 
330 inline Double_t TVector3::CosTheta() const {
331  Double_t ptot = Mag();
332  return ptot == 0.0 ? 1.0 : fZ/ptot;
333 }
334 
335 inline void TVector3::SetMag(Double_t ma) {
336  Double_t factor = Mag();
337  if (factor == 0) {
338  Warning("SetMag","zero vector can't be stretched");
339  } else {
340  factor = ma/factor;
341  SetX(fX*factor);
342  SetY(fY*factor);
343  SetZ(fZ*factor);
344  }
345 }
346 
348  Double_t p = Perp();
349  if (p != 0.0) {
350  fX *= r/p;
351  fY *= r/p;
352  }
353 }
354 
355 inline Double_t TVector3::DeltaPhi(const TVector3 & v) const {
356  return TVector2::Phi_mpi_pi(Phi()-v.Phi());
357 }
358 
359 inline Double_t TVector3::Eta() const {
360  return PseudoRapidity();
361 }
362 
363 inline Double_t TVector3::DrEtaPhi(const TVector3 & v) const{
364  return DeltaR(v);
365 }
366 
367 
369  return TVector2 (Eta(),Phi());
370 }
371 
372 inline TVector2 TVector3::XYvector() const {
373  return TVector2(fX,fY);
374 }
375 
376 #endif
TVector3 & operator=(const TVector3 &)
Definition: TVector3.h:250
Double_t y() const
Definition: TVector3.h:218
void RotateUz(const TVector3 &)
NewUzVector must be normalized !
Definition: TVector3.cxx:349
virtual ~TVector3()
Definition: TVector3.cxx:186
Double_t Phi() const
return the azimuth angle. returns phi from -pi to pi
Definition: TVector3.cxx:280
Double_t Z() const
Definition: TVector3.h:222
Double_t z() const
Definition: TVector3.h:219
void SetY(Double_t)
Definition: TVector3.h:228
void SetPerp(Double_t)
Definition: TVector3.h:347
float Float_t
Definition: RtypesCore.h:53
const char Option_t
Definition: RtypesCore.h:62
Double_t Theta() const
return the polar angle
Definition: TVector3.cxx:288
Bool_t operator!=(const TVector3 &) const
Definition: TVector3.h:261
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
Double_t operator*(const TVector3 &, const TVector3 &)
Definition: TVector3.cxx:481
const Bool_t kFALSE
Definition: Rtypes.h:92
Double_t fY
Definition: TVector3.h:190
Double_t Dot(const TVector3 &) const
Definition: TVector3.h:290
Double_t Perp() const
return the transverse component (R in cylindrical coordinate system)
Definition: TVector3.cxx:263
void RotateX(Double_t)
rotate vector around X
Definition: TVector3.cxx:307
TVector2 EtaPhiVector() const
Definition: TVector3.h:368
Double_t Mag() const
return the magnitude (rho in spherical coordinate system)
Definition: TVector3.cxx:255
Double_t Angle(const TVector3 &) const
return the angle w.r.t. another 3-vector
Definition: TVector3.cxx:239
Double_t PseudoRapidity() const
Double_t m = Mag(); return 0.5*log( (m+fZ)/(m-fZ) ); guard against Pt=0.
Definition: TVector3.cxx:370
#define ClassDef(name, id)
Definition: Rtypes.h:254
Double_t fZ
Definition: TVector3.h:190
Double_t DeltaPhi(const TVector3 &) const
Definition: TVector3.h:355
TVector3 & operator*=(Double_t)
Definition: TVector3.h:283
TVector3 operator-(const TVector3 &, const TVector3 &)
Definition: TVector3.cxx:469
TVector3 Unit() const
return unit vector parallel to this.
Definition: TVector3.cxx:296
void Rotate(Double_t, const TVector3 &)
rotate vector
Definition: TVector3.cxx:340
static Double_t Phi_mpi_pi(Double_t x)
(static function) returns phi angle in the interval [-PI,PI)
Definition: TVector2.cxx:95
void SetXYZ(Double_t x, Double_t y, Double_t z)
Definition: TVector3.h:231
TVector3 operator+(const TVector3 &, const TVector3 &)
Definition: TVector3.cxx:465
TVector3 operator-() const
Definition: TVector3.h:279
Double_t CosTheta() const
Definition: TVector3.h:330
void RotateY(Double_t)
rotate vector around Y
Definition: TVector3.cxx:318
Bool_t operator==(const TVector3 &) const
Definition: TVector3.h:257
TPaveText * pt
void SetTheta(Double_t)
Set theta keeping mag and phi constant (BaBar).
Definition: TVector3.cxx:400
ROOT::R::TRInterface & r
Definition: Object.C:4
Double_t operator[](int) const
Definition: TVector3.h:215
SVector< double, 2 > v
Definition: Dict.h:5
Double_t Perp2() const
Definition: TVector3.h:312
TVector3 & operator+=(const TVector3 &)
Definition: TVector3.h:265
Double_t X() const
Definition: TVector3.h:220
Double_t x() const
Definition: TVector3.h:217
TVector2 XYvector() const
Definition: TVector3.h:372
TVector3 & operator-=(const TVector3 &)
Definition: TVector3.h:272
Double_t Mag2() const
Definition: TVector3.h:298
Double_t Px() const
Definition: TVector3.h:223
Double_t fX
Definition: TVector3.h:190
double Double_t
Definition: RtypesCore.h:55
void SetZ(Double_t)
Definition: TVector3.h:229
Double_t Eta() const
Definition: TVector3.h:359
Double_t Pz() const
Definition: TVector3.h:225
TVector3 & Transform(const TRotation &)
transform this vector with a TRotation
Definition: TVector3.cxx:232
void RotateZ(Double_t)
rotate vector around Z
Definition: TVector3.cxx:329
void GetXYZ(Double_t *carray) const
Definition: TVector3.h:237
Mother of all ROOT objects.
Definition: TObject.h:58
void SetPhi(Double_t)
Set phi keeping mag and theta constant (BaBar).
Definition: TVector3.cxx:412
Double_t operator()(int) const
dereferencing operator const
Definition: TVector3.cxx:191
Double_t DeltaR(const TVector3 &) const
return deltaR with respect to v
Definition: TVector3.cxx:422
TVector3 Cross(const TVector3 &) const
Definition: TVector3.h:294
void SetMag(Double_t)
Definition: TVector3.h:335
Double_t Y() const
Definition: TVector3.h:221
void SetX(Double_t)
Definition: TVector3.h:227
Double_t Py() const
Definition: TVector3.h:224
void SetPtThetaPhi(Double_t pt, Double_t theta, Double_t phi)
set Pt, Theta and Phi
Definition: TVector3.cxx:390
const Bool_t kTRUE
Definition: Rtypes.h:91
Double_t DrEtaPhi(const TVector3 &) const
Definition: TVector3.h:363
TVector3 Orthogonal() const
Definition: TVector3.h:301
void Print(Option_t *option="") const
This method must be overridden when a class wants to print itself.
Definition: TVector3.cxx:496
void SetMagThetaPhi(Double_t mag, Double_t theta, Double_t phi)
setter with mag, theta, phi
Definition: TVector3.cxx:432
Double_t Pt() const
Definition: TVector3.h:315
Double_t Scalar
Definition: TVector3.h:31
void SetPtEtaPhi(Double_t pt, Double_t eta, Double_t phi)
set Pt, Eta and Phi
Definition: TVector3.cxx:382
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:904