ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TVector2.cxx
Go to the documentation of this file.
1 // @(#)root/physics:$Id$
2 // Author: Pasha Murat 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 #include "TROOT.h"
13 #include "TVector2.h"
14 #include "TClass.h"
15 #include "TMath.h"
16 
18 Double_t const kTWOPI = 2.*kPI;
19 
21 
22 ////////////////////////////////////////////////////////////////////////////////
23 ///constructor
24 
26 {
27  fX = 0.;
28  fY = 0.;
29 }
30 
31 ////////////////////////////////////////////////////////////////////////////////
32 ///constructor
33 
35 {
36  fX = v[0];
37  fY = v[1];
38 }
39 
40 ////////////////////////////////////////////////////////////////////////////////
41 ///constructor
42 
44 {
45  fX = x0;
46  fY = y0;
47 }
48 
49 ////////////////////////////////////////////////////////////////////////////////
50 
52 {
53 }
54 
55 ////////////////////////////////////////////////////////////////////////////////
56 /// return modulo of this vector
57 
59 {
60  return TMath::Sqrt(fX*fX+fY*fY);
61 }
62 
63 ////////////////////////////////////////////////////////////////////////////////
64 /// return module normalized to 1
65 
67 {
68  return (Mod2()) ? *this/Mod() : TVector2();
69 }
70 
71 ////////////////////////////////////////////////////////////////////////////////
72 /// return vector phi
73 
75 {
76  return TMath::Pi()+TMath::ATan2(-fY,-fX);
77 }
78 
79 ////////////////////////////////////////////////////////////////////////////////
80 /// (static function) returns phi angle in the interval [0,2*PI)
81 
83  if(TMath::IsNaN(x)){
84  gROOT->Error("TVector2::Phi_0_2pi","function called with NaN");
85  return x;
86  }
87  while (x >= kTWOPI) x -= kTWOPI;
88  while (x < 0.) x += kTWOPI;
89  return x;
90 }
91 
92 ////////////////////////////////////////////////////////////////////////////////
93 /// (static function) returns phi angle in the interval [-PI,PI)
94 
96  if(TMath::IsNaN(x)){
97  gROOT->Error("TVector2::Phi_mpi_pi","function called with NaN");
98  return x;
99  }
100  while (x >= kPI) x -= kTWOPI;
101  while (x < -kPI) x += kTWOPI;
102  return x;
103 }
104 
105 ////////////////////////////////////////////////////////////////////////////////
106 ///rotation by phi
107 
109 {
110  return TVector2( fX*TMath::Cos(phi)-fY*TMath::Sin(phi), fX*TMath::Sin(phi)+fY*TMath::Cos(phi) );
111 }
112 
113 ////////////////////////////////////////////////////////////////////////////////
114 ///set vector using mag and phi
115 
117 {
118  Double_t amag = TMath::Abs(mag);
119  fX = amag * TMath::Cos(phi);
120  fY = amag * TMath::Sin(phi);
121 }
122 ////////////////////////////////////////////////////////////////////////////////
123 /// Stream an object of class TVector2.
124 
125 void TVector2::Streamer(TBuffer &R__b)
126 {
127  if (R__b.IsReading()) {
128  UInt_t R__s, R__c;
129  Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
130  if (R__v > 2) {
131  R__b.ReadClassBuffer(TVector2::Class(), this, R__v, R__s, R__c);
132  return;
133  }
134  //====process old versions before automatic schema evolution
135  if (R__v < 2) TObject::Streamer(R__b);
136  R__b >> fX;
137  R__b >> fY;
138  R__b.CheckByteCount(R__s, R__c, TVector2::IsA());
139  //====end of old versions
140 
141  } else {
142  R__b.WriteClassBuffer(TVector2::Class(),this);
143  }
144 }
145 
147 {
148  //print vector parameters
149  Printf("%s %s (x,y)=(%f,%f) (rho,phi)=(%f,%f)",GetName(),GetTitle(),X(),Y(),
150  Mod(),Phi()*TMath::RadToDeg());
151 }
152 
static Double_t Phi_0_2pi(Double_t x)
(static function) returns phi angle in the interval [0,2*PI)
Definition: TVector2.cxx:82
TVector2 Unit() const
return module normalized to 1
Definition: TVector2.cxx:66
Double_t fY
Definition: TVector2.h:25
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
Double_t Mod() const
return modulo of this vector
Definition: TVector2.cxx:58
Bool_t IsReading() const
Definition: TBuffer.h:83
short Version_t
Definition: RtypesCore.h:61
const char Option_t
Definition: RtypesCore.h:62
Double_t const kTWOPI
Definition: TVector2.cxx:18
Double_t fX
Definition: TVector2.h:24
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
#define gROOT
Definition: TROOT.h:344
Double_t RadToDeg()
Definition: TMath.h:49
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
Double_t x[n]
Definition: legend1.C:17
void Class()
Definition: Class.C:29
static Double_t Phi_mpi_pi(Double_t x)
(static function) returns phi angle in the interval [-PI,PI)
Definition: TVector2.cxx:95
virtual ~TVector2()
Definition: TVector2.cxx:51
Double_t ATan2(Double_t, Double_t)
Definition: TMath.h:454
SVector< double, 2 > v
Definition: Dict.h:5
TClass * IsA() const
unsigned int UInt_t
Definition: RtypesCore.h:42
#define Printf
Definition: TGeoToOCC.h:18
ClassImp(TVector2) TVector2
constructor
Definition: TVector2.cxx:20
Double_t Cos(Double_t)
Definition: TMath.h:424
Double_t Pi()
Definition: TMath.h:44
Float_t phi
Definition: shapesAnim.C:6
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
Double_t X() const
Definition: TVector2.h:72
virtual const char * GetName() const
Returns name of object.
Definition: TObject.cxx:415
double Double_t
Definition: RtypesCore.h:55
Double_t Phi() const
return vector phi
Definition: TVector2.cxx:74
Double_t const kPI
Definition: TVector2.cxx:17
TVector2 Rotate(Double_t phi) const
rotation by phi
Definition: TVector2.cxx:108
Int_t IsNaN(Double_t x)
Definition: TMath.h:617
Double_t Sin(Double_t)
Definition: TMath.h:421
Double_t Y() const
Definition: TVector2.h:73
void SetMagPhi(Double_t mag, Double_t phi)
set vector using mag and phi
Definition: TVector2.cxx:116
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
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: TVector2.cxx:146
Double_t Mod2() const
Definition: TVector2.h:67
virtual Version_t ReadVersion(UInt_t *start=0, UInt_t *bcnt=0, const TClass *cl=0)=0