Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
LorentzVectorK.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 part of the implementation of the HepLorentzVector class:
7// Those methods which originated from ZOOM and which deal with relativistic
8// kinematic properties.
9//
10
12
13#include <cmath>
14#include <iostream>
15
16namespace CLHEP {
17
18//-******************
19// Metric flexibility
20//-******************
21
23 ZMpvMetric_t oldMetric = (metric > 0) ? TimePositive : TimeNegative;
24 if ( a1 == TimeNegative ) {
25 metric = -1.0;
26 } else {
27 metric = 1.0;
28 }
29 return oldMetric;
30}
31
33 return ( (metric > 0) ? TimePositive : TimeNegative );
34}
35
36//-********
37// plus
38// minus
39//-********
40
41double HepLorentzVector::plus (const Hep3Vector & ref) const {
42 double r = ref.mag();
43 if (r == 0) {
44 std::cerr << "HepLorentzVector::plus() - "
45 << "A zero vector used as reference to LorentzVector plus-part"
46 << std::endl;
47 return ee;
48 }
49 return ee + pp.dot(ref)/r;
50} /* plus */
51
52double HepLorentzVector::minus (const Hep3Vector & ref) const {
53 double r = ref.mag();
54 if (r == 0) {
55 std::cerr << "HepLorentzVector::minus() - "
56 << "A zero vector used as reference to LorentzVector minus-part"
57 << std::endl;
58 return ee;
59 }
60 return ee - pp.dot(ref)/r;
61} /* plus */
62
64 return HepLorentzVector (0, 0, 0, (t() < 0.0 ? -m() : m()));
65}
66
67//-********
68// beta
69// gamma
70//-********
71
72double HepLorentzVector::beta() const {
73 if (ee == 0) {
74 if (pp.mag2() == 0) {
75 return 0;
76 } else {
77 std::cerr << "HepLorentzVector::beta() - "
78 << "beta computed for HepLorentzVector with t=0 -- infinite result"
79 << std::endl;
80 return 1./ee;
81 }
82 }
83// if (restMass2() <= 0) {
84// std::cerr << "HepLorentzVector::beta() - "
85// << "beta computed for a non-timelike HepLorentzVector" << std::endl;
86// // result will make analytic sense but is physically meaningless
87// }
88 return std::sqrt (pp.mag2() / (ee*ee)) ;
89} /* beta */
90
92 double v2 = pp.mag2();
93 double t2 = ee*ee;
94 if (ee == 0) {
95 if (pp.mag2() == 0) {
96 return 1;
97 } else {
98 std::cerr << "HepLorentzVector::gamma() - "
99 << "gamma computed for HepLorentzVector with t=0 -- zero result"
100 << std::endl;
101 return 0;
102 }
103 }
104 if (t2 < v2) {
105 std::cerr << "HepLorentzVector::gamma() - "
106 << "gamma computed for a spacelike HepLorentzVector -- imaginary result"
107 << std::endl;
108 // analytic result would be imaginary.
109 return 0;
110// } else if ( t2 == v2 ) {
111// std::cerr << "HepLorentzVector::gamma() - "
112// << "gamma computed for a lightlike HepLorentzVector -- infinite result"
113// << std::endl;
114 }
115 return 1./std::sqrt(1. - v2/t2 );
116} /* gamma */
117
118
119//-***************
120// rapidity
121// pseudorapidity
122// eta
123//-***************
124
126 double z1 = pp.getZ();
127// if (std::fabs(ee) == std::fabs(z1)) {
128// std::cerr << "HepLorentzVector::rapidity() - "
129// << "rapidity for 4-vector with |E| = |Pz| -- infinite result"
130// << std::endl;
131// }
132 if (std::fabs(ee) < std::fabs(z1)) {
133 std::cerr << "HepLorentzVector::rapidity() - "
134 << "rapidity for spacelike 4-vector with |E| < |Pz| -- undefined"
135 << std::endl;
136 return 0;
137 }
138 double q = (ee + z1) / (ee - z1);
139 //-| This cannot be negative now, since both numerator
140 //-| and denominator have the same sign as ee.
141 return .5 * std::log(q);
142} /* rapidity */
143
144double HepLorentzVector::rapidity(const Hep3Vector & ref) const {
145 double r = ref.mag2();
146 if (r == 0) {
147 std::cerr << "HepLorentzVector::rapidity() - "
148 << "A zero vector used as reference to LorentzVector rapidity"
149 << std::endl;
150 return 0;
151 }
152 double vdotu = pp.dot(ref)/std::sqrt(r);
153// if (std::fabs(ee) == std::fabs(vdotu)) {
154// std::cerr << "HepLorentzVector::rapidity() - "
155// << "rapidity for 4-vector with |E| = |Pu| -- infinite result"
156// << std::endl;
157// }
158 if (std::fabs(ee) < std::fabs(vdotu)) {
159 std::cerr << "HepLorentzVector::rapidity() - "
160 << "rapidity for spacelike 4-vector with |E| < |P*ref| -- undefined "
161 << std::endl;
162 return 0;
163 }
164 double q = (ee + vdotu) / (ee - vdotu);
165 return .5 * std::log(q);
166} /* rapidity(ref) */
167
169 double v1 = pp.mag();
170// if (std::fabs(ee) == std::fabs(v1)) {
171// std::cerr << "HepLorentzVector::coLinearRapidity() - "
172// << "co-Linear rapidity for 4-vector with |E| = |P| -- infinite result"
173// << std::endl;
174// }
175 if (std::fabs(ee) < std::fabs(v1)) {
176 std::cerr << "HepLorentzVector::coLinearRapidity() - "
177 << "co-linear rapidity for spacelike 4-vector -- undefined"
178 << std::endl;
179 return 0;
180 }
181 double q = (ee + v1) / (ee - v1);
182 return .5 * std::log(q);
183} /* rapidity */
184
185//-*************
186// invariantMass
187//-*************
188
190 double m1 = invariantMass2(w);
191 if (m1 < 0) {
192 // We should find out why:
193 if ( ee * w.ee < 0 ) {
194 std::cerr << "HepLorentzVector::invariantMass() - "
195 << "invariant mass meaningless: \n"
196 << "a negative-mass input led to spacelike 4-vector sum" << std::endl;
197 return 0;
198 } else if ( (isSpacelike() && !isLightlike()) ||
199 (w.isSpacelike() && !w.isLightlike()) ) {
200 std::cerr << "HepLorentzVector::invariantMass() - "
201 << "invariant mass meaningless because of spacelike input"
202 << std::endl;
203 return 0;
204 } else {
205 // Invariant mass squared for a pair of timelike or lightlike vectors
206 // mathematically cannot be negative. If the vectors are within the
207 // tolerance of being lightlike or timelike, we can assume that prior
208 // or current roundoffs have caused the negative result, and return 0
209 // without comment.
210 return 0;
211 }
212 }
213 return (ee+w.ee >=0 ) ? std::sqrt(m1) : - std::sqrt(m1);
214} /* invariantMass */
215
216//-***************
217// findBoostToCM
218//-***************
219
221 return -boostVector();
222} /* boostToCM() */
223
225 double t1 = ee + w.ee;
226 Hep3Vector v1 = pp + w.pp;
227 if (t1 == 0) {
228 if (v1.mag2() == 0) {
229 return Hep3Vector(0,0,0);
230 } else {
231 std::cerr << "HepLorentzVector::findBoostToCM() - "
232 << "boostToCM computed for two 4-vectors with combined t=0 -- "
233 << "infinite result" << std::endl;
234 return Hep3Vector(v1*(1./t1)); // Yup, 1/0 -- that is how we return infinity
235 }
236 }
237// if (t1*t1 - v1.mag2() <= 0) {
238// std::cerr << "HepLorentzVector::findBoostToCM() - "
239// << "boostToCM computed for pair of HepLorentzVectors with non-timelike sum"
240// << std::endl;
241// // result will make analytic sense but is physically meaningless
242// }
243 return Hep3Vector(v1 * (-1./t1));
244} /* boostToCM(w) */
245
246} // namespace CLHEP
247
double getZ() const
double mag2() const
double dot(const Hep3Vector &) const
double mag() const
double invariantMass() const
Hep3Vector boostVector() const
HepLorentzVector rest4Vector() const
static ZMpvMetric_t getMetric()
double minus() const
bool isLightlike(double epsilon=tolerance) const
double invariantMass2() const
double coLinearRapidity() const
Hep3Vector findBoostToCM() const
static ZMpvMetric_t setMetric(ZMpvMetric_t met)
bool isSpacelike() const
Definition: DoubConv.h:17
@ TimePositive
Definition: LorentzVector.h:59
@ TimeNegative
Definition: LorentzVector.h:59