CLHEP 2.4.6.4
C++ Class Library for High Energy Physics
Loading...
Searching...
No Matches
CLHEP::Hep3Vector Class Reference

#include <ThreeVector.h>

Public Types

enum  {
  X =0 , Y =1 , Z =2 , NUM_COORDINATES =3 ,
  SIZE =NUM_COORDINATES
}
 

Public Member Functions

 Hep3Vector ()
 
 Hep3Vector (double x)
 
 Hep3Vector (double x, double y)
 
 Hep3Vector (double x, double y, double z)
 
 Hep3Vector (const Hep3Vector &)
 
 Hep3Vector (Hep3Vector &&)=default
 
 ~Hep3Vector ()
 
double operator() (int) const
 
double operator[] (int) const
 
doubleoperator() (int)
 
doubleoperator[] (int)
 
double x () const
 
double y () const
 
double z () const
 
void setX (double)
 
void setY (double)
 
void setZ (double)
 
void set (double x, double y, double z)
 
double phi () const
 
double theta () const
 
double cosTheta () const
 
double cos2Theta () const
 
double mag2 () const
 
double mag () const
 
void setPhi (double)
 
void setTheta (double)
 
void setMag (double)
 
double perp2 () const
 
double perp () const
 
void setPerp (double)
 
void setCylTheta (double)
 
double perp2 (const Hep3Vector &) const
 
double perp (const Hep3Vector &) const
 
Hep3Vectoroperator= (const Hep3Vector &)
 
Hep3Vectoroperator= (Hep3Vector &&)=default
 
bool operator== (const Hep3Vector &) const
 
bool operator!= (const Hep3Vector &) const
 
bool isNear (const Hep3Vector &, double epsilon=tolerance) const
 
double howNear (const Hep3Vector &v) const
 
double deltaR (const Hep3Vector &v) const
 
Hep3Vectoroperator+= (const Hep3Vector &)
 
Hep3Vectoroperator-= (const Hep3Vector &)
 
Hep3Vector operator- () const
 
Hep3Vectoroperator*= (double)
 
Hep3Vectoroperator/= (double)
 
Hep3Vector unit () const
 
Hep3Vector orthogonal () const
 
double dot (const Hep3Vector &) const
 
Hep3Vector cross (const Hep3Vector &) const
 
double angle (const Hep3Vector &) const
 
double pseudoRapidity () const
 
void setEta (double p)
 
void setCylEta (double p)
 
Hep3VectorrotateX (double)
 
Hep3VectorrotateY (double)
 
Hep3VectorrotateZ (double)
 
Hep3VectorrotateUz (const Hep3Vector &)
 
Hep3Vectorrotate (double, const Hep3Vector &)
 
Hep3Vectoroperator*= (const HepRotation &)
 
Hep3Vectortransform (const HepRotation &)
 
void setRThetaPhi (double r, double theta, double phi)
 
void setREtaPhi (double r, double eta, double phi)
 
void setRhoPhiZ (double rho, double phi, double z)
 
void setRhoPhiTheta (double rho, double phi, double theta)
 
void setRhoPhiEta (double rho, double phi, double eta)
 
double getX () const
 
double getY () const
 
double getZ () const
 
double getR () const
 
double getTheta () const
 
double getPhi () const
 
double r () const
 
double rho () const
 
double getRho () const
 
double eta () const
 
double getEta () const
 
void setR (double s)
 
void setRho (double s)
 
int compare (const Hep3Vector &v) const
 
bool operator> (const Hep3Vector &v) const
 
bool operator< (const Hep3Vector &v) const
 
bool operator>= (const Hep3Vector &v) const
 
bool operator<= (const Hep3Vector &v) const
 
double diff2 (const Hep3Vector &v) const
 
bool isParallel (const Hep3Vector &v, double epsilon=tolerance) const
 
bool isOrthogonal (const Hep3Vector &v, double epsilon=tolerance) const
 
double howParallel (const Hep3Vector &v) const
 
double howOrthogonal (const Hep3Vector &v) const
 
double beta () const
 
double gamma () const
 
double coLinearRapidity () const
 
double angle () const
 
double theta (const Hep3Vector &v2) const
 
double cosTheta (const Hep3Vector &v2) const
 
double cos2Theta (const Hep3Vector &v2) const
 
Hep3Vector project () const
 
Hep3Vector project (const Hep3Vector &v2) const
 
Hep3Vector perpPart () const
 
Hep3Vector perpPart (const Hep3Vector &v2) const
 
double rapidity () const
 
double rapidity (const Hep3Vector &v2) const
 
double eta (const Hep3Vector &v2) const
 
double polarAngle (const Hep3Vector &v2) const
 
double deltaPhi (const Hep3Vector &v2) const
 
double azimAngle (const Hep3Vector &v2) const
 
double polarAngle (const Hep3Vector &v2, const Hep3Vector &ref) const
 
double azimAngle (const Hep3Vector &v2, const Hep3Vector &ref) const
 
Hep3Vectorrotate (const Hep3Vector &axis, double delta)
 
Hep3Vectorrotate (const HepAxisAngle &ax)
 
Hep3Vectorrotate (const HepEulerAngles &e)
 
Hep3Vectorrotate (double phi, double theta, double psi)
 

Static Public Member Functions

static double setTolerance (double tol)
 
static double getTolerance ()
 

Static Public Attributes

static const int ToleranceTicks = 100
 

Protected Member Functions

void setSpherical (double r, double theta, double phi)
 
void setCylindrical (double r, double phi, double z)
 
double negativeInfinity () const
 

Protected Attributes

double data [3]
 

Static Protected Attributes

static double tolerance = Hep3Vector::ToleranceTicks * 2.22045e-16
 

Detailed Description

Author

Definition at line 36 of file ThreeVector.h.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
Enumerator
NUM_COORDINATES 
SIZE 

Definition at line 42 of file ThreeVector.h.

Constructor & Destructor Documentation

◆ Hep3Vector() [1/6]

CLHEP::Hep3Vector::Hep3Vector ( )

◆ Hep3Vector() [2/6]

CLHEP::Hep3Vector::Hep3Vector ( double  x)
explicit

◆ Hep3Vector() [3/6]

CLHEP::Hep3Vector::Hep3Vector ( double  x,
double  y 
)

◆ Hep3Vector() [4/6]

CLHEP::Hep3Vector::Hep3Vector ( double  x,
double  y,
double  z 
)

◆ Hep3Vector() [5/6]

CLHEP::Hep3Vector::Hep3Vector ( const Hep3Vector )
inline

◆ Hep3Vector() [6/6]

CLHEP::Hep3Vector::Hep3Vector ( Hep3Vector &&  )
inlinedefault

◆ ~Hep3Vector()

CLHEP::Hep3Vector::~Hep3Vector ( )
inline

Member Function Documentation

◆ angle() [1/2]

double CLHEP::Hep3Vector::angle ( ) const
inline

Referenced by polarAngle().

◆ angle() [2/2]

double CLHEP::Hep3Vector::angle ( const Hep3Vector ) const

Referenced by azimAngle(), main(), and polarAngle().

◆ azimAngle() [1/2]

double CLHEP::Hep3Vector::azimAngle ( const Hep3Vector v2) const

◆ azimAngle() [2/2]

double CLHEP::Hep3Vector::azimAngle ( const Hep3Vector v2,
const Hep3Vector ref 
) const

Definition at line 38 of file SpaceVectorD.cc.

39 {
40
41 Hep3Vector vperp ( perpPart(ref) );
42 if ( vperp.mag2() == 0 ) {
43 ZMthrowC (ZMxpvAmbiguousAngle(
44 "Cannot find azimuthal angle with reference direction parallel to "
45 "vector 1 -- will return zero"));
46 return 0;
47 }
48
49 Hep3Vector v2perp ( v2.perpPart(ref) );
50 if ( v2perp.mag2() == 0 ) {
51 ZMthrowC (ZMxpvAmbiguousAngle(
52 "Cannot find azimuthal angle with reference direction parallel to "
53 "vector 2 -- will return zero"));
54 return 0;
55 }
56
57 double ang = vperp.angle(v2perp);
58
59 // Now compute the sign of the answer: that of U*(VxV2) or
60 // the equivalent expression V*(V2xU).
61
62 if ( dot(v2.cross(ref)) >= 0 ) {
63 return ang;
64 } else {
65 return -ang;
66 }
67
68 //-| Note that if V*(V2xU) is zero, we want to return 0 or PI
69 //-| depending on whether vperp is aligned or antialigned with v2perp.
70 //-| The computed angle() expression does this properly.
71
72} /* azimAngle (v2, ref) */
#define ZMthrowC(A)
Definition: ZMxpv.h:133
double dot(const Hep3Vector &) const
Hep3Vector perpPart() const

◆ beta()

double CLHEP::Hep3Vector::beta ( ) const

Definition at line 28 of file SpaceVectorP.cc.

28 {
29 double b = std::sqrt(mag2());
30 if (b >= 1) {
31 ZMthrowA (ZMxpvTachyonic(
32 "Beta taken for Hep3Vector of at least unit length"));
33 }
34 return b;
35}
#define ZMthrowA(A)
Definition: ZMxpv.h:128
double mag2() const

Referenced by coLinearRapidity(), CLHEP::HepRotation::distance2(), CLHEP::HepRotationX::distance2(), CLHEP::HepRotationY::distance2(), and CLHEP::HepRotationZ::distance2().

◆ coLinearRapidity()

double CLHEP::Hep3Vector::coLinearRapidity ( ) const

Definition at line 66 of file SpaceVectorP.cc.

66 {
67 double b = beta();
68 if (b == 1) {
69 ZMthrowA (ZMxpvTachyonic(
70 "Co-linear Rapidity taken for Hep3Vector of unit length -- "
71 "the log should return infinity"));
72 }
73 if (b > 1) {
74 ZMthrowA (ZMxpvTachyonic(
75 "Co-linear Rapidity taken for Hep3Vector of more than unit length -- "
76 "the log would return a NAN" ));
77 }
78 // Want inverse std::tanh(b):
79 return (.5 * std::log((1+b)/(1-b)) );
80}
double beta() const
Definition: SpaceVectorP.cc:28

◆ compare()

int CLHEP::Hep3Vector::compare ( const Hep3Vector v) const

Definition at line 121 of file SpaceVector.cc.

121 {
122 if ( z() > v.z() ) {
123 return 1;
124 } else if ( z() < v.z() ) {
125 return -1;
126 } else if ( y() > v.y() ) {
127 return 1;
128 } else if ( y() < v.y() ) {
129 return -1;
130 } else if ( x() > v.x() ) {
131 return 1;
132 } else if ( x() < v.x() ) {
133 return -1;
134 } else {
135 return 0;
136 }
137} /* Compare */
double z() const
double x() const
double y() const

Referenced by CLHEP::HepLorentzVector::compare(), operator<(), operator<=(), operator>(), and operator>=().

◆ cos2Theta() [1/2]

double CLHEP::Hep3Vector::cos2Theta ( ) const
inline

◆ cos2Theta() [2/2]

double CLHEP::Hep3Vector::cos2Theta ( const Hep3Vector v2) const

Definition at line 167 of file ThreeVector.cc.

167 {
168 double arg;
169 double ptot2 = mag2();
170 double qtot2 = q.mag2();
171 if ( ptot2 == 0 || qtot2 == 0 ) {
172 arg = 1.0;
173 }else{
174 double pdq = dot(q);
175 arg = (pdq/ptot2) * (pdq/qtot2);
176 // More naive methods overflow on vectors which can be squared
177 // but can't be raised to the 4th power.
178 if(arg > 1.0) arg = 1.0;
179 }
180 return arg;
181}

◆ cosTheta() [1/2]

double CLHEP::Hep3Vector::cosTheta ( ) const
inline

Referenced by main().

◆ cosTheta() [2/2]

double CLHEP::Hep3Vector::cosTheta ( const Hep3Vector v2) const

Definition at line 154 of file ThreeVector.cc.

154 {
155 double arg;
156 double ptot2 = mag2()*q.mag2();
157 if(ptot2 <= 0) {
158 arg = 0.0;
159 }else{
160 arg = dot(q)/std::sqrt(ptot2);
161 if(arg > 1.0) arg = 1.0;
162 if(arg < -1.0) arg = -1.0;
163 }
164 return arg;
165}

◆ cross()

Hep3Vector CLHEP::Hep3Vector::cross ( const Hep3Vector ) const
inline

◆ deltaPhi()

double CLHEP::Hep3Vector::deltaPhi ( const Hep3Vector v2) const

Definition at line 138 of file ThreeVector.cc.

138 {
139 double dphi = v2.getPhi() - getPhi();
140 if ( dphi > CLHEP::pi ) {
141 dphi -= CLHEP::twopi;
142 } else if ( dphi <= -CLHEP::pi ) {
143 dphi += CLHEP::twopi;
144 }
145 return dphi;
146} /* deltaPhi */
double getPhi() const

Referenced by deltaR(), and CLHEP::HepLorentzVector::deltaR().

◆ deltaR()

double CLHEP::Hep3Vector::deltaR ( const Hep3Vector v) const

Definition at line 148 of file ThreeVector.cc.

148 {
149 double a = eta() - v.eta();
150 double b = deltaPhi(v);
151 return std::sqrt ( a*a + b*b );
152} /* deltaR */
double eta() const
double deltaPhi(const Hep3Vector &v2) const
Definition: ThreeVector.cc:138

◆ diff2()

double CLHEP::Hep3Vector::diff2 ( const Hep3Vector v) const
inline

◆ dot()

◆ eta() [1/2]

double CLHEP::Hep3Vector::eta ( ) const

Referenced by deltaR().

◆ eta() [2/2]

double CLHEP::Hep3Vector::eta ( const Hep3Vector v2) const

Definition at line 113 of file SpaceVectorP.cc.

113 {
114 // Defined as -std::log ( std::tan ( .5* theta(u) ) );
115 //
116 // Quicker is to use cosTheta:
117 // std::tan (theta/2) = std::sin(theta)/(1 + std::cos(theta))
118
119 double r1 = getR();
120 double v2r = v2.mag();
121 if ( (r1 == 0) || (v2r == 0) ) {
122 ZMthrowA (ZMxpvAmbiguousAngle(
123 "Cannot find pseudorapidity of a zero vector relative to a vector"));
124 return 0.;
125 }
126 double c = dot(v2)/(r1*v2r);
127 if ( c >= 1 ) {
128 c = 1; //-| We don't want to return NAN because of roundoff
129 ZMthrowC (ZMxpvInfinity(
130 "Pseudorapidity of vector relative to parallel vector -- "
131 "will give infinite result"));
132 // We can just go on; tangent will be 0, so
133 // std::log (tangent) will be -INFINITY, so result
134 // will be +INFINITY.
135 }
136 if ( c <= -1 ) {
137 ZMthrowC (ZMxpvInfinity(
138 "Pseudorapidity of vector relative to anti-parallel vector -- "
139 "will give negative infinite result"));
140 //-| We don't want to return NAN because of roundoff
141 return ( negativeInfinity() );
142 // If we just went on, the tangent would be NAN
143 // so return would be NAN. But the proper limit
144 // of tan is +Infinity, so the return should be
145 // -INFINITY.
146 }
147
148 double tangent = std::sqrt (1-c*c) / ( 1 + c );
149 return (- std::log (tangent));
150
151} /* eta (u) */
double getR() const
double negativeInfinity() const
Definition: SpaceVector.cc:284

◆ gamma()

double CLHEP::Hep3Vector::gamma ( ) const

Definition at line 37 of file SpaceVectorP.cc.

37 {
38 double bbeta = std::sqrt(mag2());
39 if (bbeta == 1) {
40 ZMthrowA (ZMxpvTachyonic(
41 "Gamma taken for Hep3Vector of unit magnitude -- infinite result"));
42 }
43 if (bbeta > 1) {
44 ZMthrowA (ZMxpvTachyonic(
45 "Gamma taken for Hep3Vector of more than unit magnitude -- "
46 "the sqrt function would return NAN" ));
47 }
48 return 1/std::sqrt(1-bbeta*bbeta);
49}

◆ getEta()

double CLHEP::Hep3Vector::getEta ( ) const

◆ getPhi()

double CLHEP::Hep3Vector::getPhi ( ) const
inline

◆ getR()

double CLHEP::Hep3Vector::getR ( ) const
inline

Referenced by eta(), and setEta().

◆ getRho()

double CLHEP::Hep3Vector::getRho ( ) const
inline

Referenced by setCylEta(), and setCylTheta().

◆ getTheta()

double CLHEP::Hep3Vector::getTheta ( ) const
inline

Referenced by polarAngle().

◆ getTolerance()

static double CLHEP::Hep3Vector::getTolerance ( )
inlinestatic

◆ getX()

double CLHEP::Hep3Vector::getX ( ) const
inline

◆ getY()

double CLHEP::Hep3Vector::getY ( ) const
inline

◆ getZ()

◆ howNear()

double CLHEP::Hep3Vector::howNear ( const Hep3Vector v) const

Definition at line 125 of file ThreeVector.cc.

125 {
126 // | V1 - V2 | **2 / V1 dot V2, up to 1
127 double d = (*this - v).mag2();
128 double vdv = dot(v);
129 if ( (vdv > 0) && (d < vdv) ) {
130 return std::sqrt (d/vdv);
131 } else if ( (vdv == 0) && (d == 0) ) {
132 return 0;
133 } else {
134 return 1;
135 }
136} /* howNear */

◆ howOrthogonal()

double CLHEP::Hep3Vector::howOrthogonal ( const Hep3Vector v) const

Definition at line 219 of file SpaceVector.cc.

219 {
220 // | V1 dot V2 | / | V1 x V2 |
221
222 double v1v2 = std::fabs(dot(v));
223 //-| Safe because both v1 and v2 can be squared
224 if ( v1v2 == 0 ) {
225 return 0; // Even if one or both are 0, they are considered orthogonal
226 }
227 Hep3Vector v1Xv2 ( cross(v) );
228 double abscross = v1Xv2.mag();
229 if ( v1v2 >= abscross ) {
230 return 1;
231 } else {
232 return v1v2/abscross;
233 }
234
235} /* howOrthogonal() */
Hep3Vector cross(const Hep3Vector &) const

◆ howParallel()

double CLHEP::Hep3Vector::howParallel ( const Hep3Vector v) const

Definition at line 168 of file SpaceVector.cc.

168 {
169 // | V1 x V2 | / | V1 dot V2 |
170 double v1v2 = std::fabs(dot(v));
171 if ( v1v2 == 0 ) {
172 // Zero is parallel to no other vector except for zero.
173 return ( (mag2() == 0) && (v.mag2() == 0) ) ? 0 : 1;
174 }
175 Hep3Vector v1Xv2 ( cross(v) );
176 double abscross = v1Xv2.mag();
177 if ( abscross >= v1v2 ) {
178 return 1;
179 } else {
180 return abscross/v1v2;
181 }
182} /* howParallel() */

◆ isNear()

bool CLHEP::Hep3Vector::isNear ( const Hep3Vector v,
double  epsilon = tolerance 
) const

Definition at line 120 of file ThreeVector.cc.

120 {
121 double limit = dot(v)*epsilon*epsilon;
122 return ( (*this - v).mag2() <= limit );
123} /* isNear() */

◆ isOrthogonal()

bool CLHEP::Hep3Vector::isOrthogonal ( const Hep3Vector v,
double  epsilon = tolerance 
) const

Definition at line 237 of file SpaceVector.cc.

238 {
239// | V1 x V2 | **2 <= epsilon **2 | V1 dot V2 | **2
240// V1 is *this, V2 is v
241
242 static const double TOOBIG = std::pow(2.0,507);
243 static const double SCALE = std::pow(2.0,-507);
244 double v1v2 = std::fabs(dot(v));
245 //-| Safe because both v1 and v2 can be squared
246 if ( v1v2 >= TOOBIG ) {
247 Hep3Vector sv1 ( *this * SCALE );
248 Hep3Vector sv2 ( v * SCALE );
249 Hep3Vector sv1Xsv2 = sv1.cross(sv2);
250 double x2 = sv1Xsv2.mag2();
251 double limit = epsilon*epsilon*x2;
252 double y2 = v1v2*SCALE*SCALE;
253 return ( y2*y2 <= limit );
254 }
255
256 // At this point we know v1v2 can be squared.
257
258 Hep3Vector eps_v1Xv2 ( cross(epsilon*v) );
259 if ( (std::fabs (eps_v1Xv2.x()) > TOOBIG) ||
260 (std::fabs (eps_v1Xv2.y()) > TOOBIG) ||
261 (std::fabs (eps_v1Xv2.z()) > TOOBIG) ) {
262 return true;
263 }
264
265 // At this point we know all the math we need can be done.
266
267 return ( v1v2*v1v2 <= eps_v1Xv2.mag2() );
268
269} /* isOrthogonal() */

◆ isParallel()

bool CLHEP::Hep3Vector::isParallel ( const Hep3Vector v,
double  epsilon = tolerance 
) const

Definition at line 184 of file SpaceVector.cc.

185 {
186 // | V1 x V2 | **2 <= epsilon **2 | V1 dot V2 | **2
187 // V1 is *this, V2 is v
188
189 static const double TOOBIG = std::pow(2.0,507);
190 static const double SCALE = std::pow(2.0,-507);
191 double v1v2 = std::fabs(dot(v));
192 if ( v1v2 == 0 ) {
193 return ( (mag2() == 0) && (v.mag2() == 0) );
194 }
195 if ( v1v2 >= TOOBIG ) {
196 Hep3Vector sv1 ( *this * SCALE );
197 Hep3Vector sv2 ( v * SCALE );
198 Hep3Vector sv1Xsv2 = sv1.cross(sv2);
199 double x2 = sv1Xsv2.mag2();
200 double limit = v1v2*SCALE*SCALE;
201 limit = epsilon*epsilon*limit*limit;
202 return ( x2 <= limit );
203 }
204
205 // At this point we know v1v2 can be squared.
206
207 Hep3Vector v1Xv2 ( cross(v) );
208 if ( (std::fabs (v1Xv2.x()) > TOOBIG) ||
209 (std::fabs (v1Xv2.y()) > TOOBIG) ||
210 (std::fabs (v1Xv2.z()) > TOOBIG) ) {
211 return false;
212 }
213
214 return ( (v1Xv2.mag2()) <= ((epsilon * v1v2) * (epsilon * v1v2)) );
215
216} /* isParallel() */

◆ mag()

◆ mag2()

◆ negativeInfinity()

double CLHEP::Hep3Vector::negativeInfinity ( ) const
protected

Definition at line 284 of file SpaceVector.cc.

284 {
285 // A byte-order-independent way to return -Infinity
286 struct Dib {
287 union {
288 double d;
289 unsigned char i[8];
290 } u;
291 };
292 Dib negOne;
293 Dib posTwo;
294 negOne.u.d = -1.0;
295 posTwo.u.d = 2.0;
296 Dib value;
297 int k;
298 for (k=0; k<8; k++) {
299 value.u.i[k] = negOne.u.i[k] | posTwo.u.i[k];
300 }
301 return value.u.d;
302}

Referenced by eta().

◆ operator!=()

bool CLHEP::Hep3Vector::operator!= ( const Hep3Vector ) const
inline

◆ operator()() [1/2]

double & CLHEP::Hep3Vector::operator() ( int  )
inline

◆ operator()() [2/2]

double CLHEP::Hep3Vector::operator() ( int  ) const
inline

◆ operator*=() [1/2]

Hep3Vector & CLHEP::Hep3Vector::operator*= ( const HepRotation m1)

Definition at line 17 of file ThreeVectorR.cc.

17 {
18 return *this = m1 * (*this);
19}

◆ operator*=() [2/2]

Hep3Vector & CLHEP::Hep3Vector::operator*= ( double  )
inline

Referenced by rotate().

◆ operator+=()

Hep3Vector & CLHEP::Hep3Vector::operator+= ( const Hep3Vector )
inline

◆ operator-()

Hep3Vector CLHEP::Hep3Vector::operator- ( ) const
inline

◆ operator-=()

Hep3Vector & CLHEP::Hep3Vector::operator-= ( const Hep3Vector )
inline

◆ operator/=()

Hep3Vector & CLHEP::Hep3Vector::operator/= ( double  c)

Definition at line 304 of file ThreeVector.cc.

304 {
305 if (c == 0) {
306 ZMthrowA (ZMxpvInfiniteVector(
307 "Attempt to do vector /= 0 -- "
308 "division by zero would produce infinite or NAN components"));
309 }
310 *this *= 1.0/c;
311 return *this;
312}

◆ operator<()

bool CLHEP::Hep3Vector::operator< ( const Hep3Vector v) const

Definition at line 143 of file SpaceVector.cc.

143 {
144 return (compare(v) < 0);
145}
int compare(const Hep3Vector &v) const
Definition: SpaceVector.cc:121

◆ operator<=()

bool CLHEP::Hep3Vector::operator<= ( const Hep3Vector v) const

Definition at line 149 of file SpaceVector.cc.

149 {
150 return (compare(v) <= 0);
151}

◆ operator=() [1/2]

Hep3Vector & CLHEP::Hep3Vector::operator= ( const Hep3Vector )
inline

◆ operator=() [2/2]

Hep3Vector & CLHEP::Hep3Vector::operator= ( Hep3Vector &&  )
inlinedefault

◆ operator==()

bool CLHEP::Hep3Vector::operator== ( const Hep3Vector ) const
inline

◆ operator>()

bool CLHEP::Hep3Vector::operator> ( const Hep3Vector v) const

Definition at line 140 of file SpaceVector.cc.

140 {
141 return (compare(v) > 0);
142}

◆ operator>=()

bool CLHEP::Hep3Vector::operator>= ( const Hep3Vector v) const

Definition at line 146 of file SpaceVector.cc.

146 {
147 return (compare(v) >= 0);
148}

◆ operator[]() [1/2]

double & CLHEP::Hep3Vector::operator[] ( int  )
inline

◆ operator[]() [2/2]

double CLHEP::Hep3Vector::operator[] ( int  ) const
inline

◆ orthogonal()

Hep3Vector CLHEP::Hep3Vector::orthogonal ( ) const
inline

◆ perp() [1/2]

double CLHEP::Hep3Vector::perp ( ) const
inline

Referenced by main().

◆ perp() [2/2]

double CLHEP::Hep3Vector::perp ( const Hep3Vector ) const
inline

◆ perp2() [1/2]

double CLHEP::Hep3Vector::perp2 ( ) const
inline

Referenced by main().

◆ perp2() [2/2]

double CLHEP::Hep3Vector::perp2 ( const Hep3Vector ) const
inline

◆ perpPart() [1/2]

Hep3Vector CLHEP::Hep3Vector::perpPart ( ) const
inline

Referenced by azimAngle().

◆ perpPart() [2/2]

Hep3Vector CLHEP::Hep3Vector::perpPart ( const Hep3Vector v2) const
inline

◆ phi()

double CLHEP::Hep3Vector::phi ( ) const
inline

Referenced by main().

◆ polarAngle() [1/2]

double CLHEP::Hep3Vector::polarAngle ( const Hep3Vector v2) const

Definition at line 26 of file SpaceVectorD.cc.

26 {
27 return std::fabs(v2.getTheta() - getTheta());
28} /* polarAngle */
double getTheta() const

◆ polarAngle() [2/2]

double CLHEP::Hep3Vector::polarAngle ( const Hep3Vector v2,
const Hep3Vector ref 
) const

Definition at line 30 of file SpaceVectorD.cc.

31 {
32 return std::fabs( v2.angle(ref) - angle(ref) );
33} /* polarAngle (v2, ref) */
double angle() const

◆ project() [1/2]

Hep3Vector CLHEP::Hep3Vector::project ( ) const
inline

Referenced by project().

◆ project() [2/2]

Hep3Vector CLHEP::Hep3Vector::project ( const Hep3Vector v2) const

Definition at line 86 of file SpaceVectorP.cc.

86 {
87 double mag2v2 = v2.mag2();
88 if (mag2v2 == 0) {
89 ZMthrowA (ZMxpvZeroVector(
90 "Attempt to take projection of vector against zero reference vector "));
91 return project();
92 }
93 return ( v2 * (dot(v2)/mag2v2) );
94}
Hep3Vector project() const

◆ pseudoRapidity()

double CLHEP::Hep3Vector::pseudoRapidity ( ) const

Definition at line 58 of file ThreeVector.cc.

58 {
59 double m1 = mag();
60 if ( m1== 0 ) return 0.0;
61 if ( m1== z() ) return 1.0E72;
62 if ( m1== -z() ) return -1.0E72;
63 return 0.5*std::log( (m1+z())/(m1-z()) );
64}
double mag() const

◆ r()

double CLHEP::Hep3Vector::r ( ) const
inline

◆ rapidity() [1/2]

double CLHEP::Hep3Vector::rapidity ( ) const

Definition at line 51 of file SpaceVectorP.cc.

51 {
52 if (std::fabs(z()) == 1) {
53 ZMthrowC (ZMxpvTachyonic(
54 "Rapidity in Z direction taken for Hep3Vector with |Z| = 1 -- \n"
55 "the log should return infinity"));
56 }
57 if (std::fabs(z()) > 1) {
58 ZMthrowA (ZMxpvTachyonic(
59 "Rapidity in Z direction taken for Hep3Vector with |Z| > 1 -- \n"
60 "the log would return a NAN" ));
61 }
62 // Want inverse std::tanh(z()):
63 return (.5 * std::log((1+z())/(1-z())) );
64}

◆ rapidity() [2/2]

double CLHEP::Hep3Vector::rapidity ( const Hep3Vector v2) const

Definition at line 96 of file SpaceVectorP.cc.

96 {
97 double vmag = v2.mag();
98 if ( vmag == 0 ) {
99 ZMthrowA (ZMxpvZeroVector(
100 "Rapidity taken with respect to zero vector" ));
101 return 0;
102 }
103 double z1 = dot(v2)/vmag;
104 if (std::fabs(z1) >= 1) {
105 ZMthrowA (ZMxpvTachyonic(
106 "Rapidity taken for too large a Hep3Vector "
107 "-- would return infinity or NAN"));
108 }
109 // Want inverse std::tanh(z1):
110 return (.5 * std::log((1+z1)/(1-z1)) );
111}

◆ rho()

double CLHEP::Hep3Vector::rho ( ) const
inline

◆ rotate() [1/5]

Hep3Vector & CLHEP::Hep3Vector::rotate ( const Hep3Vector axis,
double  delta 
)

Definition at line 26 of file SpaceVectorR.cc.

27 {
28 double r1 = axis.mag();
29 if ( r1 == 0 ) {
30 ZMthrowA (ZMxpvZeroVector(
31 "Attempt to rotate around a zero vector axis! "));
32 return *this;
33 }
34 double scale=1.0/r1;
35 double ux = scale*axis.getX();
36 double uy = scale*axis.getY();
37 double uz = scale*axis.getZ();
38 double cd = std::cos(ddelta);
39 double sd = std::sin(ddelta);
40 double ocd = 1 - cd;
41 double rx;
42 double ry;
43 double rz;
44
45 { double ocdux = ocd * ux;
46 rx = x() * ( cd + ocdux * ux ) +
47 y() * ( ocdux * uy - sd * uz ) +
48 z() * ( ocdux * uz + sd * uy ) ;
49 }
50
51 { double ocduy = ocd * uy;
52 ry = y() * ( cd + ocduy * uy ) +
53 z() * ( ocduy * uz - sd * ux ) +
54 x() * ( ocduy * ux + sd * uz ) ;
55 }
56
57 { double ocduz = ocd * uz;
58 rz = z() * ( cd + ocduz * uz ) +
59 x() * ( ocduz * ux - sd * uy ) +
60 y() * ( ocduz * uy + sd * ux ) ;
61 }
62
63 set(rx, ry, rz);
64 return *this;
65} /* rotate */
void set(double x, double y, double z)

◆ rotate() [2/5]

Hep3Vector & CLHEP::Hep3Vector::rotate ( const HepAxisAngle ax)

Definition at line 109 of file SpaceVectorR.cc.

109 {
110 return rotate( ax.getAxis(), ax.delta() );
111}
Hep3Vector & rotate(double, const Hep3Vector &)
Definition: ThreeVectorR.cc:25

◆ rotate() [3/5]

Hep3Vector & CLHEP::Hep3Vector::rotate ( const HepEulerAngles e)

Definition at line 113 of file SpaceVectorR.cc.

113 {
114 return rotate( ex.phi(), ex.theta(), ex.psi() );
115}

◆ rotate() [4/5]

Hep3Vector & CLHEP::Hep3Vector::rotate ( double  phi,
double  theta,
double  psi 
)

Definition at line 73 of file SpaceVectorR.cc.

75 {
76
77 double rx;
78 double ry;
79 double rz;
80
81 double sinPhi = std::sin( phi1 ), cosPhi = std::cos( phi1 );
82 double sinTheta = std::sin( theta1 ), cosTheta1 = std::cos( theta1 );
83 double sinPsi = std::sin( psi1 ), cosPsi = std::cos( psi1 );
84
85 rx = (cosPsi * cosPhi - cosTheta1 * sinPsi * sinPhi) * x() +
86 (cosPsi * sinPhi + cosTheta1 * sinPsi * cosPhi) * y() +
87 (sinPsi * sinTheta) * z() ;
88
89 ry = (- sinPsi * cosPhi - cosTheta1 * cosPsi * sinPhi) * x() +
90 (- sinPsi * sinPhi + cosTheta1 * cosPsi * cosPhi) * y() +
91 (cosPsi * sinTheta) * z() ;
92
93 rz = (sinTheta * sinPhi) * x() +
94 (- sinTheta * cosPhi) * y() +
95 (cosTheta1) * z() ;
96
97 set(rx, ry, rz);
98 return *this;
99
100} /* rotate */

◆ rotate() [5/5]

Hep3Vector & CLHEP::Hep3Vector::rotate ( double  angle1,
const Hep3Vector aaxis 
)

Definition at line 25 of file ThreeVectorR.cc.

25 {
26 HepRotation trans;
27 trans.rotate(angle1, aaxis);
28 operator*=(trans);
29 return *this;
30}
Hep3Vector & operator*=(double)

Referenced by main(), rotate(), CLHEP::HepLorentzVector::rotate(), and CLHEP::rotationOf().

◆ rotateUz()

Hep3Vector & CLHEP::Hep3Vector::rotateUz ( const Hep3Vector NewUzVector)

Definition at line 36 of file ThreeVector.cc.

36 {
37 // NewUzVector must be normalized !
38
39 double u1 = NewUzVector.x();
40 double u2 = NewUzVector.y();
41 double u3 = NewUzVector.z();
42 double up = u1*u1 + u2*u2;
43
44 if (up > 0) {
45 up = std::sqrt(up);
46 double px = (u1 * u3 * x() - u2 * y()) / up + u1 * z();
47 double py = (u2 * u3 * x() + u1 * y()) / up + u2 * z();
48 double pz = -up * x() + u3 * z();
49 set(px, py, pz);
50 } else if (u3 < 0.) {
51 setX(-x());
52 setZ(-z());
53 } // phi=0 teta=pi
54
55 return *this;
56}
void setZ(double)
void setX(double)

Referenced by CLHEP::HepLorentzVector::rotateUz().

◆ rotateX()

Hep3Vector & CLHEP::Hep3Vector::rotateX ( double  phi1)

Definition at line 90 of file ThreeVector.cc.

90 {
91 double sinphi = std::sin(phi1);
92 double cosphi = std::cos(phi1);
93 double ty = y() * cosphi - z() * sinphi;
94 double tz = z() * cosphi + y() * sinphi;
95 setY(ty);
96 setZ(tz);
97 return *this;
98} /* rotateX */
void setY(double)

Referenced by CLHEP::HepLorentzVector::rotateX(), and CLHEP::rotationXOf().

◆ rotateY()

Hep3Vector & CLHEP::Hep3Vector::rotateY ( double  phi1)

Definition at line 100 of file ThreeVector.cc.

100 {
101 double sinphi = std::sin(phi1);
102 double cosphi = std::cos(phi1);
103 double tx = x() * cosphi + z() * sinphi;
104 double tz = z() * cosphi - x() * sinphi;
105 setX(tx);
106 setZ(tz);
107 return *this;
108} /* rotateY */

Referenced by main(), CLHEP::HepLorentzVector::rotateY(), and CLHEP::rotationYOf().

◆ rotateZ()

Hep3Vector & CLHEP::Hep3Vector::rotateZ ( double  phi1)

Definition at line 110 of file ThreeVector.cc.

110 {
111 double sinphi = std::sin(phi1);
112 double cosphi = std::cos(phi1);
113 double tx = x() * cosphi - y() * sinphi;
114 double ty = y() * cosphi + x() * sinphi;
115 setX(tx);
116 setY(ty);
117 return *this;
118} /* rotateZ */

Referenced by main(), CLHEP::HepLorentzVector::rotateZ(), and CLHEP::rotationZOf().

◆ set()

◆ setCylEta()

void CLHEP::Hep3Vector::setCylEta ( double  p)

Definition at line 254 of file ThreeVector.cc.

254 {
255
256 // In cylindrical coords, set eta while keeping rho and phi fixed
257
258 double theta1 = 2 * std::atan ( std::exp (-eta1) );
259
260 //-| The remaining code is similar to setCylTheta, The reason for
261 //-| using a copy is so as to be able to change the messages in the
262 //-| ZMthrows to say eta rather than theta. Besides, we assumedly
263 //-| need not check for theta of 0 or PI.
264
265 if ( (x() == 0) && (y() == 0) ) {
266 if (z() == 0) {
267 ZMthrowC (ZMxpvZeroVector(
268 "Attempt to set cylEta of zero vector -- vector is unchanged"));
269 return;
270 }
271 if (theta1 == 0) {
272 setZ(std::fabs(z()));
273 return;
274 }
275 if (theta1 == CLHEP::pi) {
276 setZ(-std::fabs(z()));
277 return;
278 }
279 ZMthrowC (ZMxpvZeroVector(
280 "Attempt set cylindrical eta of vector along Z axis "
281 "to a non-trivial value, while keeping rho fixed -- "
282 "will return zero vector"));
283 setZ(0.0);
284 return;
285 }
286 double phi1 (getPhi());
287 double rho1 = getRho();
288 setZ(rho1 / std::tan (theta1));
289 setY(rho1 * std::sin (phi1));
290 setX(rho1 * std::cos (phi1));
291
292} /* setCylEta */
double getRho() const

◆ setCylindrical()

void CLHEP::Hep3Vector::setCylindrical ( double  r,
double  phi,
double  z 
)
protected

Definition at line 54 of file SpaceVector.cc.

57 {
58 if ( rho1 < 0 ) {
59 ZMthrowC (ZMxpvNegativeR(
60 "Cylindrical coordinates supplied with negative Rho"));
61 // No special return needed if warning is ignored.
62 }
63 setZ(z1);
64 setY(rho1 * std::sin (phi1));
65 setX(rho1 * std::cos (phi1));
66 return;
67} /* setCylindrical (r, phi, z) */

◆ setCylTheta()

void CLHEP::Hep3Vector::setCylTheta ( double  theta1)

Definition at line 209 of file ThreeVector.cc.

209 {
210
211 // In cylindrical coords, set theta while keeping rho and phi fixed
212
213 if ( (x() == 0) && (y() == 0) ) {
214 if (z() == 0) {
215 ZMthrowC (ZMxpvZeroVector(
216 "Attempt to set cylTheta of zero vector -- vector is unchanged"));
217 return;
218 }
219 if (theta1 == 0) {
220 setZ(std::fabs(z()));
221 return;
222 }
223 if (theta1 == CLHEP::pi) {
224 setZ(-std::fabs(z()));
225 return;
226 }
227 ZMthrowC (ZMxpvZeroVector(
228 "Attempt set cylindrical theta of vector along Z axis "
229 "to a non-trivial value, while keeping rho fixed -- "
230 "will return zero vector"));
231 setZ(0.0);
232 return;
233 }
234 if ( (theta1 < 0) || (theta1 > CLHEP::pi) ) {
235 ZMthrowC (ZMxpvUnusualTheta(
236 "Setting Cyl theta of a vector based on a value not in [0, PI]"));
237 // No special return needed if warning is ignored.
238 }
239 double phi1 (getPhi());
240 double rho1 = getRho();
241 if ( (theta1 == 0) || (theta1 == CLHEP::pi) ) {
242 ZMthrowC (ZMxpvInfiniteVector(
243 "Attempt to set cylindrical theta to 0 or PI "
244 "while keeping rho fixed -- infinite Z will be computed"));
245 setZ((theta1==0) ? 1.0E72 : -1.0E72);
246 return;
247 }
248 setZ(rho1 / std::tan (theta1));
249 setY(rho1 * std::sin (phi1));
250 setX(rho1 * std::cos (phi1));
251
252} /* setCylTheta */

◆ setEta()

void CLHEP::Hep3Vector::setEta ( double  p)

Definition at line 183 of file ThreeVector.cc.

183 {
184 double phi1 = 0;
185 double r1;
186 if ( (x() == 0) && (y() == 0) ) {
187 if (z() == 0) {
188 ZMthrowC (ZMxpvZeroVector(
189 "Attempt to set eta of zero vector -- vector is unchanged"));
190 return;
191 }
192 ZMthrowC (ZMxpvZeroVector(
193 "Attempt to set eta of vector along Z axis -- will use phi = 0"));
194 r1 = std::fabs(z());
195 } else {
196 r1 = getR();
197 phi1 = getPhi();
198 }
199 double tanHalfTheta = std::exp ( -eta1 );
200 double cosTheta1 =
201 (1 - tanHalfTheta*tanHalfTheta) / (1 + tanHalfTheta*tanHalfTheta);
202 double rho1 = r1*std::sqrt(1 - cosTheta1*cosTheta1);
203 setZ(r1 * cosTheta1);
204 setY(rho1 * std::sin (phi1));
205 setX(rho1 * std::cos (phi1));
206 return;
207}

◆ setMag()

void CLHEP::Hep3Vector::setMag ( double  ma)

Definition at line 23 of file ThreeVector.cc.

23 {
24 double factor = mag();
25 if (factor == 0) {
26 ZMthrowA ( ZMxpvZeroVector (
27 "Hep3Vector::setMag : zero vector can't be stretched"));
28 }else{
29 factor = ma/factor;
30 setX(x()*factor);
31 setY(y()*factor);
32 setZ(z()*factor);
33 }
34}

◆ setPerp()

void CLHEP::Hep3Vector::setPerp ( double  )
inline

◆ setPhi()

void CLHEP::Hep3Vector::setPhi ( double  )
inline

◆ setR()

void CLHEP::Hep3Vector::setR ( double  s)
inline

◆ setREtaPhi()

void CLHEP::Hep3Vector::setREtaPhi ( double  r,
double  eta,
double  phi 
)
inline

◆ setRho()

void CLHEP::Hep3Vector::setRho ( double  s)
inline

◆ setRhoPhiEta()

void CLHEP::Hep3Vector::setRhoPhiEta ( double  rho,
double  phi,
double  eta 
)

Definition at line 96 of file SpaceVector.cc.

99 {
100 if (rho1 == 0) {
101 ZMthrowC (ZMxpvZeroVector(
102 "Attempt set vector components rho, phi, eta with zero rho -- "
103 "zero vector is returned, ignoring eta and phi"));
104 set(0.0, 0.0, 0.0);
105 return;
106 }
107 double theta1 (2 * std::atan ( std::exp (-eta1) ));
108 setZ(rho1 / std::tan (theta1));
109 setY(rho1 * std::sin (phi1));
110 setX(rho1 * std::cos (phi1));
111 return;
112} /* setCyl (rho, phi, eta) */

◆ setRhoPhiTheta()

void CLHEP::Hep3Vector::setRhoPhiTheta ( double  rho,
double  phi,
double  theta 
)

Definition at line 69 of file SpaceVector.cc.

72 {
73 if (rho1 == 0) {
74 ZMthrowC (ZMxpvZeroVector(
75 "Attempt set vector components rho, phi, theta with zero rho -- "
76 "zero vector is returned, ignoring theta and phi"));
77 set(0.0, 0.0, 0.0);
78 return;
79 }
80 if ( (theta1 == 0) || (theta1 == CLHEP::pi) ) {
81 ZMthrowA (ZMxpvInfiniteVector(
82 "Attempt set cylindrical vector vector with finite rho and "
83 "theta along the Z axis: infinite Z would be computed"));
84 }
85 if ( (theta1 < 0) || (theta1 > CLHEP::pi) ) {
86 ZMthrowC (ZMxpvUnusualTheta(
87 "Rho, phi, theta set with theta not in [0, PI]"));
88 // No special return needed if warning is ignored.
89 }
90 setZ(rho1 / std::tan (theta1));
91 setY(rho1 * std::sin (phi1));
92 setX(rho1 * std::cos (phi1));
93 return;
94} /* setCyl (rho, phi, theta) */

◆ setRhoPhiZ()

void CLHEP::Hep3Vector::setRhoPhiZ ( double  rho,
double  phi,
double  z 
)
inline

◆ setRThetaPhi()

void CLHEP::Hep3Vector::setRThetaPhi ( double  r,
double  theta,
double  phi 
)
inline

◆ setSpherical()

void CLHEP::Hep3Vector::setSpherical ( double  r,
double  theta,
double  phi 
)
protected

Definition at line 33 of file SpaceVector.cc.

36 {
37 if ( r1 < 0 ) {
38 ZMthrowC (ZMxpvNegativeR(
39 "Spherical coordinates set with negative R"));
40 // No special return needed if warning is ignored.
41 }
42 if ( (theta1 < 0) || (theta1 > CLHEP::pi) ) {
43 ZMthrowC (ZMxpvUnusualTheta(
44 "Spherical coordinates set with theta not in [0, PI]"));
45 // No special return needed if warning is ignored.
46 }
47 double rho1 ( r1*std::sin(theta1));
48 setZ(r1 * std::cos(theta1));
49 setY(rho1 * std::sin (phi1));
50 setX(rho1 * std::cos (phi1));
51 return;
52} /* setSpherical (r, theta1, phi) */

◆ setTheta()

void CLHEP::Hep3Vector::setTheta ( double  )
inline

◆ setTolerance()

double CLHEP::Hep3Vector::setTolerance ( double  tol)
static

Definition at line 271 of file SpaceVector.cc.

271 {
272// Set the tolerance for Hep3Vectors to be considered near one another
273 double oldTolerance (tolerance);
274 tolerance = tol;
275 return oldTolerance;
276}
static double tolerance
Definition: ThreeVector.h:394

◆ setX()

◆ setY()

◆ setZ()

◆ theta() [1/2]

double CLHEP::Hep3Vector::theta ( ) const
inline

Referenced by main().

◆ theta() [2/2]

double CLHEP::Hep3Vector::theta ( const Hep3Vector v2) const
inline

◆ transform()

Hep3Vector & CLHEP::Hep3Vector::transform ( const HepRotation m1)

Definition at line 21 of file ThreeVectorR.cc.

21 {
22 return *this = m1 * (*this);
23}

◆ unit()

◆ x()

◆ y()

◆ z()

Member Data Documentation

◆ data

double CLHEP::Hep3Vector::data[3]
protected

Definition at line 391 of file ThreeVector.h.

◆ tolerance

double CLHEP::Hep3Vector::tolerance = Hep3Vector::ToleranceTicks * 2.22045e-16
staticprotected

Definition at line 394 of file ThreeVector.h.

Referenced by setTolerance().

◆ ToleranceTicks

const int CLHEP::Hep3Vector::ToleranceTicks = 100
static

Definition at line 295 of file ThreeVector.h.


The documentation for this class was generated from the following files: