ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TVector3.cxx
Go to the documentation of this file.
1 // @(#)root/physics:$Id$
2 // Author: Pasha Murat, Peter Malzacher 12/02/99
3 // Aug 11 1999: added Pt == 0 guard to Eta()
4 // Oct 8 1999: changed Warning to Error and
5 // return fX in Double_t & operator()
6 // Oct 20 1999: Bug fix: sign in PseudoRapidity
7 // Warning-> Error in Double_t operator()
8 
9 //______________________________________________________________________________
10 //*-*-*-*-*-*-*-*-*-*-*-*The Physics Vector package *-*-*-*-*-*-*-*-*-*-*-*
11 //*-* ========================== *
12 //*-* The Physics Vector package consists of five classes: *
13 //*-* - TVector2 *
14 //*-* - TVector3 *
15 //*-* - TRotation *
16 //*-* - TLorentzVector *
17 //*-* - TLorentzRotation *
18 //*-* It is a combination of CLHEPs Vector package written by *
19 //*-* Leif Lonnblad, Andreas Nilsson and Evgueni Tcherniaev *
20 //*-* and a ROOT package written by Pasha Murat. *
21 //*-* for CLHEP see: http://wwwinfo.cern.ch/asd/lhc++/clhep/ *
22 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
23 //BEGIN_HTML <!--
24 /* -->
25 <H2>
26 TVector3</H2>
27 <TT>TVector3</TT> is a general three vector class, which can be used for
28 the description of different vectors in 3D.
29 <H3>
30 Declaration / Access to the components</H3>
31 <TT>TVector3</TT> has been implemented as a vector of three <TT>Double_t</TT>
32 variables, representing the cartesian coordinates. By default all components
33 are initialized to zero:
34 
35 <P><TT>&nbsp; TVector3 v1;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp; //
36 v1 = (0,0,0)</TT>
37 <BR><TT>&nbsp; TVector3 v3(1,2,3); // v3 = (1,2,3)</TT>
38 <BR><TT>&nbsp; TVector3 v4(v2);&nbsp;&nbsp;&nbsp; // v4 = v2</TT>
39 
40 <P>It is also possible (but not recommended) to initialize a <TT>TVector3</TT>
41 with a <TT>Double_t</TT> or <TT>Float_t</TT> C array.
42 
43 <P>You can get the basic components either by name or by index using <TT>operator()</TT>:
44 
45 <P><TT>&nbsp; xx = v1.X();&nbsp;&nbsp;&nbsp; or&nbsp;&nbsp;&nbsp; xx =
46 v1(0);</TT>
47 <BR><TT>&nbsp; yy = v1.Y();&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;
48 yy = v1(1);</TT>
49 <BR><TT>&nbsp; zz = v1.Z();&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;
50 zz = v1(2);</TT>
51 
52 <P>The memberfunctions <TT>SetX()</TT>, <TT>SetY()</TT>, <TT>SetZ()</TT>
53 and<TT> SetXYZ()</TT> allow to set the components:
54 
55 <P><TT>&nbsp; v1.SetX(1.); v1.SetY(2.); v1.SetZ(3.);</TT>
56 <BR><TT>&nbsp; v1.SetXYZ(1.,2.,3.);</TT>
57 <BR>&nbsp;
58 <H3>
59 Noncartesian coordinates</H3>
60 To get information on the <TT>TVector3</TT> in spherical (rho,phi,theta)
61 or cylindrical (z,r,theta) coordinates, the
62 <BR>the member functions <TT>Mag()</TT> (=magnitude=rho in spherical coordinates),
63 <TT>Mag2()</TT>, <TT>Theta()</TT>, <TT>CosTheta()</TT>, <TT>Phi()</TT>,
64 <TT>Perp()</TT> (the transverse component=r in cylindrical coordinates),
65 <TT>Perp2()</TT> can be used:
66 
67 <P><TT>&nbsp; Double_t m&nbsp; = v.Mag();&nbsp;&nbsp;&nbsp; // get magnitude
68 (=rho=Sqrt(x*x+y*y+z*z)))</TT>
69 <BR><TT>&nbsp; Double_t m2 = v.Mag2();&nbsp;&nbsp; // get magnitude squared</TT>
70 <BR><TT>&nbsp; Double_t t&nbsp; = v.Theta();&nbsp; // get polar angle</TT>
71 <BR><TT>&nbsp; Double_t ct = v.CosTheta();// get cos of theta</TT>
72 <BR><TT>&nbsp; Double_t p&nbsp; = v.Phi();&nbsp;&nbsp;&nbsp; // get azimuth
73 angle</TT>
74 <BR><TT>&nbsp; Double_t pp = v.Perp();&nbsp;&nbsp; // get transverse component</TT>
75 <BR><TT>&nbsp; Double_t pp2= v.Perp2();&nbsp; // get transvers component
76 squared</TT>
77 
78 <P>It is also possible to get the transverse component with respect to
79 another vector:
80 
81 <P><TT>&nbsp; Double_t ppv1 = v.Perp(v1);</TT>
82 <BR><TT>&nbsp; Double_t pp2v1 = v.Perp2(v1);</TT>
83 
84 <P>The pseudo-rapidity ( eta=-ln (tan (theta/2)) ) can be obtained by <TT>Eta()</TT>
85 or <TT>PseudoRapidity()</TT>:
86 <BR>&nbsp;
87 <BR><TT>&nbsp; Double_t eta = v.PseudoRapidity();</TT>
88 
89 <P>There are set functions to change one of the noncartesian coordinates:
90 
91 <P><TT>&nbsp; v.SetTheta(.5); // keeping rho and phi</TT>
92 <BR><TT>&nbsp; v.SetPhi(.8);&nbsp;&nbsp; // keeping rho and theta</TT>
93 <BR><TT>&nbsp; v.SetMag(10.);&nbsp; // keeping theta and phi</TT>
94 <BR><TT>&nbsp; v.SetPerp(3.);&nbsp; // keeping z and phi</TT>
95 <BR>&nbsp;
96 <H3>
97 Arithmetic / Comparison</H3>
98 The <TT>TVector3</TT> class provides the operators to add, subtract, scale and compare
99 vectors:
100 
101 <P><TT>&nbsp; v3&nbsp; = -v1;</TT>
102 <BR><TT>&nbsp; v1&nbsp; = v2+v3;</TT>
103 <BR><TT>&nbsp; v1 += v3;</TT>
104 <BR><TT>&nbsp; v1&nbsp; = v1 - v3</TT>
105 <BR><TT>&nbsp; v1 -= v3;</TT>
106 <BR><TT>&nbsp; v1 *= 10;</TT>
107 <BR><TT>&nbsp; v1&nbsp; = 5*v2;</TT>
108 
109 <P><TT>&nbsp; if(v1==v2) {...}</TT>
110 <BR><TT>&nbsp; if(v1!=v2) {...}</TT>
111 <BR>&nbsp;
112 <H3>
113 Related Vectors</H3>
114 <TT>&nbsp; v2 = v1.Unit();&nbsp;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp; // get unit
115 vector parallel to v1</TT>
116 <BR><TT>&nbsp; v2 = v1.Orthogonal(); // get vector orthogonal to v1</TT>
117 <H3>
118 Scalar and vector products</H3>
119 <TT>&nbsp; s = v1.Dot(v2);&nbsp;&nbsp; // scalar product</TT>
120 <BR><TT>&nbsp; s = v1 * v2;&nbsp;&nbsp;&nbsp;&nbsp;&nbsp; // scalar product</TT>
121 <BR><TT>&nbsp; v = v1.Cross(v2); // vector product</TT>
122 <H3>
123 &nbsp;Angle between two vectors</H3>
124 <TT>&nbsp; Double_t a = v1.Angle(v2);</TT>
125 <H3>
126 Rotations</H3>
127 
128 <H5>
129 Rotation around axes</H5>
130 <TT>&nbsp; v.RotateX(.5);</TT>
131 <BR><TT>&nbsp; v.RotateY(TMath::Pi());</TT>
132 <BR><TT>&nbsp; v.RotateZ(angle);</TT>
133 <H5>
134 Rotation around a vector</H5>
135 <TT>&nbsp; v1.Rotate(TMath::Pi()/4, v2); // rotation around v2</TT>
136 <H5>
137 Rotation by TRotation</H5>
138 <TT>TVector3</TT> objects can be rotated by objects of the <TT>TRotation</TT>
139 class using the <TT>Transform()</TT> member functions,
140 <BR>the <TT>operator *=</TT> or the <TT>operator *</TT> of the TRotation
141 class:
142 
143 <P><TT>&nbsp; TRotation m;</TT>
144 <BR><TT>&nbsp; ...</TT>
145 <BR><TT>&nbsp; v1.transform(m);</TT>
146 <BR><TT>&nbsp; v1 = m*v1;</TT>
147 <BR><TT>&nbsp; v1 *= m; // Attention v1 = m*v1</TT>
148 <H5>
149 Transformation from rotated frame</H5>
150 <TT>&nbsp; TVector3 direction = v.Unit()</TT>
151 <BR><TT>&nbsp; v1.RotateUz(direction); // direction must be TVector3 of
152 unit length</TT>
153 
154 <P>transforms v1 from the rotated frame (z' parallel to direction, x' in
155 the theta plane and y' in the xy plane as well as perpendicular to the
156 theta plane) to the (x,y,z) frame.
157 
158 <!--*/
159 // -->END_HTML
160 //
161 
162 #include "TVector3.h"
163 
164 #include "TBuffer.h"
165 #include "TRotation.h"
166 #include "TMath.h"
167 #include "TClass.h"
168 
170 
171 ////////////////////////////////////////////////////////////////////////////////
172 
174 : fX(0.0), fY(0.0), fZ(0.0) {}
175 
177  fX(p.fX), fY(p.fY), fZ(p.fZ) {}
178 
180 : fX(xx), fY(yy), fZ(zz) {}
181 
183 : fX(x0[0]), fY(x0[1]), fZ(x0[2]) {}
184 
186 : fX(x0[0]), fY(x0[1]), fZ(x0[2]) {}
187 
189 
190 ////////////////////////////////////////////////////////////////////////////////
191 ///dereferencing operator const
192 
194  switch(i) {
195  case 0:
196  return fX;
197  case 1:
198  return fY;
199  case 2:
200  return fZ;
201  default:
202  Error("operator()(i)", "bad index (%d) returning 0",i);
203  }
204  return 0.;
205 }
206 
207 ////////////////////////////////////////////////////////////////////////////////
208 ///dereferencing operator
209 
211  switch(i) {
212  case 0:
213  return fX;
214  case 1:
215  return fY;
216  case 2:
217  return fZ;
218  default:
219  Error("operator()(i)", "bad index (%d) returning &fX",i);
220  }
221  return fX;
222 }
223 
224 ////////////////////////////////////////////////////////////////////////////////
225 ///multiplication operator
226 
228  return *this = m * (*this);
229 }
230 
231 ////////////////////////////////////////////////////////////////////////////////
232 ///transform this vector with a TRotation
233 
235  return *this = m * (*this);
236 }
237 
238 ////////////////////////////////////////////////////////////////////////////////
239 /// return the angle w.r.t. another 3-vector
240 
242 {
243  Double_t ptot2 = Mag2()*q.Mag2();
244  if(ptot2 <= 0) {
245  return 0.0;
246  } else {
247  Double_t arg = Dot(q)/TMath::Sqrt(ptot2);
248  if(arg > 1.0) arg = 1.0;
249  if(arg < -1.0) arg = -1.0;
250  return TMath::ACos(arg);
251  }
252 }
253 
254 ////////////////////////////////////////////////////////////////////////////////
255 /// return the magnitude (rho in spherical coordinate system)
256 
258 {
259  return TMath::Sqrt(Mag2());
260 }
261 
262 ////////////////////////////////////////////////////////////////////////////////
263 ///return the transverse component (R in cylindrical coordinate system)
264 
266 {
267  return TMath::Sqrt(Perp2());
268 }
269 
270 
271 ////////////////////////////////////////////////////////////////////////////////
272 ///return the transverse component (R in cylindrical coordinate system)
273 
275 {
276  return TMath::Sqrt(Perp2(p));
277 }
278 
279 ////////////////////////////////////////////////////////////////////////////////
280 ///return the azimuth angle. returns phi from -pi to pi
281 
283 {
284  return fX == 0.0 && fY == 0.0 ? 0.0 : TMath::ATan2(fY,fX);
285 }
286 
287 ////////////////////////////////////////////////////////////////////////////////
288 ///return the polar angle
289 
291 {
292  return fX == 0.0 && fY == 0.0 && fZ == 0.0 ? 0.0 : TMath::ATan2(Perp(),fZ);
293 }
294 
295 ////////////////////////////////////////////////////////////////////////////////
296 /// return unit vector parallel to this.
297 
299 {
300  Double_t tot2 = Mag2();
301  Double_t tot = (tot2 > 0) ? 1.0/TMath::Sqrt(tot2) : 1.0;
302  TVector3 p(fX*tot,fY*tot,fZ*tot);
303  return p;
304 }
305 
306 ////////////////////////////////////////////////////////////////////////////////
307 ///rotate vector around X
308 
310  Double_t s = TMath::Sin(angle);
311  Double_t c = TMath::Cos(angle);
312  Double_t yy = fY;
313  fY = c*yy - s*fZ;
314  fZ = s*yy + c*fZ;
315 }
316 
317 ////////////////////////////////////////////////////////////////////////////////
318 ///rotate vector around Y
319 
321  Double_t s = TMath::Sin(angle);
322  Double_t c = TMath::Cos(angle);
323  Double_t zz = fZ;
324  fZ = c*zz - s*fX;
325  fX = s*zz + c*fX;
326 }
327 
328 ////////////////////////////////////////////////////////////////////////////////
329 ///rotate vector around Z
330 
332  Double_t s = TMath::Sin(angle);
333  Double_t c = TMath::Cos(angle);
334  Double_t xx = fX;
335  fX = c*xx - s*fY;
336  fY = s*xx + c*fY;
337 }
338 
339 ////////////////////////////////////////////////////////////////////////////////
340 ///rotate vector
341 
342 void TVector3::Rotate(Double_t angle, const TVector3 & axis){
343  TRotation trans;
344  trans.Rotate(angle, axis);
345  operator*=(trans);
346 }
347 
348 ////////////////////////////////////////////////////////////////////////////////
349 /// NewUzVector must be normalized !
350 
351 void TVector3::RotateUz(const TVector3& NewUzVector) {
352  Double_t u1 = NewUzVector.fX;
353  Double_t u2 = NewUzVector.fY;
354  Double_t u3 = NewUzVector.fZ;
355  Double_t up = u1*u1 + u2*u2;
356 
357  if (up) {
358  up = TMath::Sqrt(up);
359  Double_t px = fX, py = fY, pz = fZ;
360  fX = (u1*u3*px - u2*py + u1*up*pz)/up;
361  fY = (u2*u3*px + u1*py + u2*up*pz)/up;
362  fZ = (u3*u3*px - px + u3*up*pz)/up;
363  } else if (u3 < 0.) { fX = -fX; fZ = -fZ; } // phi=0 teta=pi
364  else {};
365 }
366 
367 ////////////////////////////////////////////////////////////////////////////////
368 ///Double_t m = Mag();
369 ///return 0.5*log( (m+fZ)/(m-fZ) );
370 /// guard against Pt=0
371 
373  double cosTheta = CosTheta();
374  if (cosTheta*cosTheta < 1) return -0.5* TMath::Log( (1.0-cosTheta)/(1.0+cosTheta) );
375  if (fZ == 0) return 0;
376  Warning("PseudoRapidity","transvers momentum = 0! return +/- 10e10");
377  if (fZ > 0) return 10e10;
378  else return -10e10;
379 }
380 
381 ////////////////////////////////////////////////////////////////////////////////
382 ///set Pt, Eta and Phi
383 
385  Double_t apt = TMath::Abs(pt);
386  SetXYZ(apt*TMath::Cos(phi), apt*TMath::Sin(phi), apt/TMath::Tan(2.0*TMath::ATan(TMath::Exp(-eta))) );
387 }
388 
389 ////////////////////////////////////////////////////////////////////////////////
390 ///set Pt, Theta and Phi
391 
393  fX = pt * TMath::Cos(phi);
394  fY = pt * TMath::Sin(phi);
395  Double_t tanTheta = TMath::Tan(theta);
396  fZ = tanTheta ? pt / tanTheta : 0;
397 }
398 
399 ////////////////////////////////////////////////////////////////////////////////
400 /// Set theta keeping mag and phi constant (BaBar).
401 
403 {
404  Double_t ma = Mag();
405  Double_t ph = Phi();
406  SetX(ma*TMath::Sin(th)*TMath::Cos(ph));
407  SetY(ma*TMath::Sin(th)*TMath::Sin(ph));
408  SetZ(ma*TMath::Cos(th));
409 }
410 
411 ////////////////////////////////////////////////////////////////////////////////
412 /// Set phi keeping mag and theta constant (BaBar).
413 
415 {
416  Double_t xy = Perp();
417  SetX(xy*TMath::Cos(ph));
418  SetY(xy*TMath::Sin(ph));
419 }
420 
421 ////////////////////////////////////////////////////////////////////////////////
422 ///return deltaR with respect to v
423 
425 {
426  Double_t deta = Eta()-v.Eta();
427  Double_t dphi = TVector2::Phi_mpi_pi(Phi()-v.Phi());
428  return TMath::Sqrt( deta*deta+dphi*dphi );
429 }
430 
431 ////////////////////////////////////////////////////////////////////////////////
432 ///setter with mag, theta, phi
433 
435 {
436  Double_t amag = TMath::Abs(mag);
437  fX = amag * TMath::Sin(theta) * TMath::Cos(phi);
438  fY = amag * TMath::Sin(theta) * TMath::Sin(phi);
439  fZ = amag * TMath::Cos(theta);
440 }
441 
442 ////////////////////////////////////////////////////////////////////////////////
443 /// Stream an object of class TVector3.
444 
445 void TVector3::Streamer(TBuffer &R__b)
446 {
447  if (R__b.IsReading()) {
448  UInt_t R__s, R__c;
449  Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
450  if (R__v > 2) {
451  R__b.ReadClassBuffer(TVector3::Class(), this, R__v, R__s, R__c);
452  return;
453  }
454  //====process old versions before automatic schema evolution
455  if (R__v < 2) TObject::Streamer(R__b);
456  R__b >> fX;
457  R__b >> fY;
458  R__b >> fZ;
459  R__b.CheckByteCount(R__s, R__c, TVector3::IsA());
460  //====end of old versions
461 
462  } else {
463  R__b.WriteClassBuffer(TVector3::Class(),this);
464  }
465 }
466 
467 TVector3 operator + (const TVector3 & a, const TVector3 & b) {
468  return TVector3(a.X() + b.X(), a.Y() + b.Y(), a.Z() + b.Z());
469 }
470 
471 TVector3 operator - (const TVector3 & a, const TVector3 & b) {
472  return TVector3(a.X() - b.X(), a.Y() - b.Y(), a.Z() - b.Z());
473 }
474 
476  return TVector3(a*p.X(), a*p.Y(), a*p.Z());
477 }
478 
480  return TVector3(a*p.X(), a*p.Y(), a*p.Z());
481 }
482 
483 Double_t operator * (const TVector3 & a, const TVector3 & b) {
484  return a.Dot(b);
485 }
486 
487 TVector3 operator * (const TMatrix & m, const TVector3 & v ) {
488  return TVector3( m(0,0)*v.X()+m(0,1)*v.Y()+m(0,2)*v.Z(),
489  m(1,0)*v.X()+m(1,1)*v.Y()+m(1,2)*v.Z(),
490  m(2,0)*v.X()+m(2,1)*v.Y()+m(2,2)*v.Z());
491 }
492 
493 
494 //const TVector3 kXHat(1.0, 0.0, 0.0);
495 //const TVector3 kYHat(0.0, 1.0, 0.0);
496 //const TVector3 kZHat(0.0, 0.0, 1.0);
497 
499 {
500  //print vector parameters
501  Printf("%s %s (x,y,z)=(%f,%f,%f) (rho,theta,phi)=(%f,%f,%f)",GetName(),GetTitle(),X(),Y(),Z(),
503 }
void RotateUz(const TVector3 &)
NewUzVector must be normalized !
Definition: TVector3.cxx:351
virtual ~TVector3()
Definition: TVector3.cxx:188
Double_t Phi() const
return the azimuth angle. returns phi from -pi to pi
Definition: TVector3.cxx:282
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
Float_t pz
Definition: hprod.C:33
TVector3 operator*(const TVector3 &p, Double_t a)
Definition: TVector3.cxx:475
Double_t Z() const
Definition: TVector3.h:222
Bool_t IsReading() const
Definition: TBuffer.h:83
Double_t Log(Double_t x)
Definition: TMath.h:526
void SetY(Double_t)
Definition: TVector3.h:228
short Version_t
Definition: RtypesCore.h:61
float Float_t
Definition: RtypesCore.h:53
Float_t theta
Definition: shapesAnim.C:5
return c
const char Option_t
Definition: RtypesCore.h:62
Double_t Theta() const
return the polar angle
Definition: TVector3.cxx:290
Buffer base class used for serializing objects.
Definition: TBuffer.h:42
virtual Int_t CheckByteCount(UInt_t startpos, UInt_t bcnt, const TClass *clss)=0
Double_t RadToDeg()
Definition: TMath.h:49
TArc * a
Definition: textangle.C:12
Double_t fY
Definition: TVector3.h:190
Double_t Dot(const TVector3 &) const
Definition: TVector3.h:290
Float_t py
Definition: hprod.C:33
Double_t Perp() const
return the transverse component (R in cylindrical coordinate system)
Definition: TVector3.cxx:265
void RotateX(Double_t)
rotate vector around X
Definition: TVector3.cxx:309
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
Double_t Mag() const
return the magnitude (rho in spherical coordinate system)
Definition: TVector3.cxx:257
Double_t Angle(const TVector3 &) const
return the angle w.r.t. another 3-vector
Definition: TVector3.cxx:241
Double_t PseudoRapidity() const
Double_t m = Mag(); return 0.5*log( (m+fZ)/(m-fZ) ); guard against Pt=0.
Definition: TVector3.cxx:372
Double_t fZ
Definition: TVector3.h:190
void Class()
Definition: Class.C:29
TVector3 & operator*=(Double_t)
Definition: TVector3.h:283
TVector3 Unit() const
return unit vector parallel to this.
Definition: TVector3.cxx:298
void Rotate(Double_t, const TVector3 &)
rotate vector
Definition: TVector3.cxx:342
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
Double_t ATan2(Double_t, Double_t)
Definition: TMath.h:454
Double_t CosTheta() const
Definition: TVector3.h:330
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:918
TRotation & Rotate(Double_t, const TVector3 &)
Definition: TRotation.cxx:307
TVector3 operator-(const TVector3 &a, const TVector3 &b)
Definition: TVector3.cxx:471
void RotateY(Double_t)
rotate vector around Y
Definition: TVector3.cxx:320
TPaveText * pt
void SetTheta(Double_t)
Set theta keeping mag and phi constant (BaBar).
Definition: TVector3.cxx:402
SVector< double, 2 > v
Definition: Dict.h:5
XPoint xy[kMAXMK]
Definition: TGX11.cxx:122
Double_t Perp2() const
Definition: TVector3.h:312
TClass * IsA() const
unsigned int UInt_t
Definition: RtypesCore.h:42
TMarker * m
Definition: textangle.C:8
Double_t X() const
Definition: TVector3.h:220
Double_t Mag2() const
Definition: TVector3.h:298
ClassImp(TVector3) TVector3
Definition: TVector3.cxx:169
Double_t ACos(Double_t)
Definition: TMath.h:445
#define Printf
Definition: TGeoToOCC.h:18
Double_t Cos(Double_t)
Definition: TMath.h:424
Float_t phi
Definition: shapesAnim.C:6
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
Double_t Exp(Double_t x)
Definition: TMath.h:495
Double_t fX
Definition: TVector3.h:190
virtual const char * GetName() const
Returns name of object.
Definition: TObject.cxx:415
double Double_t
Definition: RtypesCore.h:55
void SetZ(Double_t)
Definition: TVector3.h:229
Double_t Eta() const
Definition: TVector3.h:359
TVector3 & Transform(const TRotation &)
transform this vector with a TRotation
Definition: TVector3.cxx:234
void RotateZ(Double_t)
rotate vector around Z
Definition: TVector3.cxx:331
Mother of all ROOT objects.
Definition: TObject.h:58
void SetPhi(Double_t)
Set phi keeping mag and theta constant (BaBar).
Definition: TVector3.cxx:414
Float_t px
Definition: hprod.C:33
Bool_t axis
Definition: geodemo.C:37
Double_t operator()(int) const
dereferencing operator const
Definition: TVector3.cxx:193
Double_t DeltaR(const TVector3 &) const
return deltaR with respect to v
Definition: TVector3.cxx:424
Double_t Sin(Double_t)
Definition: TMath.h:421
Double_t Y() const
Definition: TVector3.h:221
void SetX(Double_t)
Definition: TVector3.h:227
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
TVector3 operator+(const TVector3 &a, const TVector3 &b)
Definition: TVector3.cxx:467
void SetPtThetaPhi(Double_t pt, Double_t theta, Double_t phi)
set Pt, Theta and Phi
Definition: TVector3.cxx:392
float * q
Definition: THbookFile.cxx:87
virtual const char * GetTitle() const
Returns title of object.
Definition: TObject.cxx:459
void Print(Option_t *option="") const
This method must be overridden when a class wants to print itself.
Definition: TVector3.cxx:498
void SetMagThetaPhi(Double_t mag, Double_t theta, Double_t phi)
setter with mag, theta, phi
Definition: TVector3.cxx:434
Double_t Tan(Double_t)
Definition: TMath.h:427
virtual Version_t ReadVersion(UInt_t *start=0, UInt_t *bcnt=0, const TClass *cl=0)=0
void SetPtEtaPhi(Double_t pt, Double_t eta, Double_t phi)
set Pt, Eta and Phi
Definition: TVector3.cxx:384
Double_t ATan(Double_t)
Definition: TMath.h:451
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:904