ECCE @ EIC Software
Reference for
ECCE @ EIC
simulation and reconstruction software on GitHub
Home page
Related Pages
Modules
Namespaces
Classes
Files
External Links
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Groups
Pages
G4LorentzConvertor.cc
Go to the documentation of this file.
Or view
the newest version in sPHENIX GitHub for file G4LorentzConvertor.cc
1
//
2
// ********************************************************************
3
// * License and Disclaimer *
4
// * *
5
// * The Geant4 software is copyright of the Copyright Holders of *
6
// * the Geant4 Collaboration. It is provided under the terms and *
7
// * conditions of the Geant4 Software License, included in the file *
8
// * LICENSE and available at http://cern.ch/geant4/license . These *
9
// * include a list of copyright holders. *
10
// * *
11
// * Neither the authors of this software system, nor their employing *
12
// * institutes,nor the agencies providing financial support for this *
13
// * work make any representation or warranty, express or implied, *
14
// * regarding this software system or assume any liability for its *
15
// * use. Please see the license in the file LICENSE and URL above *
16
// * for the full disclaimer and the limitation of liability. *
17
// * *
18
// * This code implementation is the result of the scientific and *
19
// * technical work of the GEANT4 collaboration. *
20
// * By using, copying, modifying or distributing the software (or *
21
// * any work based on the software) you agree to acknowledge its *
22
// * use in resulting scientific publications, and indicate your *
23
// * acceptance of all terms of the Geant4 Software license. *
24
// ********************************************************************
25
//
26
//
27
// 20100108 Michael Kelsey -- Use G4LorentzVector internally
28
// 20100112 M. Kelsey -- Remove G4CascadeMomentum, use G4LorentzVector directly
29
// 20100308 M. Kelsey -- Bug fix in toTheTargetRestFrame: scm_momentum should
30
// be data member, not local.
31
// 20100409 M. Kelsey -- Protect std::sqrt(ga) against round-off negatives
32
// 20100519 M. Kelsey -- Add interfaces to pass G4InuclParticles directly
33
// 20100617 M. Kelsey -- Add more diagnostic messages with multiple levels
34
// 20100713 M. Kelsey -- All diagnostic messages should be verbose > 1
35
// 20110525 M. Kelsey -- Add some diagnostic messages in both ::rotate()
36
// 20110602 M. Kelsey -- Simplify some kinematics, dropping intermediate calcs
37
38
#include "
G4LorentzConvertor.hh
"
39
#include "
G4ThreeVector.hh
"
40
#include "
G4HadronicException.hh
"
41
#include "
G4InuclParticle.hh
"
42
43
44
const
G4double
G4LorentzConvertor::small
= 1.0e-10;
45
46
G4LorentzConvertor::G4LorentzConvertor
()
47
: verboseLevel(0),
v2
(0.), ecm_tot(0.), valong(0.), degenerated(
false
) {}
48
49
G4LorentzConvertor::
50
G4LorentzConvertor
(
const
G4LorentzVector
& bmom,
G4double
bmass,
51
const
G4LorentzVector
& tmom,
G4double
tmass)
52
: verboseLevel(0),
v2
(0.), ecm_tot(0.), valong(0.), degenerated(
false
) {
53
setBullet
(bmom, bmass);
54
setTarget
(tmom, tmass);
55
}
56
57
G4LorentzConvertor::
58
G4LorentzConvertor
(
const
G4InuclParticle
* bullet,
59
const
G4InuclParticle
*
target
)
60
: verboseLevel(0),
v2
(0.), ecm_tot(0.), valong(0.), degenerated(
false
) {
61
setBullet
(bullet);
62
setTarget
(target);
63
}
64
65
// Extract four-vectors from input particles
66
void
G4LorentzConvertor::setBullet
(
const
G4InuclParticle
* bullet) {
67
setBullet
(bullet->
getMomentum
());
68
}
69
70
void
G4LorentzConvertor::setTarget
(
const
G4InuclParticle
*
target
) {
71
setTarget
(target->
getMomentum
());
72
}
73
74
75
// Boost bullet and target four-vectors into desired frame
76
77
void
G4LorentzConvertor::toTheCenterOfMass
() {
78
if
(
verboseLevel
> 2)
79
G4cout
<<
" >>> G4LorentzConvertor::toTheCenterOfMass"
<<
G4endl
;
80
81
velocity
= (
target_mom
+
bullet_mom
).boostVector();
82
if
(
verboseLevel
> 3)
G4cout
<<
" boost "
<<
velocity
<<
G4endl
;
83
84
// "SCM" is reverse target momentum in the CM frame
85
scm_momentum
=
target_mom
;
86
scm_momentum
.
boost
(-
velocity
);
87
scm_momentum
.
setVect
(-
scm_momentum
.
vect
());
88
89
if
(
verboseLevel
> 3)
G4cout
<<
" pscm "
<<
scm_momentum
.
vect
() <<
G4endl
;
90
91
fillKinematics
();
92
}
93
94
void
G4LorentzConvertor::toTheTargetRestFrame
() {
95
if
(
verboseLevel
> 2)
96
G4cout
<<
" >>> G4LorentzConvertor::toTheTargetRestFrame"
<<
G4endl
;
97
98
velocity
=
target_mom
.
boostVector
();
99
if
(
verboseLevel
> 3)
G4cout
<<
" boost "
<<
velocity
<<
G4endl
;
100
101
// "SCM" is bullet momentum in the target's frame
102
scm_momentum
=
bullet_mom
;
103
scm_momentum
.
boost
(-
velocity
);
104
105
if
(
verboseLevel
> 3)
G4cout
<<
" pseudo-pscm "
<<
scm_momentum
.
vect
() <<
G4endl
;
106
107
fillKinematics
();
108
}
109
110
// Compute kinematic quantities for rotate() functions
111
112
void
G4LorentzConvertor::fillKinematics
() {
113
ecm_tot
= (
target_mom
+
bullet_mom
).
m
();
114
115
scm_direction
=
scm_momentum
.
vect
().
unit
();
116
valong
=
velocity
.
dot
(
scm_direction
);
117
118
v2
=
velocity
.
mag2
();
119
120
G4double
pvsq =
v2
-
valong
*
valong
;
// velocity perp to scm_momentum
121
if
(
verboseLevel
> 3)
G4cout
<<
" pvsq "
<< pvsq <<
G4endl
;
122
123
degenerated
= (pvsq <
small
);
124
if
(
degenerated
&&
verboseLevel
> 2)
125
G4cout
<<
" degenerated case (already along Z) "
<<
G4endl
;
126
127
if
(
verboseLevel
> 3) {
128
G4cout
<<
" v2 "
<<
v2
<<
" valong "
<< valong
129
<<
" valong*valong "
<< valong*valong <<
G4endl
;
130
}
131
}
132
133
G4LorentzVector
134
G4LorentzConvertor::backToTheLab
(
const
G4LorentzVector
&
mom
)
const
{
135
if
(
verboseLevel
> 2)
136
G4cout
<<
" >>> G4LorentzConvertor::backToTheLab"
<<
G4endl
;
137
138
if
(
verboseLevel
> 3)
139
G4cout
<<
" at rest: px "
<< mom.
x
() <<
" py "
<< mom.
y
() <<
" pz "
140
<< mom.
z
() <<
" e "
<< mom.
e
() << G4endl
141
<<
" v2 "
<<
v2
<<
G4endl
;
142
143
G4LorentzVector
mom1 =
mom
;
144
if
(
v2
>
small
) mom1.
boost
(
velocity
);
145
146
if
(
verboseLevel
> 3)
147
G4cout
<<
" at lab: px "
<< mom1.
x
() <<
" py "
<< mom1.
y
() <<
" pz "
148
<< mom1.
z
() <<
G4endl
;
149
150
return
mom1;
151
}
152
153
154
// Bullet kinematics in target rest frame (LAB frame, usually)
155
156
G4double
G4LorentzConvertor::getKinEnergyInTheTRS
()
const
{
157
if
(
verboseLevel
> 2)
158
G4cout
<<
" >>> G4LorentzConvertor::getKinEnergyInTheTRS"
<<
G4endl
;
159
160
G4LorentzVector
bmom =
bullet_mom
;
161
bmom.
boost
(-
target_mom
.
boostVector
());
162
return
bmom.
e
()-bmom.
m
();
163
}
164
165
G4double
G4LorentzConvertor::getTRSMomentum
()
const
{
166
if
(
verboseLevel
> 2)
167
G4cout
<<
" >>> G4LorentzConvertor::getTRSMomentum"
<<
G4endl
;
168
169
G4LorentzVector
bmom =
bullet_mom
;
170
bmom.
boost
(-
target_mom
.
boostVector
());
171
return
bmom.
rho
();
172
}
173
174
G4LorentzVector
G4LorentzConvertor::rotate
(
const
G4LorentzVector
&
mom
)
const
{
175
if
(
verboseLevel
> 2)
176
G4cout
<<
" >>> G4LorentzConvertor::rotate(G4LorentzVector)"
<<
G4endl
;
177
178
if
(
verboseLevel
> 3) {
179
G4cout
<<
" valong "
<<
valong
<<
" degenerated "
<<
degenerated
<< G4endl
180
<<
" before rotation: px "
<< mom.
x
() <<
" py "
<< mom.
y
()
181
<<
" pz "
<< mom.
z
() <<
G4endl
;
182
}
183
184
G4LorentzVector
mom_rot =
mom
;
185
if
(!
degenerated
) {
186
if
(
verboseLevel
> 2)
187
G4cout
<<
" rotating to align with reference z axis "
<<
G4endl
;
188
189
G4ThreeVector
vscm =
velocity
-
valong
*
scm_direction
;
190
G4ThreeVector
vxcm = scm_direction.
cross
(
velocity
);
191
192
if
(vscm.
mag
() >
small
&& vxcm.
mag
() >
small
) {
// Double check
193
if
(
verboseLevel
> 3) {
194
G4cout
<<
" reference z axis "
<< scm_direction
195
<<
" vscm "
<< vscm <<
" vxcm "
<< vxcm <<
G4endl
;
196
}
197
198
mom_rot.
setVect
(mom.
x
()*vscm.
unit
() + mom.
y
()*vxcm.
unit
() +
199
mom.
z
()*
scm_direction
);
200
}
else
{
201
if
(
verboseLevel
)
202
G4cerr
<<
">>> G4LorentzVector::rotate zero with !degenerated"
<<
G4endl
;
203
}
204
}
205
206
if
(
verboseLevel
> 3) {
207
G4cout
<<
" after rotation: px "
<< mom_rot.
x
() <<
" py "
<< mom_rot.
y
()
208
<<
" pz "
<< mom_rot.
z
() <<
G4endl
;
209
}
210
211
return
mom_rot;
212
}
213
214
G4LorentzVector
G4LorentzConvertor::rotate
(
const
G4LorentzVector
& mom1,
215
const
G4LorentzVector
&
mom
)
const
{
216
if
(
verboseLevel
> 2)
217
G4cout
<<
" >>> G4LorentzConvertor::rotate(G4LorentzVector,G4LorentzVector)"
218
<<
G4endl
;
219
220
if
(
verboseLevel
> 3) {
221
G4cout
<<
" before rotation: px "
<< mom.
x
() <<
" py "
<< mom.
y
()
222
<<
" pz "
<< mom.
z
() <<
G4endl
;
223
}
224
225
G4ThreeVector
mom1_dir = mom1.
vect
().
unit
();
226
G4double
pv =
velocity
.
dot
(mom1_dir);
227
228
G4double
vperp =
v2
- pv*pv;
// velocity perpendicular to mom1
229
if
(
verboseLevel
> 3) {
230
G4cout
<<
" vperp "
<< vperp <<
" small? "
<< (vperp <=
small
) << G4endl;
231
}
232
233
G4LorentzVector
mom_rot =
mom
;
234
235
if
(vperp >
small
) {
236
if
(
verboseLevel
> 2)
237
G4cout
<<
" rotating to align with first z axis "
<<
G4endl
;
238
239
G4ThreeVector
vmom1 =
velocity
- pv*mom1_dir;
240
G4ThreeVector
vxm1 = mom1_dir.
cross
(
velocity
);
241
242
if
(vmom1.
mag
() >
small
&& vxm1.
mag
() >
small
) {
// Double check
243
if
(
verboseLevel
> 3) {
244
G4cout
<<
" first z axis "
<< mom1_dir << G4endl
245
<<
" vmom1 "
<< vmom1 <<
" vxm1 "
<< vxm1 <<
G4endl
;
246
}
247
248
mom_rot.
setVect
(mom.
x
()*vmom1.
unit
() + mom.
y
()*vxm1.
unit
() +
249
mom.
z
()*mom1_dir );
250
}
else
{
251
if
(
verboseLevel
)
252
G4cerr
<<
">>> G4LorentzVector::rotate zero with !degenerated"
<<
G4endl
;
253
}
254
}
255
256
if
(
verboseLevel
> 3) {
257
G4cout
<<
" after rotation: px "
<< mom_rot.
x
() <<
" py "
<< mom_rot.
y
()
258
<<
" pz "
<< mom_rot.
z
() <<
G4endl
;
259
}
260
261
return
mom_rot;
262
}
263
264
G4bool
G4LorentzConvertor::reflectionNeeded
()
const
{
265
if
(
verboseLevel
> 2)
266
G4cout
<<
" >>> G4LorentzConvertor::reflectionNeeded (query)"
<<
G4endl
;
267
268
if
(
verboseLevel
> 3) {
269
G4cout
<<
" v2 = "
<<
v2
<<
" SCM z = "
<<
scm_momentum
.
z
()
270
<<
" degenerated? "
<<
degenerated
<<
G4endl
;
271
}
272
273
if
(
v2
<
small
&& !
degenerated
)
274
throw
G4HadronicException
(__FILE__, __LINE__,
"G4LorentzConvertor::reflectionNeeded - return value undefined"
);
275
276
if
(
verboseLevel
> 2) {
277
G4cout
<<
" reflection across XY is"
278
<< ((
v2
>=
small
&& (!
degenerated
||
scm_momentum
.
z
()<0.0))?
""
:
" NOT"
)
279
<<
" needed"
<<
G4endl
;
280
}
281
282
return
(
v2
>=
small
&& (!
degenerated
||
scm_momentum
.
z
()<0.0));
283
}
284
285
286
// Reporting functions for diagnostics
287
288
void
G4LorentzConvertor::printBullet
()
const
{
289
G4cout
<<
" G4LC bullet: px "
<<
bullet_mom
.
px
() <<
" py "
<<
bullet_mom
.
py
()
290
<<
" pz "
<<
bullet_mom
.
pz
() <<
" e "
<<
bullet_mom
.
e
()
291
<<
" mass "
<<
bullet_mom
.
m
() <<
G4endl
;
292
}
293
294
void
G4LorentzConvertor::printTarget
()
const
{
295
G4cout
<<
" G4LC target: px "
<<
target_mom
.
px
() <<
" py "
<<
target_mom
.
py
()
296
<<
" pz "
<<
target_mom
.
pz
() <<
" e "
<<
target_mom
.
e
()
297
<<
" mass "
<<
target_mom
.
m
() <<
G4endl
;
298
}
299
geant4
tree
geant4-10.6-release
source
processes
hadronic
models
cascade
cascade
src
G4LorentzConvertor.cc
Built by
Jin Huang
. updated:
Wed Jun 29 2022 17:25:40
using
1.8.2 with
ECCE GitHub integration