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#include "TBuffer.h"
23
24
26Double_t const kTWOPI = 2.*kPI;
27
29
30////////////////////////////////////////////////////////////////////////////////
31/// Constructor
32
34{
35 fX = 0.;
36 fY = 0.;
37}
38
39////////////////////////////////////////////////////////////////////////////////
40/// Constructor
41
43{
44 fX = v[0];
45 fY = v[1];
46}
47
48////////////////////////////////////////////////////////////////////////////////
49/// Constructor
50
52{
53 fX = x0;
54 fY = y0;
55}
56
57////////////////////////////////////////////////////////////////////////////////
58
60{
61}
62
63////////////////////////////////////////////////////////////////////////////////
64/// Return modulo of this vector
65
67{
68 return TMath::Sqrt(fX*fX+fY*fY);
69}
70
71////////////////////////////////////////////////////////////////////////////////
72/// Return module normalized to 1
73
75{
76 return (Mod2()) ? *this/Mod() : TVector2();
77}
78
79////////////////////////////////////////////////////////////////////////////////
80/// Return vector phi
81
83{
84 return TMath::Pi()+TMath::ATan2(-fY,-fX);
85}
86
87////////////////////////////////////////////////////////////////////////////////
88/// Returns phi angle in the interval [0,2*PI)
89
91 if(TMath::IsNaN(x)){
92 gROOT->Error("TVector2::Phi_0_2pi","function called with NaN");
93 return x;
94 }
95 while (x >= kTWOPI) x -= kTWOPI;
96 while (x < 0.) x += kTWOPI;
97 return x;
98}
99
100////////////////////////////////////////////////////////////////////////////////
101/// Returns phi angle in the interval [-PI,PI)
102
104 if(TMath::IsNaN(x)){
105 gROOT->Error("TVector2::Phi_mpi_pi","function called with NaN");
106 return x;
107 }
108 while (x >= kPI) x -= kTWOPI;
109 while (x < -kPI) x += kTWOPI;
110 return x;
111}
112
113////////////////////////////////////////////////////////////////////////////////
114/// Rotation by phi
115
117{
118 return TVector2( fX*TMath::Cos(phi)-fY*TMath::Sin(phi), fX*TMath::Sin(phi)+fY*TMath::Cos(phi) );
119}
120
121////////////////////////////////////////////////////////////////////////////////
122/// Set vector using mag and phi
123
125{
126 Double_t amag = TMath::Abs(mag);
127 fX = amag * TMath::Cos(phi);
128 fY = amag * TMath::Sin(phi);
129}
130////////////////////////////////////////////////////////////////////////////////
131/// Stream an object of class TVector2.
132
133void TVector2::Streamer(TBuffer &R__b)
134{
135 if (R__b.IsReading()) {
136 UInt_t R__s, R__c;
137 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
138 if (R__v > 2) {
139 R__b.ReadClassBuffer(TVector2::Class(), this, R__v, R__s, R__c);
140 return;
141 }
142 //====process old versions before automatic schema evolution
143 if (R__v < 2) TObject::Streamer(R__b);
144 R__b >> fX;
145 R__b >> fY;
146 R__b.CheckByteCount(R__s, R__c, TVector2::IsA());
147 //====end of old versions
148
149 } else {
151 }
152}
153
155{
156 //print vector parameters
157 Printf("%s %s (x,y)=(%f,%f) (rho,phi)=(%f,%f)",GetName(),GetTitle(),X(),Y(),
158 Mod(),Phi()*TMath::RadToDeg());
159}
void Class()
Definition: Class.C:29
short Version_t
Definition: RtypesCore.h:63
unsigned int UInt_t
Definition: RtypesCore.h:44
double Double_t
Definition: RtypesCore.h:57
const char Option_t
Definition: RtypesCore.h:64
#define ClassImp(name)
Definition: Rtypes.h:361
#define gROOT
Definition: TROOT.h:406
void Printf(const char *fmt,...)
Double_t const kPI
Definition: TVector2.cxx:25
Double_t const kTWOPI
Definition: TVector2.cxx:26
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:82
void SetMagPhi(Double_t mag, Double_t phi)
Set vector using mag and phi.
Definition: TVector2.cxx:124
virtual ~TVector2()
Definition: TVector2.cxx:59
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:90
Double_t Mod2() const
Definition: TVector2.h:67
TVector2()
Constructor.
Definition: TVector2.cxx:33
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:103
TVector2 Rotate(Double_t phi) const
Rotation by phi.
Definition: TVector2.cxx:116
void Print(Option_t *option="") const
This method must be overridden when a class wants to print itself.
Definition: TVector2.cxx:154
Double_t Mod() const
Return modulo of this vector.
Definition: TVector2.cxx:66
TVector2 Unit() const
Return module normalized to 1.
Definition: TVector2.cxx:74
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