Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BigBanger.cc
Go to the documentation of this file.
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// 20100114 M. Kelsey -- Remove G4CascadeMomentum, use G4LorentzVector directly
28// 20100301 M. Kelsey -- In generateBangInSCM(), restore old G4CascMom calcs.
29// for (N-1)th outgoing nucleon.
30// 20100319 M. Kelsey -- Use new generateWithRandomAngles for theta,phi stuff
31// 20100407 M. Kelsey -- Replace std::vector<> returns with data members.
32// 20100413 M. Kelsey -- Pass G4CollisionOutput by ref to ::collide()
33// 20100517 M. Kelsey -- Inherit from common base class, clean up code
34// 20100628 M. Kelsey -- Use old "bindingEnergy" fn as wrapper, add balance
35// checking after bang.
36// 20100630 M. Kelsey -- Just do simple boost for target, instead of using
37// G4LorentzConverter with dummy bullet.
38// 20100701 M. Kelsey -- Re-throw momentum list, not just angles!
39// 20100714 M. Kelsey -- Move conservation checking to base class
40// 20100726 M. Kelsey -- Move std::vector<> buffer to .hh file
41// 20100923 M. Kelsey -- Migrate to integer A and Z
42// 20110214 M. Kelsey -- Follow G4InuclParticle::Model enumerator migration
43// 20110806 M. Kelsey -- Pre-allocate buffers to reduce memory churn
44// 20110922 M. Kelsey -- Follow G4InuclParticle::print(ostream&) migration
45// 20120608 M. Kelsey -- Fix variable-name "shadowing" compiler warnings.
46// 20130622 Inherit from G4CascadeDeexciteBase, move to deExcite() interface
47// with G4Fragment
48// 20130924 M. Kelsey -- Replace std::pow with G4Pow::powN() for CPU speed
49// 20150608 M. Kelsey -- Label all while loops as terminating.
50
51#include <algorithm>
52
53#include "G4BigBanger.hh"
54#include "G4SystemOfUnits.hh"
55#include "G4CollisionOutput.hh"
56#include "G4InuclNuclei.hh"
60#include "G4Pow.hh"
61#include "Randomize.hh"
62
63using namespace G4InuclSpecialFunctions;
64
65typedef std::vector<G4InuclElementaryParticle>::iterator particleIterator;
66
68
70 G4CollisionOutput& globalOutput) {
71 if (verboseLevel) G4cout << " >>> G4BigBanger::deExcite" << G4endl;
72
73 getTargetData(target);
74 G4ThreeVector toTheLabFrame = PEX.boostVector(); // From rest to lab
75
76 // This "should" be difference between E-target and sum of m(nucleons)
77 G4double etot = (EEXS - bindingEnergy(A,Z)) * MeV/GeV; // To Bertini units
78 if (etot < 0.0) etot = 0.0;
79
80 if (verboseLevel > 2) {
81 G4cout << " BigBanger: target\n" << target
82 << "\n etot " << etot << G4endl;
83 }
84
85 if (verboseLevel > 3) {
86 G4LorentzVector PEXrest = PEX;
87 PEXrest.boost(-toTheLabFrame);
88 G4cout << " target rest frame: px " << PEXrest.px() << " py "
89 << PEXrest.py() << " pz " << PEXrest.pz() << " E " << PEXrest.e()
90 << G4endl;
91 }
92
93 generateBangInSCM(etot, A, Z);
94
95 if (verboseLevel > 2) {
96 G4cout << " particles " << particles.size() << G4endl;
97 for(G4int i = 0; i < G4int(particles.size()); i++)
98 G4cout << particles[i] << G4endl;
99 }
100
101 if (particles.empty()) { // No bang! Don't know why...
102 G4cerr << " >>> G4BigBanger unable to process fragment "
103 << target << G4endl;
104
105 // FIXME: This will violate baryon number, momentum, energy, etc.
106 return;
107 }
108
109 // convert back to Lab
110 G4LorentzVector totscm;
111 G4LorentzVector totlab;
112
113 if (verboseLevel > 2) G4cout << " BigBanger: boosting to lab" << G4endl;
114
115 particleIterator ipart;
116 for(ipart = particles.begin(); ipart != particles.end(); ipart++) {
117 G4LorentzVector mom = ipart->getMomentum();
118 if (verboseLevel > 2) totscm += mom;
119
120 mom.boost(toTheLabFrame);
121 if (verboseLevel > 2) totlab += mom;
122
123 ipart->setMomentum(mom);
124 if (verboseLevel > 2) G4cout << *ipart << G4endl;
125 }
126
127 std::sort(particles.begin(), particles.end(), G4ParticleLargerEkin());
128
129 validateOutput(target, particles); // Checks <vector> directly
130
131 if (verboseLevel > 2) {
132 G4cout << " In SCM: total outgoing momentum " << G4endl
133 << " E " << totscm.e() << " px " << totscm.x()
134 << " py " << totscm.y() << " pz " << totscm.z() << G4endl;
135 G4cout << " In Lab: mom cons " << G4endl
136 << " E " << PEX.e() - totlab.e() // PEX now includes EEXS
137 << " px " << PEX.x() - totlab.x()
138 << " py " << PEX.y() - totlab.y()
139 << " pz " << PEX.z() - totlab.z() << G4endl;
140 }
141
142 globalOutput.addOutgoingParticles(particles);
143}
144
145void G4BigBanger::generateBangInSCM(G4double etot, G4int a, G4int z) {
146 if (verboseLevel > 3) {
147 G4cout << " >>> G4BigBanger::generateBangInSCM" << G4endl;
148 }
149
150 const G4double ang_cut = 0.9999;
151 const G4int itry_max = 1000;
152
153 if (verboseLevel > 2) {
154 G4cout << " a " << a << " z " << z << G4endl;
155 }
156
157 particles.clear(); // Reset output vector before filling
158
159 if (a == 1) { // Special -- bare nucleon doesn't really "explode"
160 G4int knd = (z>0) ? 1 : 2;
161 particles.push_back(G4InuclElementaryParticle(knd)); // zero momentum
162 return;
163 }
164
165 // NOTE: If distribution fails, need to regenerate magnitudes and angles!
166 //*** generateMomentumModules(etot, a, z);
167
168 scm_momentums.reserve(a);
169 G4LorentzVector tot_mom;
170
171 G4bool bad = true;
172 G4int itry = 0;
173 while(bad && itry < itry_max) { /* Loop checking 08.06.2015 MHK */
174 itry++;
175 scm_momentums.clear();
176
177 generateMomentumModules(etot, a, z);
178 if (a == 2) {
179 // This is only a three-vector, not a four-vector
180 G4LorentzVector mom = generateWithRandomAngles(momModules[0]);
181 scm_momentums.push_back(mom);
182 scm_momentums.push_back(-mom); // Only safe since three-vector!
183 bad = false;
184 } else {
185 tot_mom *= 0.; // Easy way to reset accumulator
186
187 for(G4int i = 0; i < a-2; i++) { // All but last two are thrown
188 // This is only a three-vector, not a four-vector
189 G4LorentzVector mom = generateWithRandomAngles(momModules[i]);
190 scm_momentums.push_back(mom);
191 tot_mom += mom;
192 };
193
194 // handle last two
195 G4double tot_mod = tot_mom.rho();
196 G4double ct = -0.5*(tot_mod*tot_mod + momModules[a-2]*momModules[a-2]
197 - momModules[a-1]*momModules[a-1]) / tot_mod
198 / momModules[a-2];
199
200 if (verboseLevel > 2) G4cout << " ct last " << ct << G4endl;
201
202 if(std::fabs(ct) < ang_cut) {
203 // This is only a three-vector, not a four-vector
204 G4LorentzVector mom2 = generateWithFixedTheta(ct, momModules[a - 2]);
205
206 // rotate to the normal system
207 G4LorentzVector apr = tot_mom/tot_mod;
208 G4double a_tr = std::sqrt(apr.x()*apr.x() + apr.y()*apr.y());
209 G4LorentzVector mom;
210 mom.setX(mom2.z()*apr.x() + ( mom2.x()*apr.y() + mom2.y()*apr.z()*apr.x())/a_tr);
211 mom.setY(mom2.z()*apr.y() + (-mom2.x()*apr.x() + mom2.y()*apr.z()*apr.y())/a_tr);
212 mom.setZ(mom2.z()*apr.z() - mom2.y()*a_tr);
213
214 scm_momentums.push_back(mom);
215
216 // and the last one (again, not actually a four-vector!)
217 G4LorentzVector mom1 = -mom - tot_mom;
218
219 scm_momentums.push_back(mom1);
220 bad = false;
221 } // if (abs(ct) < ang_cut)
222 } // (a > 2)
223 } // while (bad && itry<itry_max)
224
225 if (!bad) {
226 particles.resize(a); // Use assignment to avoid temporaries
227 for(G4int i = 0; i < a; i++) {
228 G4int knd = i < z ? 1 : 2;
229
230 // Set to 0.0 the 4-th component (total energy) of the Lorentz momentum scm_momentums[i]
231 // in order to avoid very rare cases of unphysical negative (total and kinetic) energy
232 // (as reported by ATLAS with Geant4 10.6, but never reproduced in our tests).
233 // Note that these 4-vectors are actually 3-vectors, with null 4-th component in nearly
234 // all cases. After calling the method G4InuclElementaryParticle::fill (see below),
235 // the 4-th component of the momentum is set (in the method G4InuclParticle::setMomentum
236 // which in turn calls G4DynamicParticle::SetMomentum ) to the square root of the sum of
237 // the square of the mass and the square of the magnitude of the 3-momentum.
238 scm_momentums[i].setE( 0.0 );
239
240 particles[i].fill(scm_momentums[i], knd, G4InuclParticle::BigBanger);
241 };
242 };
243
244 if (verboseLevel > 2) {
245 if (itry == itry_max) G4cout << " BigBanger -> can not generate bang " << G4endl;
246 }
247
248 return;
249}
250
251void G4BigBanger::generateMomentumModules(G4double etot, G4int a, G4int z) {
252 if (verboseLevel > 3) {
253 G4cout << " >>> G4BigBanger::generateMomentumModules" << G4endl;
254 }
255
256 // Proton and neutron masses
259
260 momModules.clear(); // Reset buffer for filling
261
262 G4double xtot = 0.0;
263
264 if (a > 2) { // For "large" nuclei, energy is distributed
265 G4double promax = maxProbability(a);
266
267 momModules.resize(a, 0.); // Pre-allocate to avoid memory churn
268 for(G4int i = 0; i < a; i++) {
269 momModules[i] = generateX(a, promax);
270 xtot += momModules[i];
271
272 if (verboseLevel > 2) {
273 G4cout << " i " << i << " x " << momModules[i] << G4endl;
274 }
275 }
276 } else { // Two-body case is special, must be 50%
277 xtot = 1.;
278 momModules.push_back(0.5);
279 momModules.push_back(0.5);
280 }
281
282 for(G4int i = 0; i < a; i++) {
283 G4double mass = i < z ? mp : mn;
284
285 momModules[i] *= etot/xtot;
286 momModules[i] = std::sqrt(momModules[i] * (momModules[i] + 2.0 * mass));
287
288 if (verboseLevel > 2) {
289 G4cout << " i " << i << " pmod " << momModules[i] << G4endl;
290 }
291 };
292
293 return;
294}
295
296G4double G4BigBanger::xProbability(G4double x, G4int a) const {
297 if (verboseLevel > 3) G4cout << " >>> G4BigBanger::xProbability" << G4endl;
298
299 G4Pow* theG4Pow = G4Pow::GetInstance(); // For convenience
300
301 G4double ekpr = 0.0;
302 if(x < 1.0 || x > 0.0) {
303 ekpr = x * x;
304
305 if (a%2 == 0) { // even A
306 ekpr *= std::sqrt(1.0 - x) * theG4Pow->powN((1.0 - x), (3*a-6)/2);
307 }
308 else {
309 ekpr *= theG4Pow->powN((1.0 - x), (3*a-5)/2);
310 };
311 };
312
313 return ekpr;
314}
315
316G4double G4BigBanger::maxProbability(G4int a) const {
317 if (verboseLevel > 3) {
318 G4cout << " >>> G4BigBanger::maxProbability" << G4endl;
319 }
320
321 return xProbability(2./3./(a-1.0), a);
322}
323
324G4double G4BigBanger::generateX(G4int a, G4double promax) const {
325 if (verboseLevel > 3) G4cout << " >>> G4BigBanger::generateX" << G4endl;
326
327 const G4int itry_max = 1000;
328 G4int itry = 0;
329 G4double x;
330
331 while(itry < itry_max) { /* Loop checking 08.06.2015 MHK */
332 itry++;
333 x = G4UniformRand();
334 if(xProbability(x, a) >= promax*G4UniformRand() ) return x;
335 }
336
337 if (verboseLevel > 2) {
338 G4cout << " BigBanger -> can not generate x " << G4endl;
339 }
340
341 return maxProbability(a);
342}
std::vector< G4InuclElementaryParticle >::iterator particleIterator
CLHEP::HepLorentzVector G4LorentzVector
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
HepLorentzVector & boost(double, double, double)
virtual void deExcite(const G4Fragment &target, G4CollisionOutput &output)
G4CascadeDeexciteBase(const char *name)
void getTargetData(const G4Fragment &target)
virtual G4bool validateOutput(const G4Fragment &target, G4CollisionOutput &output)
void addOutgoingParticles(const std::vector< G4InuclElementaryParticle > &particles)
static G4double getParticleMass(G4int type)
static G4Pow * GetInstance()
Definition G4Pow.cc:41
G4double powN(G4double x, G4int n) const
Definition G4Pow.cc:162
G4LorentzVector generateWithFixedTheta(G4double ct, G4double p, G4double mass=0.)
G4double bindingEnergy(G4int A, G4int Z)
G4LorentzVector generateWithRandomAngles(G4double p, G4double mass=0.)