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