Geant4-11
ThreeVector.cc
Go to the documentation of this file.
1// -*- C++ -*-
2// ---------------------------------------------------------------------------
3//
4// This file is a part of the CLHEP - a Class Library for High Energy Physics.
5//
6// This is the implementation of the Hep3Vector class.
7//
8// See also ThreeVectorR.cc for implementation of Hep3Vector methods which
9// would couple in all the HepRotation methods.
10//
11
14
15#include <cmath>
16#include <iostream>
17
18namespace CLHEP {
19
20void Hep3Vector::setMag(double ma) {
21 double factor = mag();
22 if (factor == 0) {
23 std::cerr << "Hep3Vector::setMag() - "
24 << "zero vector can't be stretched" << std::endl;
25 }else{
26 factor = ma/factor;
27 setX(x()*factor);
28 setY(y()*factor);
29 setZ(z()*factor);
30 }
31}
32
34 // NewUzVector must be normalized !
35
36 double u1 = NewUzVector.x();
37 double u2 = NewUzVector.y();
38 double u3 = NewUzVector.z();
39 double up = u1*u1 + u2*u2;
40
41 if (up > 0) {
42 up = std::sqrt(up);
43 double px = (u1 * u3 * x() - u2 * y()) / up + u1 * z();
44 double py = (u2 * u3 * x() + u1 * y()) / up + u2 * z();
45 double pz = -up * x() + u3 * z();
46 set(px, py, pz);
47 } else if (u3 < 0.) {
48 setX(-x());
49 setZ(-z());
50 } // phi=0 teta=pi
51
52 return *this;
53}
54
56 double m1 = mag();
57 if ( m1== 0 ) return 0.0;
58 if ( m1== z() ) return 1.0E72;
59 if ( m1== -z() ) return -1.0E72;
60 return 0.5*std::log( (m1+z())/(m1-z()) );
61}
62
63std::ostream & operator<< (std::ostream & os, const Hep3Vector & v) {
64 return os << "(" << v.x() << "," << v.y() << "," << v.z() << ")";
65}
66
67void ZMinput3doubles ( std::istream & is, const char * type,
68 double & x, double & y, double & z );
69
70std::istream & operator>>(std::istream & is, Hep3Vector & v) {
71 double x, y, z;
72 ZMinput3doubles ( is, "Hep3Vector", x, y, z );
73 v.set(x, y, z);
74 return is;
75} // operator>>()
76
77const Hep3Vector HepXHat(1.0, 0.0, 0.0);
78const Hep3Vector HepYHat(0.0, 1.0, 0.0);
79const Hep3Vector HepZHat(0.0, 0.0, 1.0);
80
81//-------------------
82//
83// New methods introduced when ZOOM PhysicsVectors was merged in:
84//
85//-------------------
86
88 double sinphi = std::sin(phi1);
89 double cosphi = std::cos(phi1);
90 double ty = y() * cosphi - z() * sinphi;
91 double tz = z() * cosphi + y() * sinphi;
92 setY(ty);
93 setZ(tz);
94 return *this;
95} /* rotateX */
96
98 double sinphi = std::sin(phi1);
99 double cosphi = std::cos(phi1);
100 double tx = x() * cosphi + z() * sinphi;
101 double tz = z() * cosphi - x() * sinphi;
102 setX(tx);
103 setZ(tz);
104 return *this;
105} /* rotateY */
106
108 double sinphi = std::sin(phi1);
109 double cosphi = std::cos(phi1);
110 double tx = x() * cosphi - y() * sinphi;
111 double ty = y() * cosphi + x() * sinphi;
112 setX(tx);
113 setY(ty);
114 return *this;
115} /* rotateZ */
116
117bool Hep3Vector::isNear(const Hep3Vector & v, double epsilon) const {
118 double limit = dot(v)*epsilon*epsilon;
119 return ( (*this - v).mag2() <= limit );
120} /* isNear() */
121
122double Hep3Vector::howNear(const Hep3Vector & v ) const {
123 // | V1 - V2 | **2 / V1 dot V2, up to 1
124 double d = (*this - v).mag2();
125 double vdv = dot(v);
126 if ( (vdv > 0) && (d < vdv) ) {
127 return std::sqrt (d/vdv);
128 } else if ( (vdv == 0) && (d == 0) ) {
129 return 0;
130 } else {
131 return 1;
132 }
133} /* howNear */
134
135double Hep3Vector::deltaPhi (const Hep3Vector & v2) const {
136 double dphi = v2.getPhi() - getPhi();
137 if ( dphi > CLHEP::pi ) {
138 dphi -= CLHEP::twopi;
139 } else if ( dphi <= -CLHEP::pi ) {
140 dphi += CLHEP::twopi;
141 }
142 return dphi;
143} /* deltaPhi */
144
145double Hep3Vector::deltaR ( const Hep3Vector & v ) const {
146 double a = eta() - v.eta();
147 double b = deltaPhi(v);
148 return std::sqrt ( a*a + b*b );
149} /* deltaR */
150
151double Hep3Vector::cosTheta(const Hep3Vector & q) const {
152 double arg;
153 double ptot2 = mag2()*q.mag2();
154 if(ptot2 <= 0) {
155 arg = 0.0;
156 }else{
157 arg = dot(q)/std::sqrt(ptot2);
158 if(arg > 1.0) arg = 1.0;
159 if(arg < -1.0) arg = -1.0;
160 }
161 return arg;
162}
163
164double Hep3Vector::cos2Theta(const Hep3Vector & q) const {
165 double arg;
166 double ptot2 = mag2();
167 double qtot2 = q.mag2();
168 if ( ptot2 == 0 || qtot2 == 0 ) {
169 arg = 1.0;
170 }else{
171 double pdq = dot(q);
172 arg = (pdq/ptot2) * (pdq/qtot2);
173 // More naive methods overflow on vectors which can be squared
174 // but can't be raised to the 4th power.
175 if(arg > 1.0) arg = 1.0;
176 }
177 return arg;
178}
179
180void Hep3Vector::setEta (double eta1) {
181 double phi1 = 0;
182 double r1;
183 if ( (x() == 0) && (y() == 0) ) {
184 if (z() == 0) {
185 std::cerr << "Hep3Vector::setEta() - "
186 << "Attempt to set eta of zero vector -- vector is unchanged"
187 << std::endl;
188 return;
189 }
190 std::cerr << "Hep3Vector::setEta() - "
191 << "Attempt to set eta of vector along Z axis -- will use phi = 0"
192 << std::endl;
193 r1 = std::fabs(z());
194 } else {
195 r1 = getR();
196 phi1 = getPhi();
197 }
198 double tanHalfTheta = std::exp ( -eta1 );
199 double cosTheta1 =
200 (1 - tanHalfTheta*tanHalfTheta) / (1 + tanHalfTheta*tanHalfTheta);
201 double rho1 = r1*std::sqrt(1 - cosTheta1*cosTheta1);
202 setZ(r1 * cosTheta1);
203 setY(rho1 * std::sin (phi1));
204 setX(rho1 * std::cos (phi1));
205 return;
206}
207
208void Hep3Vector::setCylTheta (double theta1) {
209
210 // In cylindrical coords, set theta while keeping rho and phi fixed
211
212 if ( (x() == 0) && (y() == 0) ) {
213 if (z() == 0) {
214 std::cerr << "Hep3Vector::setCylTheta() - "
215 << "Attempt to set cylTheta of zero vector -- vector is unchanged"
216 << std::endl;
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 std::cerr << "Hep3Vector::setCylTheta() - "
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" << std::endl;
231 setZ(0.0);
232 return;
233 }
234 if ( (theta1 < 0) || (theta1 > CLHEP::pi) ) {
235 std::cerr << "Hep3Vector::setCylTheta() - "
236 << "Setting Cyl theta of a vector based on a value not in [0, PI]"
237 << std::endl;
238 // No special return needed if warning is ignored.
239 }
240 double phi1 (getPhi());
241 double rho1 = getRho();
242 if ( (theta1 == 0) || (theta1 == CLHEP::pi) ) {
243 std::cerr << "Hep3Vector::setCylTheta() - "
244 << "Attempt to set cylindrical theta to 0 or PI "
245 << "while keeping rho fixed -- infinite Z will be computed"
246 << std::endl;
247 setZ((theta1==0) ? 1.0E72 : -1.0E72);
248 return;
249 }
250 setZ(rho1 / std::tan (theta1));
251 setY(rho1 * std::sin (phi1));
252 setX(rho1 * std::cos (phi1));
253
254} /* setCylTheta */
255
256void Hep3Vector::setCylEta (double eta1) {
257
258 // In cylindrical coords, set eta while keeping rho and phi fixed
259
260 double theta1 = 2 * std::atan ( std::exp (-eta1) );
261
262 //-| The remaining code is similar to setCylTheta, The reason for
263 //-| using a copy is so as to be able to change the messages in the
264 //-| ZMthrows to say eta rather than theta. Besides, we assumedly
265 //-| need not check for theta of 0 or PI.
266
267 if ( (x() == 0) && (y() == 0) ) {
268 if (z() == 0) {
269 std::cerr << "Hep3Vector::setCylEta() - "
270 << "Attempt to set cylEta of zero vector -- vector is unchanged"
271 << std::endl;
272 return;
273 }
274 if (theta1 == 0) {
275 setZ(std::fabs(z()));
276 return;
277 }
278 if (theta1 == CLHEP::pi) {
279 setZ(-std::fabs(z()));
280 return;
281 }
282 std::cerr << "Hep3Vector::setCylEta() - "
283 << "Attempt set cylindrical eta of vector along Z axis "
284 << "to a non-trivial value, while keeping rho fixed -- "
285 << "will return zero vector" << std::endl;
286 setZ(0.0);
287 return;
288 }
289 double phi1 (getPhi());
290 double rho1 = getRho();
291 setZ(rho1 / std::tan (theta1));
292 setY(rho1 * std::sin (phi1));
293 setX(rho1 * std::cos (phi1));
294
295} /* setCylEta */
296
297
298Hep3Vector operator/ ( const Hep3Vector & v1, double c ) {
299// if (c == 0) {
300// std::cerr << "Hep3Vector::operator/ () - "
301// << "Attempt to divide vector by 0 -- "
302// << "will produce infinities and/or NANs" << std::endl;
303// }
304 return v1 * (1.0/c);
305} /* v / c */
306
308// if (c == 0) {
309// std::cerr << "Hep3Vector::operator/ () - "
310// << "Attempt to do vector /= 0 -- "
311// << "division by zero would produce infinite or NAN components"
312// << std::endl;
313// }
314 *this *= 1.0/c;
315 return *this;
316}
317
319
320} // namespace CLHEP
G4double epsilon(G4double density, G4double temperature)
Hep3Vector & rotateY(double)
Definition: ThreeVector.cc:97
double z() const
void setEta(double p)
Definition: ThreeVector.cc:180
Hep3Vector & rotateX(double)
Definition: ThreeVector.cc:87
double eta() const
double cos2Theta() const
static DLL_API double tolerance
Definition: ThreeVector.h:394
double x() const
void setY(double)
double mag2() const
double getR() const
Hep3Vector & rotateZ(double)
Definition: ThreeVector.cc:107
Hep3Vector & operator/=(double)
Definition: ThreeVector.cc:307
double y() const
void setCylEta(double p)
Definition: ThreeVector.cc:256
double howNear(const Hep3Vector &v) const
Definition: ThreeVector.cc:122
double dot(const Hep3Vector &) const
void setZ(double)
double mag() const
bool isNear(const Hep3Vector &, double epsilon=tolerance) const
Definition: ThreeVector.cc:117
double pseudoRapidity() const
Definition: ThreeVector.cc:55
double getRho() const
double deltaPhi(const Hep3Vector &v2) const
Definition: ThreeVector.cc:135
void set(double x, double y, double z)
void setMag(double)
Definition: ThreeVector.cc:20
double getPhi() const
double deltaR(const Hep3Vector &v) const
Definition: ThreeVector.cc:145
void setX(double)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
double cosTheta() const
void setCylTheta(double)
Definition: ThreeVector.cc:208
Definition: DoubConv.h:17
DLL_API const Hep3Vector HepZHat
Definition: ThreeVector.h:419
std::istream & operator>>(std::istream &is, HepRandom &dist)
Definition: Random.cc:223
void ZMinput3doubles(std::istream &is, const char *type, double &x, double &y, double &z)
Definition: ZMinput.cc:42
static constexpr double twopi
Definition: SystemOfUnits.h:56
DLL_API const Hep3Vector HepXHat
DLL_API const Hep3Vector HepYHat
Definition: ThreeVector.h:419
std::ostream & operator<<(std::ostream &os, const HepRandom &dist)
Definition: Random.cc:219
static constexpr double pi
Definition: SystemOfUnits.h:55
HepLorentzVector operator/(const HepLorentzVector &, double a)