Logo ROOT  
Reference Guide
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/** \class TVector2
13 \ingroup Physics
14
15TVector2 is a general two vector class, which can be used for
16the description of different vectors in 2D.
17*/
18
19#include "TROOT.h"
20#include "TVector2.h"
21#include "TMath.h"
22
24Double_t const kTWOPI = 2.*kPI;
25
27
28////////////////////////////////////////////////////////////////////////////////
29/// Constructor
30
32{
33 fX = 0.;
34 fY = 0.;
35}
36
37////////////////////////////////////////////////////////////////////////////////
38/// Constructor
39
41{
42 fX = v[0];
43 fY = v[1];
44}
45
46////////////////////////////////////////////////////////////////////////////////
47/// Constructor
48
50{
51 fX = x0;
52 fY = y0;
53}
54
55////////////////////////////////////////////////////////////////////////////////
56
58{
59}
60
61////////////////////////////////////////////////////////////////////////////////
62/// Return modulo of this vector
63
65{
66 return TMath::Sqrt(fX*fX+fY*fY);
67}
68
69////////////////////////////////////////////////////////////////////////////////
70/// Return module normalized to 1
71
73{
74 return (Mod2()) ? *this/Mod() : TVector2();
75}
76
77////////////////////////////////////////////////////////////////////////////////
78/// Return vector phi
79
81{
82 return TMath::Pi()+TMath::ATan2(-fY,-fX);
83}
84
85////////////////////////////////////////////////////////////////////////////////
86/// Returns phi angle in the interval [0,2*PI)
87
89 if(TMath::IsNaN(x)){
90 gROOT->Error("TVector2::Phi_0_2pi","function called with NaN");
91 return x;
92 }
93 while (x >= kTWOPI) x -= kTWOPI;
94 while (x < 0.) x += kTWOPI;
95 return x;
96}
97
98////////////////////////////////////////////////////////////////////////////////
99/// Returns phi angle in the interval [-PI,PI)
100
102 if(TMath::IsNaN(x)){
103 gROOT->Error("TVector2::Phi_mpi_pi","function called with NaN");
104 return x;
105 }
106 while (x >= kPI) x -= kTWOPI;
107 while (x < -kPI) x += kTWOPI;
108 return x;
109}
110
111////////////////////////////////////////////////////////////////////////////////
112/// Rotation by phi
113
115{
116 return TVector2( fX*TMath::Cos(phi)-fY*TMath::Sin(phi), fX*TMath::Sin(phi)+fY*TMath::Cos(phi) );
117}
118
119////////////////////////////////////////////////////////////////////////////////
120/// Set vector using mag and phi
121
123{
124 Double_t amag = TMath::Abs(mag);
125 fX = amag * TMath::Cos(phi);
126 fY = amag * TMath::Sin(phi);
127}
128////////////////////////////////////////////////////////////////////////////////
129/// Stream an object of class TVector2.
130
131void TVector2::Streamer(TBuffer &R__b)
132{
133 if (R__b.IsReading()) {
134 UInt_t R__s, R__c;
135 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
136 if (R__v > 2) {
137 R__b.ReadClassBuffer(TVector2::Class(), this, R__v, R__s, R__c);
138 return;
139 }
140 //====process old versions before automatic schema evolution
141 if (R__v < 2) TObject::Streamer(R__b);
142 R__b >> fX;
143 R__b >> fY;
144 R__b.CheckByteCount(R__s, R__c, TVector2::IsA());
145 //====end of old versions
146
147 } else {
149 }
150}
151
153{
154 //print vector parameters
155 Printf("%s %s (x,y)=(%f,%f) (rho,phi)=(%f,%f)",GetName(),GetTitle(),X(),Y(),
156 Mod(),Phi()*TMath::RadToDeg());
157}
void Class()
Definition: Class.C:29
short Version_t
Definition: RtypesCore.h:61
unsigned int UInt_t
Definition: RtypesCore.h:42
double Double_t
Definition: RtypesCore.h:55
const char Option_t
Definition: RtypesCore.h:62
#define ClassImp(name)
Definition: Rtypes.h:365
#define gROOT
Definition: TROOT.h:415
void Printf(const char *fmt,...)
Double_t const kPI
Definition: TVector2.cxx:23
Double_t const kTWOPI
Definition: TVector2.cxx:24
Buffer base class used for serializing objects.
Definition: TBuffer.h:42
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
virtual Version_t ReadVersion(UInt_t *start=0, UInt_t *bcnt=0, const TClass *cl=0)=0
virtual Int_t CheckByteCount(UInt_t startpos, UInt_t bcnt, const TClass *clss)=0
Bool_t IsReading() const
Definition: TBuffer.h:85
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
virtual const char * GetName() const
Returns name of object.
Definition: TObject.cxx:357
virtual const char * GetTitle() const
Returns title of object.
Definition: TObject.cxx:401
TVector2 is a general two vector class, which can be used for the description of different vectors in...
Definition: TVector2.h:18
Double_t Phi() const
Return vector phi.
Definition: TVector2.cxx:80
void SetMagPhi(Double_t mag, Double_t phi)
Set vector using mag and phi.
Definition: TVector2.cxx:122
virtual ~TVector2()
Definition: TVector2.cxx:57
Double_t Y() const
Definition: TVector2.h:73
static Double_t Phi_0_2pi(Double_t x)
Returns phi angle in the interval [0,2*PI)
Definition: TVector2.cxx:88
Double_t Mod2() const
Definition: TVector2.h:67
TVector2()
Constructor.
Definition: TVector2.cxx:31
Double_t fX
Definition: TVector2.h:24
Double_t X() const
Definition: TVector2.h:72
Double_t fY
Definition: TVector2.h:25
static Double_t Phi_mpi_pi(Double_t x)
Returns phi angle in the interval [-PI,PI)
Definition: TVector2.cxx:101
TVector2 Rotate(Double_t phi) const
Rotation by phi.
Definition: TVector2.cxx:114
void Print(Option_t *option="") const
This method must be overridden when a class wants to print itself.
Definition: TVector2.cxx:152
Double_t Mod() const
Return modulo of this vector.
Definition: TVector2.cxx:64
TVector2 Unit() const
Return module normalized to 1.
Definition: TVector2.cxx:72
Double_t x[n]
Definition: legend1.C:17
Bool_t IsNaN(Double_t x)
Definition: TMath.h:882
Double_t ATan2(Double_t y, Double_t x)
Definition: TMath.h:669
Double_t Sqrt(Double_t x)
Definition: TMath.h:681
Double_t Cos(Double_t)
Definition: TMath.h:631
constexpr Double_t Pi()
Definition: TMath.h:38
Double_t Sin(Double_t)
Definition: TMath.h:627
constexpr Double_t RadToDeg()
Conversion from radian to degree:
Definition: TMath.h:74
Short_t Abs(Short_t d)
Definition: TMathBase.h:120