Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NucleiModel.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// 20100112 M. Kelsey -- Remove G4CascadeMomentum, use G4LorentzVector directly
29// 20100114 M. Kelsey -- Use G4ThreeVector for position
30// 20100309 M. Kelsey -- Use new generateWithRandomAngles for theta,phi stuff;
31// eliminate some unnecessary std::pow(), move ctor here
32// 20100407 M. Kelsey -- Replace std::vector<>::resize(0) with ::clear().
33// Move output vectors from generateParticleFate() and
34// ::generateInteractionPartners() to be data members; return via
35// const-ref instead of by value.
36// 20100413 M. Kelsey -- Pass G4CollisionOutput by ref to ::collide()
37// 20100418 M. Kelsey -- Reference output particle lists via const-ref
38// 20100421 M. Kelsey -- Replace hardwired p/n masses with G4PartDef's
39// 20100517 M. Kelsey -- Use G4CascadeINterpolator for cross-section
40// calculations. use G4Cascade*Channel for total xsec in pi-N
41// N-N channels. Move absorptionCrossSection() from SpecialFunc.
42// 20100610 M. Kelsey -- Replace another random-angle code block; add some
43// diagnostic output for partner-list production.
44// 20100617 M. Kelsey -- Replace preprocessor flag CHC_CHECK with
45// G4CASCADE_DEBUG_CHARGE
46// 20100620 M. Kelsey -- Improve error message on empty partners list, add
47// four-momentum checking after EPCollider
48// 20100621 M. Kelsey -- In boundaryTransition() account for momentum transfer
49// to secondary by boosting into recoil nucleus "rest" frame.
50// Don't need temporaries to create dummy "partners" for list.
51// 20100622 M. Kelsey -- Restore use of "bindingEnergy()" function name, which
52// is now a wrapper for G4NucleiProperties::GetBindingEnergy().
53// 20100623 M. Kelsey -- Eliminate some temporaries terminating partner-list,
54// discard recoil boost for now. Initialize all data
55// members in ctors. Allow generateModel() to be called
56// mutliple times, by clearing vectors each time through;
57// avoid extra work by returning if A and Z are same as
58// before.
59// 20100628 M. Kelsey -- Two momentum-recoil bugs; don't subtract energies!
60// 20100715 M. Kelsey -- Make G4InuclNuclei in generateModel(), use for
61// balance checking (only verbose>2) in generateParticleFate().
62// 20100721 M. Kelsey -- Use new G4CASCADE_CHECK_ECONS for conservation checks
63// 20100723 M. Kelsey -- Move G4CollisionOutput buffer to .hh for reuse
64// 20100726 M. Kelsey -- Preallocate arrays with number_of_zones dimension.
65// 20100902 M. Kelsey -- Remove resize(3) directives from qdeutron/acsecs
66// 20100907 M. Kelsey -- Limit interaction targets based on current nucleon
67// counts (e.g., no pp if protonNumberCurrent < 2!)
68// 20100914 M. Kelsey -- Migrate to integer A and Z
69// 20100928 M. Kelsey -- worthToPropagate() use nuclear potential for hadrons.
70// 20101005 M. Kelsey -- Annotate hardwired numbers for upcoming rationalizing,
71// move hardwired numbers out to static data members, factor out
72// all the model construction pieces to separate functions, fix
73// wrong value for 4/3 pi.
74// 20101019 M. Kelsey -- CoVerity reports: unitialized constructor, dtor leak
75// 20101020 M. Kelsey -- Bug fixes to refactoring changes (5 Oct). Back out
76// worthToPropagate() changes for better regression testing.
77// 20101020 M. Kelsey -- Re-activate worthToPropagate() changes.
78// 20101119 M. Kelsey -- Hide "negative path" and "no partners" messages in
79// verbosity.
80// 20110218 M. Kelsey -- Add crossSectionUnits and radiusUnits scale factors,
81// use "theoretical" numbers for radii etc., multipled by scale
82// factor; set scale factors using environment variables
83// 20110303 M. Kelsey -- Add comments why using fabs() with B.E. differences?
84// 20110321 M. Kelsey -- Replace strtof() with strtod() for envvar conversion
85// 20110321 M. Kelsey -- Use fm and fm^2 as default units, Per D. Wright
86// (NOTE: Restored from original 20110318 commit)
87// 20110324 D. Wright -- Implement trailing effect
88// 20110324 M. Kelsey -- Move ::reset() here, as it has more code.
89// 20110519 M. Kelsey -- Used "rho" after assignment, instead of recomputing
90// 20110525 M. Kelsey -- Revert scale factor changes (undo 20110321 changes)
91// 20110617 M. Kelsey -- Apply scale factor to trailing-effect radius, make
92// latter runtime adjustable (G4NUCMODEL_RAD_TRAILING)
93// 20110720 M. Kelsey -- Follow interface change for cross-section tables,
94// eliminating switch blocks.
95// 20110806 M. Kelsey -- Reduce memory churn by pre-allocating buffers
96// 20110823 M. Kelsey -- Remove local cross-section tables entirely
97// 20110825 M. Kelsey -- Add comments regarding Fermi momentum scale, set of
98// "best guess" parameter values
99// 20110831 M. Kelsey -- Make "best guess" parameters the defaults
100// 20110922 M. Kelsey -- Follow migrations G4InuclParticle::print(ostream&)
101// and G4CascadParticle::print(ostream&)
102// 20111018 M. Kelsey -- Correct kaon potential to be positive, not negative
103// 20111107 M. Kelsey -- *** REVERT TO OLD NON-PHYSICAL PARAMETERS FOR 9.5 ***
104// 20120306 M. Kelsey -- For incident projectile (CParticle incoming to
105// nucleus from outside) don't use pw cut; force interaction.
106// 20120307 M. Kelsey -- Compute zone volumes during setup, not during each
107// interaction. Include photons in absorption by dibaryons.
108// Discard unused code in totalCrossSection(). Encapsulate
109// interaction path calculations in generateInteractionLength().
110// 20120308 M. Kelsey -- Add binned photon absorption cross-section.
111// 20120320 M. Kelsey -- Bug fix: fill zone_volumes for single-zone case
112// 20120321 M. Kelsey -- Add check in isProjectile() for single-zone case.
113// 20120608 M. Kelsey -- Fix variable-name "shadowing" compiler warnings.
114// 20120822 M. Kelsey -- Move envvars to G4CascadeParameters.
115
116#include <numeric>
117
118#include "G4NucleiModel.hh"
119#include "G4PhysicalConstants.hh"
120#include "G4SystemOfUnits.hh"
121#include "G4CascadeChannel.hh"
125#include "G4CascadeParameters.hh"
126#include "G4CollisionOutput.hh"
128#include "G4HadTmpUtil.hh"
129#include "G4InuclNuclei.hh"
132#include "G4LorentzConvertor.hh"
133#include "G4Neutron.hh"
136#include "G4Proton.hh"
137
138using namespace G4InuclParticleNames;
139using namespace G4InuclSpecialFunctions;
140
141typedef std::vector<G4InuclElementaryParticle>::iterator particleIterator;
142
143// For the best approximation to a physical-units model, set the following:
144// setenv G4NUCMODEL_XSEC_SCALE 0.1
145// setenv G4NUCMODEL_RAD_SCALE 1.0
146// setenv G4NUCMODEL_RAD_2PAR 1
147// setenv G4NUCMODEL_RAD_SMALL 1.992
148// setenv G4NUCMODEL_RAD_ALPHA 0.84
149// setenv G4NUCMODEL_FERMI_SCALE 0.685
150// setenv G4NUCMODEL_RAD_TRAILING 0.70
151
152// Scaling factors for radii and cross-sections, currently different!
153const G4double G4NucleiModel::crossSectionUnits = G4CascadeParameters::xsecScale();
154const G4double G4NucleiModel::radiusUnits = G4CascadeParameters::radiusScale();
155const G4double G4NucleiModel::skinDepth = 0.611207*radiusUnits;
156
157// One- vs. two-parameter nuclear radius based on envvar
158// ==> radius = radiusScale*cbrt(A) + radiusScale2/cbrt(A)
159
160const G4double G4NucleiModel::radiusScale =
161 (G4CascadeParameters::useTwoParam() ? 1.16 : 1.2) * radiusUnits;
162const G4double G4NucleiModel::radiusScale2 =
163 (G4CascadeParameters::useTwoParam() ? -1.3456 : 0.) * radiusUnits;
164
165// NOTE: Old code used R_small = 8.0 (~2.83*units), and R_alpha = 0.7*R_small
166// Published data suggests R_small ~ 1.992 fm, R_alpha = 0.84*R_small
167
168const G4double G4NucleiModel::radiusForSmall = G4CascadeParameters::radiusSmall();
169
170const G4double G4NucleiModel::radScaleAlpha = G4CascadeParameters::radiusAlpha();
171
172// Scale factor relating Fermi momentum to density of states, units GeV.fm
173// NOTE: Old code has 0.685*units GeV.fm, literature suggests 0.470 GeV.fm,
174// but this gives too small momentum; old value gives <P_F> ~ 270 MeV
175const G4double G4NucleiModel::fermiMomentum = G4CascadeParameters::fermiScale();
176
177// Effective radius (0.87 to 1.2 fm) of nucleon, for trailing effect
178const G4double G4NucleiModel::R_nucleon = G4CascadeParameters::radiusTrailing();
179
180// Zone boundaries as fraction of nuclear radius (from outside in)
181const G4double G4NucleiModel::alfa3[3] = { 0.7, 0.3, 0.01 };
182const G4double G4NucleiModel::alfa6[6] = { 0.9, 0.6, 0.4, 0.2, 0.1, 0.05 };
183
184// Flat nuclear potentials for mesons and hyperons (GeV)
185const G4double G4NucleiModel::pion_vp = 0.007;
186const G4double G4NucleiModel::pion_vp_small = 0.007;
187const G4double G4NucleiModel::kaon_vp = 0.015;
188const G4double G4NucleiModel::hyperon_vp = 0.030;
189
190const G4double G4NucleiModel::piTimes4thirds = pi*4./3.;
191
192
193// Constructors
195 : verboseLevel(0), nuclei_radius(0.), nuclei_volume(0.), number_of_zones(0),
196 A(0), Z(0), theNucleus(0), neutronNumber(0), protonNumber(0),
197 neutronNumberCurrent(0), protonNumberCurrent(0), current_nucl1(0),
198 current_nucl2(0) {}
199
201 : verboseLevel(0), nuclei_radius(0.), nuclei_volume(0.), number_of_zones(0),
202 A(0), Z(0), theNucleus(0), neutronNumber(0), protonNumber(0),
203 neutronNumberCurrent(0), protonNumberCurrent(0), current_nucl1(0),
204 current_nucl2(0) {
205 generateModel(a,z);
206}
207
209 : verboseLevel(0), nuclei_radius(0.), nuclei_volume(0.), number_of_zones(0),
210 A(0), Z(0), theNucleus(0), neutronNumber(0), protonNumber(0),
211 neutronNumberCurrent(0), protonNumberCurrent(0), current_nucl1(0),
212 current_nucl2(0) {
214}
215
217 delete theNucleus;
218 theNucleus = 0;
219}
220
221
222// Initialize model state for new cascade
223
224void G4NucleiModel::reset(G4int nHitNeutrons, G4int nHitProtons,
225 const std::vector<G4ThreeVector>* hitPoints) {
226 neutronNumberCurrent = neutronNumber - nHitNeutrons;
227 protonNumberCurrent = protonNumber - nHitProtons;
228
229 // zero or copy collision point array for trailing effect
230 if (!hitPoints || !hitPoints->empty()) collisionPts.clear();
231 else collisionPts = *hitPoints;
232}
233
234
235// Generate nuclear model parameters for given nucleus
236
238 generateModel(nuclei->getA(), nuclei->getZ());
239}
240
242 if (verboseLevel) {
243 G4cout << " >>> G4NucleiModel::generateModel A " << a << " Z " << z
244 << G4endl;
245 }
246
247 // If model already built, just return; otherwise intialize everything
248 if (a == A && z == Z) {
249 if (verboseLevel > 1) G4cout << " model already generated" << z << G4endl;
250 reset();
251 return;
252 }
253
254 A = a;
255 Z = z;
256 delete theNucleus;
257 theNucleus = new G4InuclNuclei(A,Z); // For conservation checking
258
259 neutronNumber = A - Z;
260 protonNumber = Z;
261 reset();
262
263 if (verboseLevel > 3) {
264 G4cout << " crossSectionUnits = " << crossSectionUnits << G4endl
265 << " radiusUnits = " << radiusUnits << G4endl
266 << " skinDepth = " << skinDepth << G4endl
267 << " radiusScale = " << radiusScale << G4endl
268 << " radiusScale2 = " << radiusScale2 << G4endl
269 << " radiusForSmall = " << radiusForSmall << G4endl
270 << " radScaleAlpha = " << radScaleAlpha << G4endl
271 << " fermiMomentum = " << fermiMomentum << G4endl
272 << " piTimes4thirds = " << piTimes4thirds << G4endl;
273 }
274
275 G4double nuclearRadius; // Nuclear radius computed from A
276 if (A>4) nuclearRadius = radiusScale*G4cbrt(A) + radiusScale2/G4cbrt(A);
277 else nuclearRadius = radiusForSmall * (A==4 ? radScaleAlpha : 1.);
278
279 // This will be used to pre-allocate lots of arrays below
280 number_of_zones = (A < 5) ? 1 : (A < 100) ? 3 : 6;
281
282 // Clear all parameters arrays for reloading
283 binding_energies.clear();
284 nucleon_densities.clear();
285 zone_potentials.clear();
286 fermi_momenta.clear();
287 zone_radii.clear();
288 zone_volumes.clear();
289
290 fillBindingEnergies();
291 fillZoneRadii(nuclearRadius);
292
293 G4double tot_vol = fillZoneVolumes(nuclearRadius); // Woods-Saxon integral
294
295 fillPotentials(proton, tot_vol); // Protons
296 fillPotentials(neutron, tot_vol); // Neutrons
297
298 // Additional flat zone potentials for other hadrons
299 const std::vector<G4double> vp(number_of_zones, (A>4)?pion_vp:pion_vp_small);
300 const std::vector<G4double> kp(number_of_zones, kaon_vp);
301 const std::vector<G4double> hp(number_of_zones, hyperon_vp);
302
303 zone_potentials.push_back(vp);
304 zone_potentials.push_back(kp);
305 zone_potentials.push_back(hp);
306
307 nuclei_radius = zone_radii.back();
308 nuclei_volume = std::accumulate(zone_volumes.begin(),zone_volumes.end(),0.);
309
310 if (verboseLevel > 3) printModel();
311}
312
313
314// Load binding energy array for current nucleus
315
316void G4NucleiModel::fillBindingEnergies() {
317 if (verboseLevel > 1)
318 G4cout << " >>> G4NucleiModel::fillBindingEnergies" << G4endl;
319
320 G4double dm = bindingEnergy(A,Z);
321
322 // Binding energy differences for proton and neutron loss, respectively
323 // FIXME: Why is fabs() used here instead of the signed difference?
324 binding_energies.push_back(std::fabs(bindingEnergy(A-1,Z-1)-dm)/GeV);
325 binding_energies.push_back(std::fabs(bindingEnergy(A-1,Z)-dm)/GeV);
326}
327
328// Load zone boundary radius array for current nucleus
329
330void G4NucleiModel::fillZoneRadii(G4double nuclearRadius) {
331 if (verboseLevel > 1)
332 G4cout << " >>> G4NucleiModel::fillZoneRadii" << G4endl;
333
334 G4double skinRatio = nuclearRadius/skinDepth;
335 G4double skinDecay = std::exp(-skinRatio);
336
337 if (A < 5) { // Light ions treated as simple balls
338 zone_radii.push_back(nuclearRadius);
339 ur[0] = 0.;
340 ur[1] = 1.;
341 } else if (A < 12) { // Small nuclei have Gaussian potential
342 G4double rSq = nuclearRadius * nuclearRadius;
343 G4double gaussRadius = std::sqrt(rSq * (1.0 - 1.0/A) + 6.4);
344
345 ur[0] = 0.0;
346 for (G4int i = 0; i < number_of_zones; i++) {
347 G4double y = std::sqrt(-std::log(alfa3[i]));
348 zone_radii.push_back(gaussRadius * y);
349 ur[i+1] = y;
350 }
351 } else if (A < 100) { // Intermediate nuclei have three-zone W-S
352 ur[0] = -skinRatio;
353 for (G4int i = 0; i < number_of_zones; i++) {
354 G4double y = std::log((1.0 + skinDecay)/alfa3[i] - 1.0);
355 zone_radii.push_back(nuclearRadius + skinDepth * y);
356 ur[i+1] = y;
357 }
358 } else { // Heavy nuclei have six-zone W-S
359 ur[0] = -skinRatio;
360 for (G4int i = 0; i < number_of_zones; i++) {
361 G4double y = std::log((1.0 + skinDecay)/alfa6[i] - 1.0);
362 zone_radii.push_back(nuclearRadius + skinDepth * y);
363 ur[i+1] = y;
364 }
365 }
366}
367
368// Compute zone-by-zone density integrals
369
370G4double G4NucleiModel::fillZoneVolumes(G4double nuclearRadius) {
371 if (verboseLevel > 1)
372 G4cout << " >>> G4NucleiModel::fillZoneVolumes" << G4endl;
373
374 G4double tot_vol = 0.; // Return value omitting 4pi/3 for sphere!
375
376 if (A < 5) { // Light ions are treated as simple balls
377 v[0] = v1[0] = 1.;
378 tot_vol = zone_radii[0]*zone_radii[0]*zone_radii[0];
379 zone_volumes.push_back(tot_vol*piTimes4thirds); // True volume of sphere
380
381 return tot_vol;
382 }
383
384 PotentialType usePotential = (A < 12) ? Gaussian : WoodsSaxon;
385
386 for (G4int i = 0; i < number_of_zones; i++) {
387 if (usePotential == WoodsSaxon) {
388 v[i] = zoneIntegralWoodsSaxon(ur[i], ur[i+1], nuclearRadius);
389 } else {
390 v[i] = zoneIntegralGaussian(ur[i], ur[i+1], nuclearRadius);
391 }
392
393 tot_vol += v[i];
394
395 v1[i] = zone_radii[i]*zone_radii[i]*zone_radii[i];
396 if (i > 0) v1[i] -= zone_radii[i-1]*zone_radii[i-1]*zone_radii[i-1];
397
398 zone_volumes.push_back(v1[i]*piTimes4thirds); // True volume of shell
399 }
400
401 return tot_vol; // Sum of zone integrals, not geometric volume
402}
403
404// Load potentials for nucleons (using G4InuclParticle codes)
405void G4NucleiModel::fillPotentials(G4int type, G4double tot_vol) {
406 if (verboseLevel > 1)
407 G4cout << " >>> G4NucleiModel::fillZoneVolumes(" << type << ")" << G4endl;
408
409 if (type != proton && type != neutron) return;
410
412
413 // FIXME: This is the fabs() binding energy difference, not signed
414 const G4double dm = binding_energies[type-1];
415
416 rod.clear(); rod.reserve(number_of_zones);
417 pf.clear(); pf.reserve(number_of_zones);
418 vz.clear(); vz.reserve(number_of_zones);
419
420 G4int nNucleons = (type==proton) ? protonNumber : neutronNumber;
421 G4double dd0 = nNucleons / tot_vol / piTimes4thirds;
422
423 for (G4int i = 0; i < number_of_zones; i++) {
424 G4double rd = dd0 * v[i] / v1[i];
425 rod.push_back(rd);
426 G4double pff = fermiMomentum * G4cbrt(rd);
427 pf.push_back(pff);
428 vz.push_back(0.5 * pff * pff / mass + dm);
429 }
430
431 nucleon_densities.push_back(rod);
432 fermi_momenta.push_back(pf);
433 zone_potentials.push_back(vz);
434}
435
436// Zone integral of Woods-Saxon density function
437G4double G4NucleiModel::zoneIntegralWoodsSaxon(G4double r1, G4double r2,
438 G4double nuclearRadius) const {
439 if (verboseLevel > 1) {
440 G4cout << " >>> G4NucleiModel::zoneIntegralWoodsSaxon" << G4endl;
441 }
442
443 const G4double epsilon = 1.0e-3;
444 const G4int itry_max = 1000;
445
446 G4double skinRatio = nuclearRadius / skinDepth;
447
448 G4double d2 = 2.0 * skinRatio;
449 G4double dr = r2 - r1;
450 G4double fr1 = r1 * (r1 + d2) / (1.0 + std::exp(r1));
451 G4double fr2 = r2 * (r2 + d2) / (1.0 + std::exp(r2));
452 G4double fi = (fr1 + fr2) / 2.;
453 G4double fun1 = fi * dr;
454 G4double fun;
455 G4double jc = 1;
456 G4double dr1 = dr;
457 G4int itry = 0;
458
459 while (itry < itry_max) {
460 dr /= 2.;
461 itry++;
462
463 G4double r = r1 - dr;
464 fi = 0.0;
465 G4int jc1 = G4int(std::pow(2.0, jc - 1) + 0.1);
466
467 for (G4int i = 0; i < jc1; i++) {
468 r += dr1;
469 fi += r * (r + d2) / (1.0 + std::exp(r));
470 }
471
472 fun = 0.5 * fun1 + fi * dr;
473
474 if (std::fabs((fun - fun1) / fun) <= epsilon) break;
475
476 jc++;
477 dr1 = dr;
478 fun1 = fun;
479 } // while (itry < itry_max)
480
481 if (verboseLevel > 2 && itry == itry_max)
482 G4cout << " zoneIntegralWoodsSaxon-> n iter " << itry_max << G4endl;
483
484 G4double skinDepth3 = skinDepth*skinDepth*skinDepth;
485
486 return skinDepth3 * (fun + skinRatio*skinRatio*std::log((1.0 + std::exp(-r1)) / (1.0 + std::exp(-r2))));
487}
488
489
490// Zone integral of Gaussian density function
491G4double G4NucleiModel::zoneIntegralGaussian(G4double r1, G4double r2,
492 G4double nucRad) const {
493 if (verboseLevel > 1) {
494 G4cout << " >>> G4NucleiModel::zoneIntegralGaussian" << G4endl;
495 }
496
497 G4double gaussRadius = std::sqrt(nucRad*nucRad * (1.0 - 1.0/A) + 6.4);
498
499 const G4double epsilon = 1.0e-3;
500 const G4int itry_max = 1000;
501
502 G4double dr = r2 - r1;
503 G4double fr1 = r1 * r1 * std::exp(-r1 * r1);
504 G4double fr2 = r2 * r2 * std::exp(-r2 * r2);
505 G4double fi = (fr1 + fr2) / 2.;
506 G4double fun1 = fi * dr;
507 G4double fun;
508 G4double jc = 1;
509 G4double dr1 = dr;
510 G4int itry = 0;
511
512 while (itry < itry_max) {
513 dr /= 2.;
514 itry++;
515 G4double r = r1 - dr;
516 fi = 0.0;
517 G4int jc1 = int(std::pow(2.0, jc - 1) + 0.1);
518
519 for (G4int i = 0; i < jc1; i++) {
520 r += dr1;
521 fi += r * r * std::exp(-r * r);
522 }
523
524 fun = 0.5 * fun1 + fi * dr;
525
526 if (std::fabs((fun - fun1) / fun) <= epsilon) break;
527
528 jc++;
529 dr1 = dr;
530 fun1 = fun;
531 } // while (itry < itry_max)
532
533 if (verboseLevel > 2 && itry == itry_max)
534 G4cerr << " zoneIntegralGaussian-> n iter " << itry_max << G4endl;
535
536 return gaussRadius*gaussRadius*gaussRadius * fun;
537}
538
539
541 if (verboseLevel > 1) {
542 G4cout << " >>> G4NucleiModel::printModel" << G4endl;
543 }
544
545 G4cout << " nuclei model for A " << A << " Z " << Z << G4endl
546 << " proton binding energy " << binding_energies[0]
547 << " neutron binding energy " << binding_energies[1] << G4endl
548 << " Nuclei radius " << nuclei_radius << " volume " << nuclei_volume
549 << " number of zones " << number_of_zones << G4endl;
550
551 for (G4int i = 0; i < number_of_zones; i++)
552 G4cout << " zone " << i+1 << " radius " << zone_radii[i]
553 << " volume " << zone_volumes[i] << G4endl
554 << " protons: density " << getDensity(1,i) << " PF " <<
555 getFermiMomentum(1,i) << " VP " << getPotential(1,i) << G4endl
556 << " neutrons: density " << getDensity(2,i) << " PF " <<
557 getFermiMomentum(2,i) << " VP " << getPotential(2,i) << G4endl
558 << " pions: VP " << getPotential(3,i) << G4endl;
559}
560
561
563 G4double ekin = 0.0;
564
565 if (ip < 3 && izone < number_of_zones) { // ip for proton/neutron only
566 G4double pfermi = fermi_momenta[ip - 1][izone];
568
569 ekin = std::sqrt(pfermi*pfermi + mass*mass) - mass;
570 }
571 return ekin;
572}
573
574
577 G4double pmod = getFermiMomentum(type, zone) * G4cbrt(inuclRndm());
579
580 return generateWithRandomAngles(pmod, mass);
581}
582
583
586 if (verboseLevel > 1) {
587 G4cout << " >>> G4NucleiModel::generateNucleon" << G4endl;
588 }
589
591 return G4InuclElementaryParticle(mom, type);
592}
593
594
596G4NucleiModel::generateQuasiDeutron(G4int type1, G4int type2,
597 G4int zone) const {
598 if (verboseLevel > 1) {
599 G4cout << " >>> G4NucleiModel::generateQuasiDeutron" << G4endl;
600 }
601
602 // Quasideuteron consists of an unbound but associated nucleon pair
603
604 // FIXME: Why generate two separate nucleon momenta (randomly!) and
605 // add them, instead of just throwing a net momentum for the
606 // dinulceon state? And why do I have to capture the two
607 // return values into local variables?
608 G4LorentzVector mom1 = generateNucleonMomentum(type1, zone);
609 G4LorentzVector mom2 = generateNucleonMomentum(type2, zone);
610 G4LorentzVector dmom = mom1+mom2;
611
612 G4int dtype = 0;
613 if (type1*type2 == pro*pro) dtype = 111;
614 else if (type1*type2 == pro*neu) dtype = 112;
615 else if (type1*type2 == neu*neu) dtype = 122;
616
617 return G4InuclElementaryParticle(dmom, dtype);
618}
619
620
621void
622G4NucleiModel::generateInteractionPartners(G4CascadParticle& cparticle) {
623 if (verboseLevel > 1) {
624 G4cout << " >>> G4NucleiModel::generateInteractionPartners" << G4endl;
625 }
626
627 const G4double small = 1.0e-10; // Cutoff for identifying 0.
628
629 thePartners.clear(); // Reset buffer for next cycle
630
631 G4int ptype = cparticle.getParticle().type();
632 G4int zone = cparticle.getCurrentZone();
633
634 G4double r_in; // Destination radius if inbound
635 G4double r_out; // Destination radius if outbound
636
637 if (zone == number_of_zones) {
638 r_in = nuclei_radius;
639 r_out = 0.0;
640 } else if (zone == 0) { // particle is outside core
641 r_in = 0.0;
642 r_out = zone_radii[0];
643 } else {
644 r_in = zone_radii[zone - 1];
645 r_out = zone_radii[zone];
646 }
647
648 G4double path = cparticle.getPathToTheNextZone(r_in, r_out);
649
650 if (verboseLevel > 2) {
651 if (isProjectile(cparticle)) G4cout << " incident particle" << G4endl;
652 G4cout << " r_in " << r_in << " r_out " << r_out << " path " << path
653 << G4endl;
654 }
655
656 if (path < -small) { // something wrong
657 if (verboseLevel)
658 G4cerr << " generateInteractionPartners-> negative path length" << G4endl;
659 return;
660 }
661
662 if (std::fabs(path) < small) { // just on the boundary
663 if (verboseLevel > 3)
664 G4cout << " generateInteractionPartners-> zero path" << G4endl;
665
666 thePartners.push_back(partner()); // Dummy list terminator with zero path
667 return;
668 }
669
670 // Initialize frame-transformations for current particle
671 dummy_convertor.setBullet(cparticle.getParticle());
672
673 G4double invmfp = 0.; // Buffers for interaction probability
674 G4double spath = 0.;
675 for (G4int ip = 1; ip < 3; ip++) {
676 // Only process nucleons which remain active in target
677 if (ip==proton && protonNumberCurrent < 1) continue;
678 if (ip==neutron && neutronNumberCurrent < 1) continue;
679
680 // All nucleons are assumed to be at rest when colliding
681 G4InuclElementaryParticle particle = generateNucleon(ip, zone);
682 invmfp = inverseMeanFreePath(cparticle, particle);
683 spath = generateInteractionLength(cparticle, path, invmfp);
684
685 if (spath < path) {
686 if (verboseLevel > 3) {
687 G4cout << " adding partner[" << thePartners.size() << "]: "
688 << particle << G4endl;
689 }
690 thePartners.push_back(partner(particle, spath));
691 }
692 } // for (G4int ip...
693
694 if (verboseLevel > 2) {
695 G4cout << " after nucleons " << thePartners.size() << " path " << path << G4endl;
696 }
697
698 // Absorption possible for pions or photons interacting with dibaryons
699 if (cparticle.getParticle().pion() || cparticle.getParticle().isPhoton()) {
700 if (verboseLevel > 2) {
701 G4cout << " trying quasi-deuterons with bullet: "
702 << cparticle.getParticle() << G4endl;
703 }
704
705 // Initialize buffers for quasi-deuteron results
706 qdeutrons.clear();
707 acsecs.clear();
708
709 G4double tot_invmfp = 0.0; // Total inv. mean-free-path for all QDs
710
711 // Proton-proton state interacts with pi- or neutrals
712 if (protonNumberCurrent >= 2 && ptype != pip) {
713 G4InuclElementaryParticle ppd = generateQuasiDeutron(pro, pro, zone);
714 if (verboseLevel > 2)
715 G4cout << " ptype=" << ptype << " using pp target\n" << ppd << G4endl;
716
717 invmfp = inverseMeanFreePath(cparticle, ppd);
718 tot_invmfp += invmfp;
719 acsecs.push_back(invmfp);
720 qdeutrons.push_back(ppd);
721 }
722
723 // Proton-neutron state interacts with any pion type or photon
724 if (protonNumberCurrent >= 1 && neutronNumberCurrent >= 1) {
725 G4InuclElementaryParticle npd = generateQuasiDeutron(pro, neu, zone);
726 if (verboseLevel > 2)
727 G4cout << " ptype=" << ptype << " using np target\n" << npd << G4endl;
728
729 invmfp = inverseMeanFreePath(cparticle, npd);
730 tot_invmfp += invmfp;
731 acsecs.push_back(invmfp);
732 qdeutrons.push_back(npd);
733 }
734
735 // Neutron-neutron state interacts with pi+ or neutrals
736 if (neutronNumberCurrent >= 2 && ptype != pim) {
737 G4InuclElementaryParticle nnd = generateQuasiDeutron(neu, neu, zone);
738 if (verboseLevel > 2)
739 G4cout << " ptype=" << ptype << " using nn target\n" << nnd << G4endl;
740
741 invmfp = inverseMeanFreePath(cparticle, nnd);
742 tot_invmfp += invmfp;
743 acsecs.push_back(invmfp);
744 qdeutrons.push_back(nnd);
745 }
746
747 // Select quasideuteron interaction from non-zero cross-section choices
748 if (verboseLevel > 2) {
749 for (size_t i=0; i<qdeutrons.size(); i++) {
750 G4cout << " acsecs[" << qdeutrons[i].getDefinition()->GetParticleName()
751 << "] " << acsecs[i];
752 }
753 G4cout << G4endl;
754 }
755
756 if (tot_invmfp > small) { // Must have absorption cross-section
757 G4double apath = generateInteractionLength(cparticle, path, tot_invmfp);
758
759 if (apath < path) { // choose the qdeutron
760 G4double sl = inuclRndm() * tot_invmfp;
761 G4double as = 0.0;
762
763 for (size_t i = 0; i < qdeutrons.size(); i++) {
764 as += acsecs[i];
765 if (sl < as) {
766 if (verboseLevel > 2) G4cout << " deut type " << i << G4endl;
767 thePartners.push_back(partner(qdeutrons[i], apath));
768 break;
769 }
770 } // for (qdeutrons...
771 } // if (apath < path)
772 } // if (tot_invmfp > small)
773 } // pion or photon
774
775 if (verboseLevel > 2) {
776 G4cout << " after deuterons " << thePartners.size() << " partners"
777 << G4endl;
778 }
779
780 if (thePartners.size() > 1) { // Sort list by path length
781 std::sort(thePartners.begin(), thePartners.end(), sortPartners);
782 }
783
784 G4InuclElementaryParticle particle; // Total path at end of list
785 thePartners.push_back(partner(particle, path));
786}
787
788
791 G4ElementaryParticleCollider* theEPCollider,
792 std::vector<G4CascadParticle>& outgoing_cparticles) {
793 if (verboseLevel > 1)
794 G4cout << " >>> G4NucleiModel::generateParticleFate" << G4endl;
795
796 if (verboseLevel > 2) G4cout << " cparticle: " << cparticle << G4endl;
797
798 // Create four-vector checking
799#ifdef G4CASCADE_CHECK_ECONS
800 G4CascadeCheckBalance balance(0.005, 0.01, "G4NucleiModel"); // Second arg is in GeV
801 balance.setVerboseLevel(verboseLevel);
802#endif
803
804 outgoing_cparticles.clear(); // Clear return buffer for this event
805 generateInteractionPartners(cparticle); // Fills "thePartners" data
806
807 if (thePartners.empty()) { // smth. is wrong -> needs special treatment
808 if (verboseLevel)
809 G4cerr << " generateParticleFate-> got empty interaction-partners list "
810 << G4endl;
811 return;
812 }
813
814 G4int npart = thePartners.size(); // Last item is a total-path placeholder
815
816 if (npart == 1) { // cparticle is on the next zone entry
817 cparticle.propagateAlongThePath(thePartners[0].second);
818 cparticle.incrementCurrentPath(thePartners[0].second);
819 boundaryTransition(cparticle);
820 outgoing_cparticles.push_back(cparticle);
821
822 if (verboseLevel > 2) G4cout << " next zone \n" << cparticle << G4endl;
823 } else { // there are possible interactions
824 if (verboseLevel > 1)
825 G4cout << " processing " << npart-1 << " possible interactions" << G4endl;
826
827 G4ThreeVector old_position = cparticle.getPosition();
828 G4InuclElementaryParticle& bullet = cparticle.getParticle();
829 G4bool no_interaction = true;
830 G4int zone = cparticle.getCurrentZone();
831
832 for (G4int i=0; i<npart-1; i++) { // Last item is a total-path placeholder
833 if (i > 0) cparticle.updatePosition(old_position);
834
835 G4InuclElementaryParticle& target = thePartners[i].first;
836
837 if (verboseLevel > 3) {
838 if (target.quasi_deutron()) G4cout << " try absorption: ";
839 G4cout << " target " << target.type() << " bullet " << bullet.type()
840 << G4endl;
841 }
842
843 EPCoutput.reset();
844 theEPCollider->collide(&bullet, &target, EPCoutput);
845
846 // If collision failed, exit loop over partners
847 if (EPCoutput.numberOfOutgoingParticles() == 0) break;
848
849 if (verboseLevel > 2) {
850 EPCoutput.printCollisionOutput();
851#ifdef G4CASCADE_CHECK_ECONS
852 balance.collide(&bullet, &target, EPCoutput);
853 balance.okay(); // Do checks, but ignore result
854#endif
855 }
856
857 // Get list of outgoing particles for evaluation
858 std::vector<G4InuclElementaryParticle>& outgoing_particles =
859 EPCoutput.getOutgoingParticles();
860
861 if (!passFermi(outgoing_particles, zone)) continue; // Interaction fails
862
863 // Trailing effect: reject interaction at previously hit nucleon
864 cparticle.propagateAlongThePath(thePartners[i].second);
865 G4ThreeVector new_position = cparticle.getPosition();
866
867 if (!passTrailing(new_position)) continue;
868 collisionPts.push_back(new_position);
869
870 // Sort particles according to beta (fastest first)
871 std::sort(outgoing_particles.begin(), outgoing_particles.end(),
873
874 if (verboseLevel > 2)
875 G4cout << " adding " << outgoing_particles.size()
876 << " output particles" << G4endl;
877
878 // NOTE: Embedded temporary is optimized away (no copying gets done)
879 for (G4int ip = 0; ip < G4int(outgoing_particles.size()); ip++) {
880 outgoing_cparticles.push_back(G4CascadParticle(outgoing_particles[ip],
881 new_position, zone,
882 0.0, 0));
883 }
884
885 no_interaction = false;
886 current_nucl1 = 0;
887 current_nucl2 = 0;
888
889#ifdef G4CASCADE_DEBUG_CHARGE
890 {
891 G4double out_charge = 0.0;
892
893 for (G4int ip = 0; ip < G4int(outgoing_particles.size()); ip++)
894 out_charge += outgoing_particles[ip].getCharge();
895
896 G4cout << " multiplicity " << outgoing_particles.size()
897 << " bul type " << bullet.type()
898 << " targ type " << target.type()
899 << "\n initial charge "
900 << bullet.getCharge() + target.getCharge()
901 << " out charge " << out_charge << G4endl;
902 }
903#endif
904
905 if (verboseLevel > 2)
906 G4cout << " partner type " << target.type() << G4endl;
907
908 if (target.nucleon()) {
909 current_nucl1 = target.type();
910 } else {
911 if (verboseLevel > 2) G4cout << " good absorption " << G4endl;
912 current_nucl1 = (target.type() - 100) / 10;
913 current_nucl2 = target.type() - 100 - 10 * current_nucl1;
914 }
915
916 if (current_nucl1 == 1) {
917 if (verboseLevel > 3) G4cout << " decrement proton count" << G4endl;
918 protonNumberCurrent--;
919 } else {
920 if (verboseLevel > 3) G4cout << " decrement neutron count" << G4endl;
921 neutronNumberCurrent--;
922 }
923
924 if (current_nucl2 == 1) {
925 if (verboseLevel > 3) G4cout << " decrement proton count" << G4endl;
926 protonNumberCurrent--;
927 } else if (current_nucl2 == 2) {
928 if (verboseLevel > 3) G4cout << " decrement neutron count" << G4endl;
929 neutronNumberCurrent--;
930 }
931
932 break;
933 } // loop over partners
934
935 if (no_interaction) { // still no interactions
936 if (verboseLevel > 1) G4cout << " no interaction " << G4endl;
937
938 // For conservation checking (below), get particle before updating
939 static G4InuclElementaryParticle prescatCP; // Avoid memory churn
940 prescatCP = cparticle.getParticle();
941
942 // Last "partner" is just a total-path placeholder
943 cparticle.updatePosition(old_position);
944 cparticle.propagateAlongThePath(thePartners[npart-1].second);
945 cparticle.incrementCurrentPath(thePartners[npart-1].second);
946 boundaryTransition(cparticle);
947 outgoing_cparticles.push_back(cparticle);
948
949 // Check conservation for simple scattering (ignore target nucleus!)
950#ifdef G4CASCADE_CHECK_ECONS
951 if (verboseLevel > 2) {
952 balance.collide(&prescatCP, 0, outgoing_cparticles);
953 balance.okay(); // Report violations, but don't act on them
954 }
955#endif
956 } // if (no_interaction)
957 } // if (npart == 1) [else]
958
959 return;
960}
961
962G4bool G4NucleiModel::passFermi(const std::vector<G4InuclElementaryParticle>& particles,
963 G4int zone) {
964 if (verboseLevel > 1) {
965 G4cout << " >>> G4NucleiModel::passFermi" << G4endl;
966 }
967
968 // Only check Fermi momenta for nucleons
969 for (G4int i = 0; i < G4int(particles.size()); i++) {
970 if (!particles[i].nucleon()) continue;
971
972 G4int type = particles[i].type();
973 G4double mom = particles[i].getMomModule();
974 G4double pfermi = fermi_momenta[type-1][zone];
975
976 if (verboseLevel > 2)
977 G4cout << " type " << type << " p " << mom << " pf " << pfermi << G4endl;
978
979 if (mom < pfermi) {
980 if (verboseLevel > 2) G4cout << " rejected by Fermi" << G4endl;
981 return false;
982 }
983 }
984 return true;
985}
986
987
988// Test here for trailing effect: loop over all previous collision
989// locations and test for d > R_nucleon
990
991G4bool G4NucleiModel::passTrailing(const G4ThreeVector& hit_position) {
992 if (verboseLevel > 1)
993 G4cout << " >>> G4NucleiModel::passTrailing " << hit_position << G4endl;
994
995 G4double dist;
996 for (G4int i = 0; i < G4int(collisionPts.size() ); i++) {
997 dist = (collisionPts[i] - hit_position).mag();
998 if (verboseLevel > 2) G4cout << " dist " << dist << G4endl;
999 if (dist < R_nucleon) {
1000 if (verboseLevel > 2) G4cout << " rejected by Trailing" << G4endl;
1001 return false;
1002 }
1003 }
1004 return true; // New point far enough away to be used
1005}
1006
1007
1008void G4NucleiModel::boundaryTransition(G4CascadParticle& cparticle) {
1009 if (verboseLevel > 1) {
1010 G4cout << " >>> G4NucleiModel::boundaryTransition" << G4endl;
1011 }
1012
1013 G4int zone = cparticle.getCurrentZone();
1014
1015 if (cparticle.movingInsideNuclei() && zone == 0) {
1016 G4cerr << " boundaryTransition-> in zone 0 " << G4endl;
1017 return;
1018 }
1019
1020 G4LorentzVector mom = cparticle.getMomentum();
1021 G4ThreeVector pos = cparticle.getPosition();
1022
1023 G4int type = cparticle.getParticle().type();
1024
1025 G4double r = pos.mag();
1026 G4double pr = pos.dot(mom.vect()) / r;
1027
1028 G4int next_zone = cparticle.movingInsideNuclei() ? zone - 1 : zone + 1;
1029
1030 G4double dv = getPotential(type,zone) - getPotential(type, next_zone);
1031 if (verboseLevel > 3) {
1032 G4cout << "Potentials for type " << type << " = "
1033 << getPotential(type,zone) << " , "
1034 << getPotential(type,next_zone) << G4endl;
1035 }
1036
1037 G4double qv = dv * dv - 2.0 * dv * mom.e() + pr * pr;
1038
1039 G4double p1r;
1040
1041 if (verboseLevel > 3) {
1042 G4cout << " type " << type << " zone " << zone << " next " << next_zone
1043 << " qv " << qv << " dv " << dv << G4endl;
1044 }
1045
1046 if (qv <= 0.0) { // reflection
1047 if (verboseLevel > 3) G4cout << " reflects off boundary" << G4endl;
1048 p1r = -pr;
1049 cparticle.incrementReflectionCounter();
1050 } else { // transition
1051 if (verboseLevel > 3) G4cout << " passes thru boundary" << G4endl;
1052 p1r = std::sqrt(qv);
1053 if(pr < 0.0) p1r = -p1r;
1054 cparticle.updateZone(next_zone);
1055 cparticle.resetReflection();
1056 }
1057
1058 G4double prr = (p1r - pr) / r;
1059
1060 if (verboseLevel > 3) {
1061 G4cout << " prr " << prr << " delta px " << prr*pos.x() << " py "
1062 << prr*pos.y() << " pz " << prr*pos.z() << " mag "
1063 << std::fabs(prr*r) << G4endl;
1064 }
1065
1066 mom.setVect(mom.vect() + pos*prr);
1067 cparticle.updateParticleMomentum(mom);
1068}
1069
1070
1072 if (verboseLevel > 2) {
1073 G4cout << " isProjectile(cparticle): zone "
1074 << cparticle.getCurrentZone() << " of " << number_of_zones
1075 << " path " << cparticle.getCurrentPath()
1076 << " movingInside " << cparticle.movingInsideNuclei()
1077 << " nReflections " << cparticle.getNumberOfReflections()
1078 << G4endl;
1079 }
1080
1081 return (cparticle.getCurrentZone() == number_of_zones-1 && // At surface
1082 cparticle.getCurrentPath() == 1000. && // Input primary
1083 cparticle.getNumberOfReflections() <= 0 && // first pass
1084 (cparticle.movingInsideNuclei()||number_of_zones==1) // inbound
1085 );
1086}
1087
1089 if (verboseLevel > 1) {
1090 G4cout << " >>> G4NucleiModel::worthToPropagate" << G4endl;
1091 }
1092
1093 const G4double ekin_scale = 2.0;
1094
1095 G4bool worth = true;
1096
1097 if (cparticle.reflectedNow()) { // Just reflected -- keep going?
1098 G4int zone = cparticle.getCurrentZone();
1099 G4int ip = cparticle.getParticle().type();
1100
1101 // NOTE: Temporarily backing out use of potential for non-nucleons
1102 G4double ekin_cut = (cparticle.getParticle().nucleon()) ?
1103 getFermiKinetic(ip, zone) : 0.; //*** getPotential(ip, zone);
1104
1105 worth = cparticle.getParticle().getKineticEnergy()/ekin_scale > ekin_cut;
1106
1107 if (verboseLevel > 3) {
1108 G4cout << " type=" << ip
1109 << " ekin=" << cparticle.getParticle().getKineticEnergy()
1110 << " potential=" << ekin_cut
1111 << " : worth? " << worth << G4endl;
1112 }
1113 }
1114
1115 return worth;
1116}
1117
1118
1119G4double G4NucleiModel::getRatio(G4int ip) const {
1120 if (verboseLevel > 3) {
1121 G4cout << " >>> G4NucleiModel::getRatio " << ip << G4endl;
1122 }
1123
1124 switch (ip) {
1125 case proton: return G4double(protonNumberCurrent)/G4double(protonNumber);
1126 case neutron: return G4double(neutronNumberCurrent)/G4double(neutronNumber);
1127 case diproton: return getRatio(proton)*getRatio(proton);
1128 case unboundPN: return getRatio(proton)*getRatio(neutron);
1129 case dineutron: return getRatio(neutron)*getRatio(neutron);
1130 default: return 0.;
1131 }
1132
1133 return 0.;
1134}
1135
1136G4double G4NucleiModel::getCurrentDensity(G4int ip, G4int izone) const {
1137 const G4double pn_spec = 1.0; // Scale factor for pn vs. pp/nn
1138 //const G4double pn_spec = 0.5;
1139
1140 G4double dens = 0.;
1141
1142 if (ip < 100) dens = getDensity(ip,izone);
1143 else { // For dibaryons, remove extra 1/volume term in density product
1144 switch (ip) {
1145 case diproton:
1146 dens = getDensity(proton,izone) * getDensity(proton,izone);
1147 break;
1148 case unboundPN:
1149 dens = getDensity(proton,izone) * getDensity(neutron,izone) * pn_spec;
1150 break;
1151 case dineutron:
1152 dens = getDensity(neutron,izone) * getDensity(neutron,izone);
1153 break;
1154 default: dens = 0.;
1155 }
1156 dens *= getVolume(izone);
1157 }
1158
1159 return getRatio(ip) * dens;
1160}
1161
1162
1165 if (verboseLevel > 1) {
1166 G4cout << " >>> G4NucleiModel::initializeCascad(particle)" << G4endl;
1167 }
1168
1169 const G4double large = 1000.0;
1170
1171 // FIXME: Previous version generated random sin(theta), then used -cos(theta)
1172 // Using generateWithRandomAngles changes result!
1173 // G4ThreeVector pos = generateWithRandomAngles(nuclei_radius).vect();
1174 G4double costh = std::sqrt(1.0 - inuclRndm());
1175 G4ThreeVector pos = generateWithFixedTheta(-costh, nuclei_radius);
1176
1177 G4CascadParticle cpart(*particle, pos, number_of_zones, large, 0);
1178
1179 if (verboseLevel > 2) G4cout << cpart << G4endl;
1180
1181 return cpart;
1182}
1183
1185 G4InuclNuclei* target,
1186 modelLists& output) {
1187 if (verboseLevel) {
1188 G4cout << " >>> G4NucleiModel::initializeCascad(bullet,target,output)"
1189 << G4endl;
1190 }
1191
1192 const G4double large = 1000.0;
1193 const G4double max_a_for_cascad = 5.0;
1194 const G4double ekin_cut = 2.0;
1195 const G4double small_ekin = 1.0e-6;
1196 const G4double r_large2for3 = 62.0;
1197 const G4double r0forAeq3 = 3.92;
1198 const G4double s3max = 6.5;
1199 const G4double r_large2for4 = 69.14;
1200 const G4double r0forAeq4 = 4.16;
1201 const G4double s4max = 7.0;
1202 const G4int itry_max = 100;
1203
1204 // Convenient references to input buffer contents
1205 std::vector<G4CascadParticle>& casparticles = output.first;
1206 std::vector<G4InuclElementaryParticle>& particles = output.second;
1207
1208 casparticles.clear(); // Reset input buffer to fill this event
1209 particles.clear();
1210
1211 // first decide whether it will be cascad or compound final nuclei
1212 G4int ab = bullet->getA();
1213 G4int zb = bullet->getZ();
1214 G4int at = target->getA();
1215 G4int zt = target->getZ();
1216
1217 G4double massb = bullet->getMass(); // For creating LorentzVectors below
1218
1219 if (ab < max_a_for_cascad) {
1220
1221 G4double benb = bindingEnergy(ab,zb)/GeV / G4double(ab);
1222 G4double bent = bindingEnergy(at,zt)/GeV / G4double(at);
1223 G4double ben = benb < bent ? bent : benb;
1224
1225 if (bullet->getKineticEnergy()/ab > ekin_cut*ben) {
1226 G4int itryg = 0;
1227
1228 while (casparticles.size() == 0 && itryg < itry_max) {
1229 itryg++;
1230 particles.clear();
1231
1232 // nucleons coordinates and momenta in nuclei rest frame
1233 coordinates.clear();
1234 momentums.clear();
1235
1236 if (ab < 3) { // deuteron, simplest case
1237 G4double r = 2.214 - 3.4208 * std::log(1.0 - 0.981 * inuclRndm());
1239 coordinates.push_back(coord1);
1240 coordinates.push_back(-coord1);
1241
1242 G4double p = 0.0;
1243 G4bool bad = true;
1244 G4int itry = 0;
1245
1246 while (bad && itry < itry_max) {
1247 itry++;
1248 p = 456.0 * inuclRndm();
1249
1250 if (p * p / (p * p + 2079.36) / (p * p + 2079.36) > 1.2023e-4 * inuclRndm() &&
1251 p * r > 312.0) bad = false;
1252 }
1253
1254 if (itry == itry_max)
1255 if (verboseLevel > 2){
1256 G4cout << " deutron bullet generation-> itry = " << itry_max << G4endl;
1257 }
1258
1259 p = 0.0005 * p;
1260
1261 if (verboseLevel > 2){
1262 G4cout << " p nuc " << p << G4endl;
1263 }
1264
1266
1267 momentums.push_back(mom);
1268 mom.setVect(-mom.vect());
1269 momentums.push_back(-mom);
1270 } else {
1271 G4ThreeVector coord1;
1272
1273 G4bool badco = true;
1274
1275 G4int itry = 0;
1276
1277 if (ab == 3) {
1278 while (badco && itry < itry_max) {
1279 if (itry > 0) coordinates.clear();
1280 itry++;
1281 G4int i(0);
1282
1283 for (i = 0; i < 2; i++) {
1284 G4int itry1 = 0;
1285 G4double ss, u, rho;
1286 G4double fmax = std::exp(-0.5) / std::sqrt(0.5);
1287
1288 while (itry1 < itry_max) {
1289 itry1++;
1290 ss = -std::log(inuclRndm());
1291 u = fmax * inuclRndm();
1292 rho = std::sqrt(ss) * std::exp(-ss);
1293
1294 if (rho > u && ss < s3max) {
1295 ss = r0forAeq3 * std::sqrt(ss);
1296 coord1 = generateWithRandomAngles(ss).vect();
1297 coordinates.push_back(coord1);
1298
1299 if (verboseLevel > 2){
1300 G4cout << " i " << i << " r " << coord1.mag() << G4endl;
1301 }
1302 break;
1303 }
1304 }
1305
1306 if (itry1 == itry_max) { // bad case
1307 coord1.set(10000.,10000.,10000.);
1308 coordinates.push_back(coord1);
1309 break;
1310 }
1311 }
1312
1313 coord1 = -coordinates[0] - coordinates[1];
1314 if (verboseLevel > 2) {
1315 G4cout << " 3 r " << coord1.mag() << G4endl;
1316 }
1317
1318 coordinates.push_back(coord1);
1319
1320 G4bool large_dist = false;
1321
1322 for (i = 0; i < 2; i++) {
1323 for (G4int j = i+1; j < 3; j++) {
1324 G4double r2 = (coordinates[i]-coordinates[j]).mag2();
1325
1326 if (verboseLevel > 2) {
1327 G4cout << " i " << i << " j " << j << " r2 " << r2 << G4endl;
1328 }
1329
1330 if (r2 > r_large2for3) {
1331 large_dist = true;
1332
1333 break;
1334 }
1335 }
1336
1337 if (large_dist) break;
1338 }
1339
1340 if(!large_dist) badco = false;
1341
1342 }
1343
1344 } else { // a >= 4
1345 G4double b = 3./(ab - 2.0);
1346 G4double b1 = 1.0 - b / 2.0;
1347 G4double u = b1 + std::sqrt(b1 * b1 + b);
1348 G4double fmax = (1.0 + u/b) * u * std::exp(-u);
1349
1350 while (badco && itry < itry_max) {
1351
1352 if (itry > 0) coordinates.clear();
1353 itry++;
1354 G4int i(0);
1355
1356 for (i = 0; i < ab-1; i++) {
1357 G4int itry1 = 0;
1358 G4double ss;
1359
1360 while (itry1 < itry_max) {
1361 itry1++;
1362 ss = -std::log(inuclRndm());
1363 u = fmax * inuclRndm();
1364
1365 if (std::sqrt(ss) * std::exp(-ss) * (1.0 + ss/b) > u
1366 && ss < s4max) {
1367 ss = r0forAeq4 * std::sqrt(ss);
1368 coord1 = generateWithRandomAngles(ss).vect();
1369 coordinates.push_back(coord1);
1370
1371 if (verboseLevel > 2) {
1372 G4cout << " i " << i << " r " << coord1.mag() << G4endl;
1373 }
1374
1375 break;
1376 }
1377 }
1378
1379 if (itry1 == itry_max) { // bad case
1380 coord1.set(10000.,10000.,10000.);
1381 coordinates.push_back(coord1);
1382 break;
1383 }
1384 }
1385
1386 coord1 *= 0.0; // Cheap way to reset
1387 for(G4int j = 0; j < ab -1; j++) coord1 -= coordinates[j];
1388
1389 coordinates.push_back(coord1);
1390
1391 if (verboseLevel > 2){
1392 G4cout << " last r " << coord1.mag() << G4endl;
1393 }
1394
1395 G4bool large_dist = false;
1396
1397 for (i = 0; i < ab-1; i++) {
1398 for (G4int j = i+1; j < ab; j++) {
1399
1400 G4double r2 = (coordinates[i]-coordinates[j]).mag2();
1401
1402 if (verboseLevel > 2){
1403 G4cout << " i " << i << " j " << j << " r2 " << r2 << G4endl;
1404 }
1405
1406 if (r2 > r_large2for4) {
1407 large_dist = true;
1408
1409 break;
1410 }
1411 }
1412
1413 if (large_dist) break;
1414 }
1415
1416 if (!large_dist) badco = false;
1417 }
1418 }
1419
1420 if(badco) {
1421 G4cout << " can not generate the nucleons coordinates for a "
1422 << ab << G4endl;
1423
1424 casparticles.clear(); // Return empty buffer on error
1425 particles.clear();
1426 return;
1427
1428 } else { // momentums
1429 G4double p, u, x;
1430 G4LorentzVector mom;
1431 //G4bool badp = True;
1432
1433 for (G4int i = 0; i < ab - 1; i++) {
1434 G4int itry2 = 0;
1435
1436 while(itry2 < itry_max) {
1437 itry2++;
1438 u = -std::log(0.879853 - 0.8798502 * inuclRndm());
1439 x = u * std::exp(-u);
1440
1441 if(x > inuclRndm()) {
1442 p = std::sqrt(0.01953 * u);
1443 mom = generateWithRandomAngles(p, massb);
1444 momentums.push_back(mom);
1445
1446 break;
1447 }
1448 } // while (itry2 < itry_max)
1449
1450 if(itry2 == itry_max) {
1451 G4cout << " can not generate proper momentum for a "
1452 << ab << G4endl;
1453
1454 casparticles.clear(); // Return empty buffer on error
1455 particles.clear();
1456 return;
1457 }
1458 } // for (i=0 ...
1459 // last momentum
1460
1461 mom *= 0.; // Cheap way to reset
1462 mom.setE(bullet->getEnergy()+target->getEnergy());
1463
1464 for(G4int j=0; j< ab-1; j++) mom -= momentums[j];
1465
1466 momentums.push_back(mom);
1467 }
1468 }
1469
1470 // Coordinates and momenta at rest are generated, now back to the lab
1471 G4double rb = 0.0;
1472 G4int i(0);
1473
1474 for(i = 0; i < G4int(coordinates.size()); i++) {
1475 G4double rp = coordinates[i].mag2();
1476
1477 if(rp > rb) rb = rp;
1478 }
1479
1480 // nuclei i.p. as a whole
1481 G4double s1 = std::sqrt(inuclRndm());
1482 G4double phi = randomPHI();
1483 G4double rz = (nuclei_radius + rb) * s1;
1484 G4ThreeVector global_pos(rz*std::cos(phi), rz*std::sin(phi),
1485 -(nuclei_radius+rb)*std::sqrt(1.0-s1*s1));
1486
1487 for (i = 0; i < G4int(coordinates.size()); i++) {
1488 coordinates[i] += global_pos;
1489 }
1490
1491 // all nucleons at rest
1492 raw_particles.clear();
1493
1494 for (G4int ipa = 0; ipa < ab; ipa++) {
1495 G4int knd = ipa < zb ? 1 : 2;
1496 raw_particles.push_back(G4InuclElementaryParticle(momentums[ipa], knd));
1497 }
1498
1499 G4InuclElementaryParticle dummy(small_ekin, 1);
1500 G4LorentzConvertor toTheBulletRestFrame(&dummy, bullet);
1501 toTheBulletRestFrame.toTheTargetRestFrame();
1502
1503 particleIterator ipart;
1504
1505 for (ipart = raw_particles.begin(); ipart != raw_particles.end(); ipart++) {
1506 ipart->setMomentum(toTheBulletRestFrame.backToTheLab(ipart->getMomentum()));
1507 }
1508
1509 // fill cascad particles and outgoing particles
1510
1511 for(G4int ip = 0; ip < G4int(raw_particles.size()); ip++) {
1512 G4LorentzVector mom = raw_particles[ip].getMomentum();
1513 G4double pmod = mom.rho();
1514 G4double t0 = -mom.vect().dot(coordinates[ip]) / pmod;
1515 G4double det = t0 * t0 + nuclei_radius * nuclei_radius
1516 - coordinates[ip].mag2();
1517 G4double tr = -1.0;
1518
1519 if(det > 0.0) {
1520 G4double t1 = t0 + std::sqrt(det);
1521 G4double t2 = t0 - std::sqrt(det);
1522
1523 if(std::fabs(t1) <= std::fabs(t2)) {
1524 if(t1 > 0.0) {
1525 if(coordinates[ip].z() + mom.z() * t1 / pmod <= 0.0) tr = t1;
1526 }
1527
1528 if(tr < 0.0 && t2 > 0.0) {
1529
1530 if(coordinates[ip].z() + mom.z() * t2 / pmod <= 0.0) tr = t2;
1531 }
1532
1533 } else {
1534 if(t2 > 0.0) {
1535
1536 if(coordinates[ip].z() + mom.z() * t2 / pmod <= 0.0) tr = t2;
1537 }
1538
1539 if(tr < 0.0 && t1 > 0.0) {
1540 if(coordinates[ip].z() + mom.z() * t1 / pmod <= 0.0) tr = t1;
1541 }
1542 }
1543
1544 }
1545
1546 if(tr >= 0.0) { // cascad particle
1547 coordinates[ip] += mom.vect()*tr / pmod;
1548 casparticles.push_back(G4CascadParticle(raw_particles[ip],
1549 coordinates[ip],
1550 number_of_zones, large, 0));
1551
1552 } else {
1553 particles.push_back(raw_particles[ip]);
1554 }
1555 }
1556 }
1557
1558 if(casparticles.size() == 0) {
1559 particles.clear();
1560
1561 G4cout << " can not generate proper distribution for " << itry_max
1562 << " steps " << G4endl;
1563 }
1564 }
1565 }
1566
1567 if(verboseLevel > 2){
1568 G4int ip(0);
1569
1570 G4cout << " cascad particles: " << casparticles.size() << G4endl;
1571 for(ip = 0; ip < G4int(casparticles.size()); ip++)
1572 G4cout << casparticles[ip] << G4endl;
1573
1574 G4cout << " outgoing particles: " << particles.size() << G4endl;
1575 for(ip = 0; ip < G4int(particles.size()); ip++)
1576 G4cout << particles[ip] << G4endl;
1577 }
1578
1579 return; // Buffer has been filled
1580}
1581
1582
1583// Compute mean free path for inclusive interaction of projectile and target
1584G4double
1585G4NucleiModel::inverseMeanFreePath(const G4CascadParticle& cparticle,
1586 const G4InuclElementaryParticle& target) {
1587 G4int ptype = cparticle.getParticle().type();
1588 G4int ip = target.type();
1589 G4int zone = cparticle.getCurrentZone();
1590
1591 // Configure frame transformation to get kinetic energy for lookups
1592 dummy_convertor.setBullet(cparticle.getParticle());
1593 dummy_convertor.setTarget(&target);
1594 G4double ekin = dummy_convertor.getKinEnergyInTheTRS();
1595
1596 // Get cross section for interacting with target (dibaryons are absorptive)
1597 G4double csec = (ip < 100) ? totalCrossSection(ekin, ptype*ip)
1598 : absorptionCrossSection(ekin, ptype);
1599
1600 if(verboseLevel > 2) {
1601 G4cout << " ip " << ip << " ekin " << ekin << " csec " << csec << G4endl;
1602 }
1603
1604 if (csec <= 0.) return 0.; // No interaction, avoid unnecessary work
1605
1606 // Get nuclear density and compute mean free path
1607 return csec * getCurrentDensity(ip, zone);
1608}
1609
1610// Throw random distance for interaction of particle using path and MFP
1611
1612G4double
1613G4NucleiModel::generateInteractionLength(const G4CascadParticle& cparticle,
1614 G4double path, G4double invmfp) const {
1615 // Delay interactions of newly formed secondaries (minimum int. length)
1616 const G4double young_cut = std::sqrt(10.0) * 0.25;
1617 const G4double huge_num = 50.0; // Argument to exponential
1618 const G4double small = 1.0e-10;
1619
1620 if (invmfp < small) return path; // No interaction, avoid unnecessary work
1621
1622 G4double pw = -path * invmfp;
1623 if (pw < -huge_num) pw = -huge_num;
1624 pw = 1.0 - std::exp(pw);
1625
1626 if (verboseLevel > 2)
1627 G4cout << " mfp " << 1./invmfp << " pw " << pw << G4endl;
1628
1629 G4double spath = path; // Buffer for return value
1630
1631 // Primary particle(s) should always interact at least once
1632 if (isProjectile(cparticle) || (inuclRndm() < pw)) {
1633 spath = -std::log(1.0 - pw * inuclRndm()) / invmfp;
1634 if (cparticle.young(young_cut, spath)) spath = path;
1635
1636 if (verboseLevel > 2)
1637 G4cout << " spath " << spath << " path " << path << G4endl;
1638 }
1639
1640 return spath;
1641}
1642
1643
1644// Parametrized cross section for pion and photon absorption
1645
1647 if (type != pionPlus && type != pionMinus && type != pionZero &&
1648 type != photon) {
1649 G4cerr << " absorptionCrossSection only valid for incident pions" << G4endl;
1650 return 0.;
1651 }
1652
1653 G4double csec = 0.;
1654
1655 // Pion absorption is parametrized for low vs. medium energy
1656 if (type == pionPlus || type == pionMinus || type == pionZero) {
1657 if (ke < 0.3) csec = (0.1106 / std::sqrt(ke) - 0.8
1658 + 0.08 / ((ke-0.123)*(ke-0.123) + 0.0056) );
1659 else if (ke < 1.0) csec = 3.6735 * (1.0-ke)*(1.0-ke);
1660 }
1661
1662 // Photon cross-section is binned from parametrization by Mi. Kossov
1663 if (type == photon) {
1664 static const G4double gammaQDscale = G4CascadeParameters::gammaQDScale();
1665 static const G4double kebins[] =
1666 { 0.0, 0.01, 0.013, 0.018, 0.024, 0.032, 0.042, 0.056, 0.075, 0.1,
1667 0.13, 0.18, 0.24, 0.32, 0.42, 0.56, 0.75, 1.0, 1.3, 1.8,
1668 2.4, 3.2, 4.2, 5.6, 7.5, 10.0, 13.0, 18.0, 24.0, 32.0 };
1669 static const G4double gammaD[] =
1670 { 2.8, 1.3, 0.89, 0.56, 0.38, 0.27, 0.19, 0.14, 0.098,
1671 0.071, 0.054, 0.0003, 0.0007, 0.0027, 0.0014, 0.001, 0.0012, 0.0005,
1672 0.0003, 0.0002,0.0002, 0.0002, 0.0002, 0.0002, 0.0001, 0.0001, 0.0001,
1673 0.0001, 0.0001, 0.0001 };
1674 static G4CascadeInterpolator<30> interp(kebins);
1675 csec = interp.interpolate(ke, gammaD) * gammaQDscale;
1676 }
1677
1678 if (csec < 0.0) csec = 0.0;
1679
1680 if (verboseLevel > 2) {
1681 G4cout << " ekin " << ke << " abs. csec " << csec << " mb" << G4endl;
1682 }
1683
1684 return crossSectionUnits * csec;
1685}
1686
1687
1689{
1690 // All scattering cross-sections are available from tables
1691 const G4CascadeChannel* xsecTable = G4CascadeChannelTables::GetTable(rtype);
1692 if (!xsecTable) {
1693 G4cerr << " unknown collison type = " << rtype << G4endl;
1694 return 0.;
1695 }
1696
1697 return (crossSectionUnits * xsecTable->getCrossSection(ke));
1698}
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:60
@ photon
@ neutron
std::vector< G4InuclElementaryParticle >::iterator particleIterator
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
double z() const
double x() const
double y() const
double dot(const Hep3Vector &) const
double mag() const
void set(double x, double y, double z)
Hep3Vector vect() const
void setVect(const Hep3Vector &)
G4bool reflectedNow() const
G4double getCurrentPath() const
void updateZone(G4int izone)
void incrementCurrentPath(G4double npath)
const G4InuclElementaryParticle & getParticle() const
void updateParticleMomentum(const G4LorentzVector &mom)
void updatePosition(const G4ThreeVector &pos)
G4bool movingInsideNuclei() const
void propagateAlongThePath(G4double path)
G4int getNumberOfReflections() const
G4LorentzVector getMomentum() const
G4bool young(G4double young_path_cut, G4double cpath) const
void incrementReflectionCounter()
G4int getCurrentZone() const
G4double getPathToTheNextZone(G4double rz_in, G4double rz_out)
const G4ThreeVector & getPosition() const
static const G4CascadeChannel * GetTable(G4int initialState)
virtual G4double getCrossSection(double ke) const =0
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &output)
G4double interpolate(const G4double x, const G4double(&yb)[nBins]) const
static G4double xsecScale()
static G4double radiusScale()
static G4bool useTwoParam()
static G4double gammaQDScale()
static G4double radiusAlpha()
static G4double radiusSmall()
static G4double fermiScale()
static G4double radiusTrailing()
G4int numberOfOutgoingParticles() const
void printCollisionOutput(std::ostream &os=G4cout) const
const std::vector< G4InuclElementaryParticle > & getOutgoingParticles() const
void collide(G4InuclParticle *bullet, G4InuclParticle *target, G4CollisionOutput &output)
static G4double getParticleMass(G4int type)
G4int getZ() const
G4int getA() const
G4double getKineticEnergy() const
G4double getMass() const
G4double getEnergy() const
G4double getCharge() const
void setBullet(const G4InuclParticle *bullet)
G4LorentzVector backToTheLab(const G4LorentzVector &mom) const
G4double getKinEnergyInTheTRS() const
void setTarget(const G4InuclParticle *target)
G4double getDensity(G4int ip, G4int izone) const
G4CascadParticle initializeCascad(G4InuclElementaryParticle *particle)
G4bool worthToPropagate(const G4CascadParticle &cparticle) const
void generateParticleFate(G4CascadParticle &cparticle, G4ElementaryParticleCollider *theEPCollider, std::vector< G4CascadParticle > &cascade)
void printModel() const
G4InuclElementaryParticle generateNucleon(G4int type, G4int zone) const
G4double getVolume(G4int izone) const
G4double totalCrossSection(G4double ke, G4int rtype) const
G4double getPotential(G4int ip, G4int izone) const
G4double getFermiMomentum(G4int ip, G4int izone) const
void reset(G4int nHitNeutrons=0, G4int nHitProtons=0, const std::vector< G4ThreeVector > *hitPoints=0)
G4double absorptionCrossSection(G4double e, G4int type) const
void generateModel(G4InuclNuclei *nuclei)
std::pair< std::vector< G4CascadParticle >, std::vector< G4InuclElementaryParticle > > modelLists
G4double getFermiKinetic(G4int ip, G4int izone) const
G4LorentzVector generateNucleonMomentum(G4int type, G4int zone) const
G4bool isProjectile(const G4CascadParticle &cparticle) const
virtual void setVerboseLevel(G4int verbose=0)
const G4double pi
G4LorentzVector generateWithFixedTheta(G4double ct, G4double p, G4double mass=0.)
G4double bindingEnergy(G4int A, G4int Z)
G4LorentzVector generateWithRandomAngles(G4double p, G4double mass=0.)