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