Geant4-11
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//
27// For the best approximation to a physical-units model, set the following:
28// setenv G4NUCMODEL_XSEC_SCALE 0.1
29// setenv G4NUCMODEL_RAD_SCALE 1.0
30// setenv G4NUCMODEL_RAD_2PAR 1
31// setenv G4NUCMODEL_RAD_SMALL 1.992
32// setenv G4NUCMODEL_RAD_ALPHA 0.84
33// setenv G4NUCMODEL_FERMI_SCALE 0.685
34// setenv G4NUCMODEL_RAD_TRAILING 0.70
35//
36// 20100112 M. Kelsey -- Remove G4CascadeMomentum, use G4LorentzVector directly
37// 20100114 M. Kelsey -- Use G4ThreeVector for position
38// 20100309 M. Kelsey -- Use new generateWithRandomAngles for theta,phi stuff;
39// eliminate some unnecessary std::pow(), move ctor here
40// 20100407 M. Kelsey -- Replace std::vector<>::resize(0) with ::clear().
41// Move output vectors from generateParticleFate() and
42// ::generateInteractionPartners() to be data members; return via
43// const-ref instead of by value.
44// 20100413 M. Kelsey -- Pass G4CollisionOutput by ref to ::collide()
45// 20100418 M. Kelsey -- Reference output particle lists via const-ref
46// 20100421 M. Kelsey -- Replace hardwired p/n masses with G4PartDef's
47// 20100517 M. Kelsey -- Use G4CascadeINterpolator for cross-section
48// calculations. use G4Cascade*Channel for total xsec in pi-N
49// N-N channels. Move absorptionCrossSection() from SpecialFunc.
50// 20100610 M. Kelsey -- Replace another random-angle code block; add some
51// diagnostic output for partner-list production.
52// 20100617 M. Kelsey -- Replace preprocessor flag CHC_CHECK with
53// G4CASCADE_DEBUG_CHARGE
54// 20100620 M. Kelsey -- Improve error message on empty partners list, add
55// four-momentum checking after EPCollider
56// 20100621 M. Kelsey -- In boundaryTransition() account for momentum transfer
57// to secondary by boosting into recoil nucleus "rest" frame.
58// Don't need temporaries to create dummy "partners" for list.
59// 20100622 M. Kelsey -- Restore use of "bindingEnergy()" function name, which
60// is now a wrapper for G4NucleiProperties::GetBindingEnergy().
61// 20100623 M. Kelsey -- Eliminate some temporaries terminating partner-list,
62// discard recoil boost for now. Initialize all data
63// members in ctors. Allow generateModel() to be called
64// mutliple times, by clearing vectors each time through;
65// avoid extra work by returning if A and Z are same as
66// before.
67// 20100628 M. Kelsey -- Two momentum-recoil bugs; don't subtract energies!
68// 20100715 M. Kelsey -- Make G4InuclNuclei in generateModel(), use for
69// balance checking (only verbose>2) in generateParticleFate().
70// 20100721 M. Kelsey -- Use new G4CASCADE_CHECK_ECONS for conservation checks
71// 20100723 M. Kelsey -- Move G4CollisionOutput buffer to .hh for reuse
72// 20100726 M. Kelsey -- Preallocate arrays with number_of_zones dimension.
73// 20100902 M. Kelsey -- Remove resize(3) directives from qdeutron/acsecs
74// 20100907 M. Kelsey -- Limit interaction targets based on current nucleon
75// counts (e.g., no pp if protonNumberCurrent < 2!)
76// 20100914 M. Kelsey -- Migrate to integer A and Z
77// 20100928 M. Kelsey -- worthToPropagate() use nuclear potential for hadrons.
78// 20101005 M. Kelsey -- Annotate hardwired numbers for upcoming rationalizing,
79// move hardwired numbers out to static data members, factor out
80// all the model construction pieces to separate functions, fix
81// wrong value for 4/3 pi.
82// 20101019 M. Kelsey -- CoVerity reports: unitialized constructor, dtor leak
83// 20101020 M. Kelsey -- Bug fixes to refactoring changes (5 Oct). Back out
84// worthToPropagate() changes for better regression testing.
85// 20101020 M. Kelsey -- Re-activate worthToPropagate() changes.
86// 20101119 M. Kelsey -- Hide "negative path" and "no partners" messages in
87// verbosity.
88// 20110218 M. Kelsey -- Add crossSectionUnits and radiusUnits scale factors,
89// use "theoretical" numbers for radii etc., multipled by scale
90// factor; set scale factors using environment variables
91// 20110303 M. Kelsey -- Add comments why using fabs() with B.E. differences?
92// 20110321 M. Kelsey -- Replace strtof() with strtod() for envvar conversion
93// 20110321 M. Kelsey -- Use fm and fm^2 as default units, Per D. Wright
94// (NOTE: Restored from original 20110318 commit)
95// 20110324 D. Wright -- Implement trailing effect
96// 20110324 M. Kelsey -- Move ::reset() here, as it has more code.
97// 20110519 M. Kelsey -- Used "rho" after assignment, instead of recomputing
98// 20110525 M. Kelsey -- Revert scale factor changes (undo 20110321 changes)
99// 20110617 M. Kelsey -- Apply scale factor to trailing-effect radius, make
100// latter runtime adjustable (G4NUCMODEL_RAD_TRAILING)
101// 20110720 M. Kelsey -- Follow interface change for cross-section tables,
102// eliminating switch blocks.
103// 20110806 M. Kelsey -- Reduce memory churn by pre-allocating buffers
104// 20110823 M. Kelsey -- Remove local cross-section tables entirely
105// 20110825 M. Kelsey -- Add comments regarding Fermi momentum scale, set of
106// "best guess" parameter values
107// 20110831 M. Kelsey -- Make "best guess" parameters the defaults
108// 20110922 M. Kelsey -- Follow migrations G4InuclParticle::print(ostream&)
109// and G4CascadParticle::print(ostream&)
110// 20111018 M. Kelsey -- Correct kaon potential to be positive, not negative
111// 20111107 M. Kelsey -- *** REVERT TO OLD NON-PHYSICAL PARAMETERS FOR 9.5 ***
112// 20120306 M. Kelsey -- For incident projectile (CParticle incoming to
113// nucleus from outside) don't use pw cut; force interaction.
114// 20120307 M. Kelsey -- Compute zone volumes during setup, not during each
115// interaction. Include photons in absorption by dibaryons.
116// Discard unused code in totalCrossSection(). Encapsulate
117// interaction path calculations in generateInteractionLength().
118// 20120308 M. Kelsey -- Add binned photon absorption cross-section.
119// 20120320 M. Kelsey -- Bug fix: fill zone_volumes for single-zone case
120// 20120321 M. Kelsey -- Add check in isProjectile() for single-zone case.
121// 20120608 M. Kelsey -- Fix variable-name "shadowing" compiler warnings.
122// 20120822 M. Kelsey -- Move envvars to G4CascadeParameters.
123// 20121205 M. Kelsey -- Bug fix in generateParticleFate(): daughters should
124// have generation count incremented from parent. This allows
125// isProjectile() to simply check generation number == 0.
126// 20130221 M. Kelsey -- For incident photons, move to random point along
127// trajectory through nucleus for first (forced) interaction.
128// 20130223 M. Kelsey -- Fix bugs in weighting and selection of traj points
129// 20130302 M. Kelsey -- Replace "isProjectile()" with "forceFirst()" in
130// path-length calculation.
131// 20130314 M. Kelsey -- Restore null initializer and if-block for _TLS_.
132// 20130508 D. Wright -- Support muon-nuclear interactions.
133// 20130510 M. Kelsey -- Check for zero-momentum in choosePointAlongTraj().
134// 20130511 M. Kelsey -- Move choosePointAlongTraj() to initializeCascad(),
135// instead of after generateIP(). Change "spath<path" to <= to
136// handle at-rest particles. Skip mu-/neutron interactions.
137// Hide "zone 0 transition" error message behind verbosity.
138// 20130523 M. Kelsey -- Bug fix: Initialize dummy_convertor in inverseMFP.
139// For capture-at-rest, replace exterior nucleus with outer zone.
140// 20130524 M. Kelsey -- Move "large" and "small" cutoffs to static const.
141// Check zone argument to inverseMFP() to ensure within range.
142// 20130611 M. Kelsey -- Undo "spath<=path" change (20130511), replace with
143// explicit check on spath==0.
144// 20130619 A. Ribon -- Fixed reproducibility problem in the method
145// generateParticleFate
146// 20130627 M. Kelsey -- Check "path==0.", rather than spath.
147// 20130628 M. Kelsey -- Print deuteron type code, rather than array index,
148// Extend useQuasiDeuteron() to check good absorption
149// 20130701 M. Kelsey -- Don't average 1/MFP for total interaction; just sum;
150// inverseMFP() returns "large" value instead of input path.
151// 20130806 M. Kelsey -- Move static G4InuclEP's to file scope to avoid
152// thread collisions
153// 20130924 M. Kelsey -- Use G4Log, G4Exp, G4Pow for CPU speedup
154// 20130925 M. Kelsey -- Eliminate unnecessary use of pow() in integrals
155// 20131001 M. Kelsey -- Move QDinterp object to data member, thread local
156// 20140116 M. Kelsey -- Move all envvar parameters to const data members
157// 20141001 M. Kelsey -- Change sign of "dv" in boundaryTransition
158// 20150608 M. Kelsey -- Label all while loops as terminating.
159// 20150622 M. Kelsey -- Use G4AutoDelete for _TLS_ buffers.
160// 20180227 A. Ribon -- Replaced obsolete std::bind2nd with std::bind
161
162#include "G4NucleiModel.hh"
163#include "G4AutoDelete.hh"
164#include "G4CascadeChannel.hh"
168#include "G4CascadeParameters.hh"
169#include "G4CollisionOutput.hh"
171#include "G4Exp.hh"
172#include "Randomize.hh"
173#include "G4HadTmpUtil.hh"
174#include "G4InuclNuclei.hh"
177#include "G4Log.hh"
178#include "G4LorentzConvertor.hh"
179#include "G4Neutron.hh"
182#include "G4PhysicalConstants.hh"
183#include "G4Proton.hh"
184#include "G4SystemOfUnits.hh"
185#include <vector>
186#include <algorithm>
187#include <functional>
188#include <numeric>
189
190using namespace G4InuclParticleNames;
191using namespace G4InuclSpecialFunctions;
192
193typedef std::vector<G4InuclElementaryParticle>::iterator particleIterator;
194
195// Cut-off limits for extreme values in calculations
196const G4double G4NucleiModel::small = 1.0e-9;
197const G4double G4NucleiModel::large = 1000.;
198
199// For convenience in computing densities
201
202// Zone boundaries as fraction of nuclear radius (from outside in)
203const G4double G4NucleiModel::alfa3[3] = { 0.7, 0.3, 0.01 };
204const G4double G4NucleiModel::alfa6[6] = { 0.9, 0.6, 0.4, 0.2, 0.1, 0.05 };
205
206// Flat nuclear potentials for mesons and hyperons (GeV)
207const G4double G4NucleiModel::pion_vp = 0.007;
209const G4double G4NucleiModel::kaon_vp = 0.015;
211
212namespace {
213 // Quasideuteron absorption cross section is taken to be the
214 // deuteron photo-disintegration cross section (gam + D -> p + n)
215 // Values taken from smooth curve drawn through data from 2.4 - 500 MeV,
216 // plus angle integration of JLab data from 0.5 - 3 GeV
217 // M. Mirazita et al., Phys. Rev. C 70, 014005 (2004)
218 // Points above 3 GeV are extrapolated.
219 const G4double kebins[] = {
220 0.0, 0.0024, 0.0032, 0.0042, 0.0056, 0.0075, 0.01, 0.024, 0.075, 0.1,
221 0.13, 0.18, 0.24, 0.32, 0.42, 0.56, 0.75, 1.0, 1.3, 1.8,
222 2.4, 3.2, 4.2, 5.6, 7.5, 10.0, 13.0, 18.0, 24.0, 32.0 };
223
224 const G4double gammaQDxsec[30] = {
225 0.0, 0.7, 2.0, 2.2, 2.1, 1.8, 1.3, 0.4, 0.098, 0.071,
226 0.055, 0.055, 0.065, 0.045, 0.017, 0.007, 2.37e-3, 6.14e-4, 1.72e-4, 4.2e-5,
227 1.05e-5, 3.0e-6, 7.0e-7, 1.3e-7, 2.3e-8, 3.2e-9, 4.9e-10, 0.0, 0.0, 0.0 };
228}
229
230
231// Constructors and destructor
233 : verboseLevel(0), nuclei_radius(0.), nuclei_volume(0.), number_of_zones(0),
234 A(0), Z(0), theNucleus(0), neutronNumber(0), protonNumber(0),
235 neutronNumberCurrent(0), protonNumberCurrent(0), current_nucl1(0),
236 current_nucl2(0), gammaQDinterp(kebins),
237 crossSectionUnits(G4CascadeParameters::xsecScale()),
238 radiusUnits(G4CascadeParameters::radiusScale()),
239 skinDepth(0.611207*radiusUnits),
240 radiusScale((G4CascadeParameters::useTwoParam()?1.16:1.2)*radiusUnits),
241 radiusScale2((G4CascadeParameters::useTwoParam()?-1.3456:0.)*radiusUnits),
242 radiusForSmall(G4CascadeParameters::radiusSmall()),
243 radScaleAlpha(G4CascadeParameters::radiusAlpha()),
244 fermiMomentum(G4CascadeParameters::fermiScale()),
245 R_nucleon(G4CascadeParameters::radiusTrailing()),
246 gammaQDscale(G4CascadeParameters::gammaQDScale()),
247 potentialThickness(1.0),
248 neutronEP(neutron), protonEP(proton) {}
249
251 : verboseLevel(0), nuclei_radius(0.), nuclei_volume(0.), number_of_zones(0),
252 A(0), Z(0), theNucleus(0), neutronNumber(0), protonNumber(0),
253 neutronNumberCurrent(0), protonNumberCurrent(0), current_nucl1(0),
254 current_nucl2(0), gammaQDinterp(kebins),
255 crossSectionUnits(G4CascadeParameters::xsecScale()),
256 radiusUnits(G4CascadeParameters::radiusScale()),
257 skinDepth(0.611207*radiusUnits),
258 radiusScale((G4CascadeParameters::useTwoParam()?1.16:1.2)*radiusUnits),
259 radiusScale2((G4CascadeParameters::useTwoParam()?-1.3456:0.)*radiusUnits),
260 radiusForSmall(G4CascadeParameters::radiusSmall()),
261 radScaleAlpha(G4CascadeParameters::radiusAlpha()),
262 fermiMomentum(G4CascadeParameters::fermiScale()),
263 R_nucleon(G4CascadeParameters::radiusTrailing()),
264 gammaQDscale(G4CascadeParameters::gammaQDScale()),
265 potentialThickness(1.0),
266 neutronEP(neutron), protonEP(proton) {
267 generateModel(a,z);
268}
269
271 : verboseLevel(0), nuclei_radius(0.), nuclei_volume(0.), number_of_zones(0),
272 A(0), Z(0), theNucleus(0), neutronNumber(0), protonNumber(0),
273 neutronNumberCurrent(0), protonNumberCurrent(0), current_nucl1(0),
274 current_nucl2(0), gammaQDinterp(kebins),
275 crossSectionUnits(G4CascadeParameters::xsecScale()),
276 radiusUnits(G4CascadeParameters::radiusScale()),
277 skinDepth(0.611207*radiusUnits),
278 radiusScale((G4CascadeParameters::useTwoParam()?1.16:1.2)*radiusUnits),
279 radiusScale2((G4CascadeParameters::useTwoParam()?-1.3456:0.)*radiusUnits),
280 radiusForSmall(G4CascadeParameters::radiusSmall()),
281 radScaleAlpha(G4CascadeParameters::radiusAlpha()),
282 fermiMomentum(G4CascadeParameters::fermiScale()),
283 R_nucleon(G4CascadeParameters::radiusTrailing()),
284 gammaQDscale(G4CascadeParameters::gammaQDScale()),
285 potentialThickness(1.0),
286 neutronEP(neutron), protonEP(proton) {
288}
289
291 delete theNucleus;
292 theNucleus = 0;
293}
294
295
296// Initialize model state for new cascade
297
298void G4NucleiModel::reset(G4int nHitNeutrons, G4int nHitProtons,
299 const std::vector<G4ThreeVector>* hitPoints) {
300 neutronNumberCurrent = neutronNumber - nHitNeutrons;
301 protonNumberCurrent = protonNumber - nHitProtons;
302
303 // zero or copy collision point array for trailing effect
304 if (!hitPoints || !hitPoints->empty()) collisionPts.clear();
305 else collisionPts = *hitPoints;
306}
307
308
309// Generate nuclear model parameters for given nucleus
310
312 generateModel(nuclei->getA(), nuclei->getZ());
313}
314
316 if (verboseLevel) {
317 G4cout << " >>> G4NucleiModel::generateModel A " << a << " Z " << z
318 << G4endl;
319 }
320
321 // If model already built, just return; otherwise initialize everything
322 if (a == A && z == Z) {
323 if (verboseLevel > 1) G4cout << " model already generated" << z << G4endl;
324 reset();
325 return;
326 }
327
328 A = a;
329 Z = z;
330 delete theNucleus;
331 theNucleus = new G4InuclNuclei(A,Z); // For conservation checking
332
333 neutronNumber = A - Z;
334 protonNumber = Z;
335 reset();
336
337 if (verboseLevel > 3) {
338 G4cout << " crossSectionUnits = " << crossSectionUnits << G4endl
339 << " radiusUnits = " << radiusUnits << G4endl
340 << " skinDepth = " << skinDepth << G4endl
341 << " radiusScale = " << radiusScale << G4endl
342 << " radiusScale2 = " << radiusScale2 << G4endl
343 << " radiusForSmall = " << radiusForSmall << G4endl
344 << " radScaleAlpha = " << radScaleAlpha << G4endl
345 << " fermiMomentum = " << fermiMomentum << G4endl
346 << " piTimes4thirds = " << piTimes4thirds << G4endl;
347 }
348
349 G4double nuclearRadius; // Nuclear radius computed from A
350 if (A>4) nuclearRadius = radiusScale*G4cbrt(A) + radiusScale2/G4cbrt(A);
351 else nuclearRadius = radiusForSmall * (A==4 ? radScaleAlpha : 1.);
352
353 // This will be used to pre-allocate lots of arrays below
354 number_of_zones = (A < 5) ? 1 : (A < 100) ? 3 : 6;
355
356 // Clear all parameters arrays for reloading
357 binding_energies.clear();
358 nucleon_densities.clear();
359 zone_potentials.clear();
360 fermi_momenta.clear();
361 zone_radii.clear();
362 zone_volumes.clear();
363
365 fillZoneRadii(nuclearRadius);
366
367 G4double tot_vol = fillZoneVolumes(nuclearRadius); // Woods-Saxon integral
368
369 fillPotentials(proton, tot_vol); // Protons
370 fillPotentials(neutron, tot_vol); // Neutrons
371
372 // Additional flat zone potentials for other hadrons
373 const std::vector<G4double> vp(number_of_zones, (A>4)?pion_vp:pion_vp_small);
374 const std::vector<G4double> kp(number_of_zones, kaon_vp);
375 const std::vector<G4double> hp(number_of_zones, hyperon_vp);
376
377 zone_potentials.push_back(vp);
378 zone_potentials.push_back(kp);
379 zone_potentials.push_back(hp);
380
381 nuclei_radius = zone_radii.back();
382 nuclei_volume = std::accumulate(zone_volumes.begin(),zone_volumes.end(),0.);
383
384 if (verboseLevel > 3) printModel();
385}
386
387
388// Load binding energy array for current nucleus
389
391 if (verboseLevel > 1)
392 G4cout << " >>> G4NucleiModel::fillBindingEnergies" << G4endl;
393
395
396 // Binding energy differences for proton and neutron loss, respectively
397 // FIXME: Why is fabs() used here instead of the signed difference?
398 binding_energies.push_back(std::fabs(bindingEnergy(A-1,Z-1)-dm)/GeV);
399 binding_energies.push_back(std::fabs(bindingEnergy(A-1,Z)-dm)/GeV);
400}
401
402// Load zone boundary radius array for current nucleus
403
405 if (verboseLevel > 1)
406 G4cout << " >>> G4NucleiModel::fillZoneRadii" << G4endl;
407
408 G4double skinRatio = nuclearRadius/skinDepth;
409 G4double skinDecay = G4Exp(-skinRatio);
410
411 if (A < 5) { // Light ions treated as simple balls
412 zone_radii.push_back(nuclearRadius);
413 ur[0] = 0.;
414 ur[1] = 1.;
415 } else if (A < 12) { // Small nuclei have Gaussian potential
416 G4double rSq = nuclearRadius * nuclearRadius;
417 G4double gaussRadius = std::sqrt(rSq * (1.0 - 1.0/A) + 6.4);
418
419 ur[0] = 0.0;
420 for (G4int i = 0; i < number_of_zones; i++) {
421 G4double y = std::sqrt(-G4Log(alfa3[i]));
422 zone_radii.push_back(gaussRadius * y);
423 ur[i+1] = y;
424 }
425 } else if (A < 100) { // Intermediate nuclei have three-zone W-S
426 ur[0] = -skinRatio;
427 for (G4int i = 0; i < number_of_zones; i++) {
428 G4double y = G4Log((1.0 + skinDecay)/alfa3[i] - 1.0);
429 zone_radii.push_back(nuclearRadius + skinDepth * y);
430 ur[i+1] = y;
431 }
432 } else { // Heavy nuclei have six-zone W-S
433 ur[0] = -skinRatio;
434 for (G4int i = 0; i < number_of_zones; i++) {
435 G4double y = G4Log((1.0 + skinDecay)/alfa6[i] - 1.0);
436 zone_radii.push_back(nuclearRadius + skinDepth * y);
437 ur[i+1] = y;
438 }
439 }
440}
441
442// Compute zone-by-zone density integrals
443
445 if (verboseLevel > 1)
446 G4cout << " >>> G4NucleiModel::fillZoneVolumes" << G4endl;
447
448 G4double tot_vol = 0.; // Return value omitting 4pi/3 for sphere!
449
450 if (A < 5) { // Light ions are treated as simple balls
451 v[0] = v1[0] = 1.;
452 tot_vol = zone_radii[0]*zone_radii[0]*zone_radii[0];
453 zone_volumes.push_back(tot_vol*piTimes4thirds); // True volume of sphere
454
455 return tot_vol;
456 }
457
458 PotentialType usePotential = (A < 12) ? Gaussian : WoodsSaxon;
459
460 for (G4int i = 0; i < number_of_zones; i++) {
461 if (usePotential == WoodsSaxon) {
462 v[i] = zoneIntegralWoodsSaxon(ur[i], ur[i+1], nuclearRadius);
463 } else {
464 v[i] = zoneIntegralGaussian(ur[i], ur[i+1], nuclearRadius);
465 }
466
467 tot_vol += v[i];
468
469 v1[i] = zone_radii[i]*zone_radii[i]*zone_radii[i];
470 if (i > 0) v1[i] -= zone_radii[i-1]*zone_radii[i-1]*zone_radii[i-1];
471
472 zone_volumes.push_back(v1[i]*piTimes4thirds); // True volume of shell
473 }
474
475 return tot_vol; // Sum of zone integrals, not geometric volume
476}
477
478// Load potentials for nucleons (using G4InuclParticle codes)
480 if (verboseLevel > 1)
481 G4cout << " >>> G4NucleiModel::fillZoneVolumes(" << type << ")" << G4endl;
482
483 if (type != proton && type != neutron) return;
484
486
487 // FIXME: This is the fabs() binding energy difference, not signed
488 const G4double dm = binding_energies[type-1];
489
490 rod.clear(); rod.reserve(number_of_zones);
491 pf.clear(); pf.reserve(number_of_zones);
492 vz.clear(); vz.reserve(number_of_zones);
493
494 G4int nNucleons = (type==proton) ? protonNumber : neutronNumber;
495 G4double dd0 = nNucleons / tot_vol / piTimes4thirds;
496
497 for (G4int i = 0; i < number_of_zones; i++) {
498 G4double rd = dd0 * v[i] / v1[i];
499 rod.push_back(rd);
500 G4double pff = fermiMomentum * G4cbrt(rd);
501 pf.push_back(pff);
502 vz.push_back(0.5 * pff * pff / mass + dm);
503 }
504
505 nucleon_densities.push_back(rod);
506 fermi_momenta.push_back(pf);
507 zone_potentials.push_back(vz);
508}
509
510// Zone integral of Woods-Saxon density function
512 G4double nuclearRadius) const {
513 if (verboseLevel > 1) {
514 G4cout << " >>> G4NucleiModel::zoneIntegralWoodsSaxon" << G4endl;
515 }
516
517 const G4double epsilon = 1.0e-3;
518 const G4int itry_max = 1000;
519
520 G4double skinRatio = nuclearRadius / skinDepth;
521
522 G4double d2 = 2.0 * skinRatio;
523 G4double dr = r2 - r1;
524 G4double fr1 = r1 * (r1 + d2) / (1.0 + G4Exp(r1));
525 G4double fr2 = r2 * (r2 + d2) / (1.0 + G4Exp(r2));
526 G4double fi = (fr1 + fr2) / 2.;
527 G4double fun1 = fi * dr;
528 G4double fun;
529 G4int jc = 1;
530 G4double dr1 = dr;
531 G4int itry = 0;
532
533 while (itry < itry_max) { /* Loop checking 08.06.2015 MHK */
534 dr /= 2.;
535 itry++;
536
537 G4double r = r1 - dr;
538 fi = 0.0;
539
540 for (G4int i = 0; i < jc; i++) {
541 r += dr1;
542 fi += r * (r + d2) / (1.0 + G4Exp(r));
543 }
544
545 fun = 0.5 * fun1 + fi * dr;
546
547 if (std::fabs((fun - fun1) / fun) <= epsilon) break;
548
549 jc *= 2;
550 dr1 = dr;
551 fun1 = fun;
552 } // while (itry < itry_max)
553
554 if (verboseLevel > 2 && itry == itry_max)
555 G4cout << " zoneIntegralWoodsSaxon-> n iter " << itry_max << G4endl;
556
558
559 return skinDepth3 * (fun + skinRatio*skinRatio*G4Log((1.0 + G4Exp(-r1)) / (1.0 + G4Exp(-r2))));
560}
561
562
563// Zone integral of Gaussian density function
565 G4double nucRad) const {
566 if (verboseLevel > 1) {
567 G4cout << " >>> G4NucleiModel::zoneIntegralGaussian" << G4endl;
568 }
569
570 G4double gaussRadius = std::sqrt(nucRad*nucRad * (1.0 - 1.0/A) + 6.4);
571
572 const G4double epsilon = 1.0e-3;
573 const G4int itry_max = 1000;
574
575 G4double dr = r2 - r1;
576 G4double fr1 = r1 * r1 * G4Exp(-r1 * r1);
577 G4double fr2 = r2 * r2 * G4Exp(-r2 * r2);
578 G4double fi = (fr1 + fr2) / 2.;
579 G4double fun1 = fi * dr;
580 G4double fun;
581 G4int jc = 1;
582 G4double dr1 = dr;
583 G4int itry = 0;
584
585 while (itry < itry_max) { /* Loop checking 08.06.2015 MHK */
586 dr /= 2.;
587 itry++;
588 G4double r = r1 - dr;
589 fi = 0.0;
590
591 for (G4int i = 0; i < jc; i++) {
592 r += dr1;
593 fi += r * r * G4Exp(-r * r);
594 }
595
596 fun = 0.5 * fun1 + fi * dr;
597
598 if (std::fabs((fun - fun1) / fun) <= epsilon) break;
599
600 jc *= 2;
601 dr1 = dr;
602 fun1 = fun;
603 } // while (itry < itry_max)
604
605 if (verboseLevel > 2 && itry == itry_max)
606 G4cerr << " zoneIntegralGaussian-> n iter " << itry_max << G4endl;
607
608 return gaussRadius*gaussRadius*gaussRadius * fun;
609}
610
611
613 if (verboseLevel > 1) {
614 G4cout << " >>> G4NucleiModel::printModel" << G4endl;
615 }
616
617 G4cout << " nuclei model for A " << A << " Z " << Z << G4endl
618 << " proton binding energy " << binding_energies[0]
619 << " neutron binding energy " << binding_energies[1] << G4endl
620 << " Nuclei radius " << nuclei_radius << " volume " << nuclei_volume
621 << " number of zones " << number_of_zones << G4endl;
622
623 for (G4int i = 0; i < number_of_zones; i++)
624 G4cout << " zone " << i+1 << " radius " << zone_radii[i]
625 << " volume " << zone_volumes[i] << G4endl
626 << " protons: density " << getDensity(1,i) << " PF " <<
627 getFermiMomentum(1,i) << " VP " << getPotential(1,i) << G4endl
628 << " neutrons: density " << getDensity(2,i) << " PF " <<
629 getFermiMomentum(2,i) << " VP " << getPotential(2,i) << G4endl
630 << " pions: VP " << getPotential(3,i) << G4endl;
631}
632
633
635 G4double ekin = 0.0;
636
637 if (ip < 3 && izone < number_of_zones) { // ip for proton/neutron only
638 G4double pfermi = fermi_momenta[ip - 1][izone];
640
641 ekin = std::sqrt(pfermi*pfermi + mass*mass) - mass;
642 }
643 return ekin;
644}
645
646
649 G4double pmod = getFermiMomentum(type, zone) * G4cbrt(inuclRndm());
651
652 return generateWithRandomAngles(pmod, mass);
653}
654
655
658 if (verboseLevel > 1) {
659 G4cout << " >>> G4NucleiModel::generateNucleon" << G4endl;
660 }
661
663 return G4InuclElementaryParticle(mom, type);
664}
665
666
669 G4int zone) const {
670 if (verboseLevel > 1) {
671 G4cout << " >>> G4NucleiModel::generateQuasiDeuteron" << G4endl;
672 }
673
674 // Quasideuteron consists of an unbound but associated nucleon pair
675
676 // FIXME: Why generate two separate nucleon momenta (randomly!) and
677 // add them, instead of just throwing a net momentum for the
678 // dinulceon state? And why do I have to capture the two
679 // return values into local variables?
680 G4LorentzVector mom1 = generateNucleonMomentum(type1, zone);
681 G4LorentzVector mom2 = generateNucleonMomentum(type2, zone);
682 G4LorentzVector dmom = mom1+mom2;
683
684 G4int dtype = 0;
685 if (type1*type2 == pro*pro) dtype = 111;
686 else if (type1*type2 == pro*neu) dtype = 112;
687 else if (type1*type2 == neu*neu) dtype = 122;
688
689 return G4InuclElementaryParticle(dmom, dtype);
690}
691
692
693void
695 if (verboseLevel > 1) {
696 G4cout << " >>> G4NucleiModel::generateInteractionPartners" << G4endl;
697 }
698
699 thePartners.clear(); // Reset buffer for next cycle
700
701 G4int ptype = cparticle.getParticle().type();
702 G4int zone = cparticle.getCurrentZone();
703
704 G4double r_in; // Destination radius if inbound
705 G4double r_out; // Destination radius if outbound
706
707 if (zone == number_of_zones) {
708 r_in = nuclei_radius;
709 r_out = 0.0;
710 } else if (zone == 0) { // particle is outside core
711 r_in = 0.0;
712 r_out = zone_radii[0];
713 } else {
714 r_in = zone_radii[zone - 1];
715 r_out = zone_radii[zone];
716 }
717
718 G4double path = cparticle.getPathToTheNextZone(r_in, r_out);
719
720 if (verboseLevel > 2) {
721 if (isProjectile(cparticle)) G4cout << " incident particle: ";
722 G4cout << " r_in " << r_in << " r_out " << r_out << " path " << path
723 << G4endl;
724 }
725
726 if (path < -small) { // something wrong
727 if (verboseLevel)
728 G4cerr << " generateInteractionPartners-> negative path length" << G4endl;
729 return;
730 }
731
732 if (std::fabs(path) < small) { // Not moving, or just at boundary
733 if (cparticle.getMomentum().vect().mag() > small) {
734 if (verboseLevel > 3)
735 G4cout << " generateInteractionPartners-> zero path" << G4endl;
736
737 thePartners.push_back(partner()); // Dummy list terminator with zero path
738 return;
739 }
740
741 if (zone >= number_of_zones) // Place captured-at-rest in nucleus
742 zone = number_of_zones-1;
743 }
744
745 G4double invmfp = 0.; // Buffers for interaction probability
746 G4double spath = 0.;
747 for (G4int ip = 1; ip < 3; ip++) {
748 // Only process nucleons which remain active in target
749 if (ip==proton && protonNumberCurrent < 1) continue;
750 if (ip==neutron && neutronNumberCurrent < 1) continue;
751 if (ip==neutron && ptype==muonMinus) continue; // mu-/n forbidden
752
753 // All nucleons are assumed to be at rest when colliding
754 G4InuclElementaryParticle particle = generateNucleon(ip, zone);
755 invmfp = inverseMeanFreePath(cparticle, particle);
756 spath = generateInteractionLength(cparticle, path, invmfp);
757
758 if (path<small || spath < path) {
759 if (verboseLevel > 3) {
760 G4cout << " adding partner[" << thePartners.size() << "]: "
761 << particle << G4endl;
762 }
763 thePartners.push_back(partner(particle, spath));
764 }
765 } // for (G4int ip...
766
767 if (verboseLevel > 2) {
768 G4cout << " after nucleons " << thePartners.size() << " path " << path << G4endl;
769 }
770
771 // Absorption possible for pions or photons interacting with dibaryons
772 if (useQuasiDeuteron(cparticle.getParticle().type())) {
773 if (verboseLevel > 2) {
774 G4cout << " trying quasi-deuterons with bullet: "
775 << cparticle.getParticle() << G4endl;
776 }
777
778 // Initialize buffers for quasi-deuteron results
779 qdeutrons.clear();
780 acsecs.clear();
781
782 G4double tot_invmfp = 0.0; // Total inv. mean-free-path for all QDs
783
784 // Proton-proton state interacts with pi-, mu- or neutrals
785 if (protonNumberCurrent >= 2 && ptype != pip) {
787 if (verboseLevel > 2)
788 G4cout << " ptype=" << ptype << " using pp target\n" << ppd << G4endl;
789
790 invmfp = inverseMeanFreePath(cparticle, ppd);
791 tot_invmfp += invmfp;
792 acsecs.push_back(invmfp);
793 qdeutrons.push_back(ppd);
794 }
795
796 // Proton-neutron state interacts with any pion type or photon
797 if (protonNumberCurrent >= 1 && neutronNumberCurrent >= 1) {
799 if (verboseLevel > 2)
800 G4cout << " ptype=" << ptype << " using np target\n" << npd << G4endl;
801
802 invmfp = inverseMeanFreePath(cparticle, npd);
803 tot_invmfp += invmfp;
804 acsecs.push_back(invmfp);
805 qdeutrons.push_back(npd);
806 }
807
808 // Neutron-neutron state interacts with pi+ or neutrals
809 if (neutronNumberCurrent >= 2 && ptype != pim && ptype != mum) {
811 if (verboseLevel > 2)
812 G4cout << " ptype=" << ptype << " using nn target\n" << nnd << G4endl;
813
814 invmfp = inverseMeanFreePath(cparticle, nnd);
815 tot_invmfp += invmfp;
816 acsecs.push_back(invmfp);
817 qdeutrons.push_back(nnd);
818 }
819
820 // Select quasideuteron interaction from non-zero cross-section choices
821 if (verboseLevel > 2) {
822 for (size_t i=0; i<qdeutrons.size(); i++) {
823 G4cout << " acsecs[" << qdeutrons[i].getDefinition()->GetParticleName()
824 << "] " << acsecs[i];
825 }
826 G4cout << G4endl;
827 }
828
829 if (tot_invmfp > small) { // Must have absorption cross-section
830 G4double apath = generateInteractionLength(cparticle, path, tot_invmfp);
831
832 if (path<small || apath < path) { // choose the qdeutron
833 G4double sl = inuclRndm() * tot_invmfp;
834 G4double as = 0.0;
835
836 for (size_t i = 0; i < qdeutrons.size(); i++) {
837 as += acsecs[i];
838 if (sl < as) {
839 if (verboseLevel > 2)
840 G4cout << " deut type " << qdeutrons[i] << G4endl;
841
842 thePartners.push_back(partner(qdeutrons[i], apath));
843 break;
844 }
845 } // for (qdeutrons...
846 } // if (apath < path)
847 } // if (tot_invmfp > small)
848 } // if (useQuasiDeuteron) [pion, muon or photon]
849
850 if (verboseLevel > 2) {
851 G4cout << " after deuterons " << thePartners.size() << " partners"
852 << G4endl;
853 }
854
855 if (thePartners.size() > 1) { // Sort list by path length
856 std::sort(thePartners.begin(), thePartners.end(), sortPartners);
857 }
858
859 G4InuclElementaryParticle particle; // Total path at end of list
860 thePartners.push_back(partner(particle, path));
861}
862
863
866 G4ElementaryParticleCollider* theEPCollider,
867 std::vector<G4CascadParticle>& outgoing_cparticles) {
868 if (verboseLevel > 1)
869 G4cout << " >>> G4NucleiModel::generateParticleFate" << G4endl;
870
871 if (verboseLevel > 2) G4cout << " cparticle: " << cparticle << G4endl;
872
873 // Create four-vector checking
874#ifdef G4CASCADE_CHECK_ECONS
875 G4CascadeCheckBalance balance(0.005, 0.01, "G4NucleiModel"); // Second arg is in GeV
877#endif
878
879 outgoing_cparticles.clear(); // Clear return buffer for this event
880 generateInteractionPartners(cparticle); // Fills "thePartners" data
881
882 if (thePartners.empty()) { // smth. is wrong -> needs special treatment
883 if (verboseLevel)
884 G4cerr << " generateParticleFate-> got empty interaction-partners list "
885 << G4endl;
886 return;
887 }
888
889 G4int npart = thePartners.size(); // Last item is a total-path placeholder
890
891 if (npart == 1) { // cparticle is on the next zone entry
892 if (verboseLevel > 1)
893 G4cout << " no interactions; moving to next zone" << G4endl;
894
897 boundaryTransition(cparticle);
898 outgoing_cparticles.push_back(cparticle);
899
900 if (verboseLevel > 2) G4cout << " next zone \n" << cparticle << G4endl;
901
902 //A.R. 19-Jun-2013: Fixed rare cases of non-reproducibility.
903 current_nucl1 = 0;
904 current_nucl2 = 0;
905
906 return;
907 } // if (npart == 1)
908
909 if (verboseLevel > 1)
910 G4cout << " processing " << npart-1 << " possible interactions" << G4endl;
911
912 G4ThreeVector old_position = cparticle.getPosition();
913 G4InuclElementaryParticle& bullet = cparticle.getParticle();
914 G4bool no_interaction = true;
915 G4int zone = cparticle.getCurrentZone();
916
917 for (G4int i=0; i<npart-1; i++) { // Last item is a total-path placeholder
918 if (i > 0) cparticle.updatePosition(old_position);
919
920 G4InuclElementaryParticle& target = thePartners[i].first;
921
922 if (verboseLevel > 3) {
923 if (target.quasi_deutron()) G4cout << " try absorption: ";
924 G4cout << " target " << target.type() << " bullet " << bullet.type()
925 << G4endl;
926 }
927
929 // Pass current (A,Z) configuration for possible recoils
931 theEPCollider->setNucleusState(massNumberCurrent, protonNumberCurrent);
932 theEPCollider->collide(&bullet, &target, EPCoutput);
933
934 // If collision failed, exit loop over partners
935 if (EPCoutput.numberOfOutgoingParticles() == 0) break;
936
937 if (verboseLevel > 2) {
939#ifdef G4CASCADE_CHECK_ECONS
940 balance.collide(&bullet, &target, EPCoutput);
941 balance.okay(); // Do checks, but ignore result
942#endif
943 }
944
945 // Get list of outgoing particles for evaluation
946 std::vector<G4InuclElementaryParticle>& outgoing_particles =
948
949 if (!passFermi(outgoing_particles, zone)) continue; // Interaction fails
950
951 // Trailing effect: reject interaction at previously hit nucleon
953 const G4ThreeVector& new_position = cparticle.getPosition();
954
955 if (!passTrailing(new_position)) continue;
956 collisionPts.push_back(new_position);
957
958 // Sort particles according to beta (fastest first)
959 std::sort(outgoing_particles.begin(), outgoing_particles.end(),
961
962 if (verboseLevel > 2)
963 G4cout << " adding " << outgoing_particles.size()
964 << " output particles" << G4endl;
965
966 // NOTE: Embedded temporary is optimized away (no copying gets done)
967 G4int nextGen = cparticle.getGeneration()+1;
968 for (G4int ip = 0; ip < G4int(outgoing_particles.size()); ip++) {
969 outgoing_cparticles.push_back(G4CascadParticle(outgoing_particles[ip],
970 new_position, zone,
971 0.0, nextGen));
972 }
973
974 no_interaction = false;
975 current_nucl1 = 0;
976 current_nucl2 = 0;
977
978#ifdef G4CASCADE_DEBUG_CHARGE
979 {
980 G4double out_charge = 0.0;
981
982 for (G4int ip = 0; ip < G4int(outgoing_particles.size()); ip++)
983 out_charge += outgoing_particles[ip].getCharge();
984
985 G4cout << " multiplicity " << outgoing_particles.size()
986 << " bul type " << bullet.type()
987 << " targ type " << target.type()
988 << "\n initial charge "
989 << bullet.getCharge() + target.getCharge()
990 << " out charge " << out_charge << G4endl;
991 }
992#endif
993
994 if (verboseLevel > 2)
995 G4cout << " partner type " << target.type() << G4endl;
996
997 if (target.nucleon()) {
998 current_nucl1 = target.type();
999 } else {
1000 if (verboseLevel > 2) G4cout << " good absorption " << G4endl;
1001 current_nucl1 = (target.type() - 100) / 10;
1002 current_nucl2 = target.type() - 100 - 10 * current_nucl1;
1003 }
1004
1005 if (current_nucl1 == 1) {
1006 if (verboseLevel > 3) G4cout << " decrement proton count" << G4endl;
1008 } else {
1009 if (verboseLevel > 3) G4cout << " decrement neutron count" << G4endl;
1011 }
1012
1013 if (current_nucl2 == 1) {
1014 if (verboseLevel > 3) G4cout << " decrement proton count" << G4endl;
1016 } else if (current_nucl2 == 2) {
1017 if (verboseLevel > 3) G4cout << " decrement neutron count" << G4endl;
1019 }
1020
1021 break;
1022 } // loop over partners
1023
1024 if (no_interaction) { // still no interactions
1025 if (verboseLevel > 1) G4cout << " no interaction " << G4endl;
1026
1027 // For conservation checking (below), get particle before updating
1028 static G4ThreadLocal G4InuclElementaryParticle *prescatCP_G4MT_TLS_ = 0;
1029 if (!prescatCP_G4MT_TLS_) {
1030 prescatCP_G4MT_TLS_ = new G4InuclElementaryParticle;
1031 G4AutoDelete::Register(prescatCP_G4MT_TLS_);
1032 }
1033 G4InuclElementaryParticle &prescatCP = *prescatCP_G4MT_TLS_; // Avoid memory churn
1034 prescatCP = cparticle.getParticle();
1035
1036 // Last "partner" is just a total-path placeholder
1037 cparticle.updatePosition(old_position);
1038 cparticle.propagateAlongThePath(thePartners[npart-1].second);
1039 cparticle.incrementCurrentPath(thePartners[npart-1].second);
1040 boundaryTransition(cparticle);
1041 outgoing_cparticles.push_back(cparticle);
1042
1043 // Check conservation for simple scattering (ignore target nucleus!)
1044#ifdef G4CASCADE_CHECK_ECONS
1045 if (verboseLevel > 2) {
1046 balance.collide(&prescatCP, 0, outgoing_cparticles);
1047 balance.okay(); // Report violations, but don't act on them
1048 }
1049#endif
1050 } // if (no_interaction)
1051
1052 return;
1053}
1054
1055// Test if particle is suitable for absorption with dibaryons
1056
1058 if (qdtype==pn || qdtype==0) // All absorptive particles
1059 return (ptype==pi0 || ptype==pip || ptype==pim || ptype==gam || ptype==mum);
1060 else if (qdtype==pp) // Negative or neutral only
1061 return (ptype==pi0 || ptype==pim || ptype==gam || ptype==mum);
1062 else if (qdtype==nn) // Positive or neutral only
1063 return (ptype==pi0 || ptype==pip || ptype==gam);
1064
1065 return false; // Not a quasideuteron target
1066}
1067
1068
1069G4bool G4NucleiModel::passFermi(const std::vector<G4InuclElementaryParticle>& particles,
1070 G4int zone) {
1071 if (verboseLevel > 1) {
1072 G4cout << " >>> G4NucleiModel::passFermi" << G4endl;
1073 }
1074
1075 // Only check Fermi momenta for nucleons
1076 for (G4int i = 0; i < G4int(particles.size()); i++) {
1077 if (!particles[i].nucleon()) continue;
1078
1079 G4int type = particles[i].type();
1080 G4double mom = particles[i].getMomModule();
1081 G4double pfermi = fermi_momenta[type-1][zone];
1082
1083 if (verboseLevel > 2)
1084 G4cout << " type " << type << " p " << mom << " pf " << pfermi << G4endl;
1085
1086 if (mom < pfermi) {
1087 if (verboseLevel > 2) G4cout << " rejected by Fermi" << G4endl;
1088 return false;
1089 }
1090 }
1091 return true;
1092}
1093
1094
1095// Test here for trailing effect: loop over all previous collision
1096// locations and test for d > R_nucleon
1097
1099 if (verboseLevel > 1)
1100 G4cout << " >>> G4NucleiModel::passTrailing " << hit_position << G4endl;
1101
1102 G4double dist;
1103 for (G4int i = 0; i < G4int(collisionPts.size() ); i++) {
1104 dist = (collisionPts[i] - hit_position).mag();
1105 if (verboseLevel > 2) G4cout << " dist " << dist << G4endl;
1106 if (dist < R_nucleon) {
1107 if (verboseLevel > 2) G4cout << " rejected by Trailing" << G4endl;
1108 return false;
1109 }
1110 }
1111 return true; // New point far enough away to be used
1112}
1113
1114
1116 if (verboseLevel > 1) {
1117 G4cout << " >>> G4NucleiModel::boundaryTransition" << G4endl;
1118 }
1119
1120 G4int zone = cparticle.getCurrentZone();
1121
1122 if (cparticle.movingInsideNuclei() && zone == 0) {
1123 if (verboseLevel) G4cerr << " boundaryTransition-> in zone 0 " << G4endl;
1124 return;
1125 }
1126
1127 G4LorentzVector mom = cparticle.getMomentum();
1128 G4ThreeVector pos = cparticle.getPosition();
1129
1130 G4int type = cparticle.getParticle().type();
1131
1132 G4double r = pos.mag();
1133 G4double p = mom.vect().mag(); // NAT
1134 G4double pr = pos.dot(mom.vect()) / r;
1135 G4double pperp2 = p*p - pr*pr; // NAT
1136
1137 G4int next_zone = cparticle.movingInsideNuclei() ? zone - 1 : zone + 1;
1138
1139 // NOTE: dv is the height of the wall seen by the particle
1140 G4double dv = getPotential(type,next_zone) - getPotential(type,zone);
1141 if (verboseLevel > 3) {
1142 G4cout << "Potentials for type " << type << " = "
1143 << getPotential(type,zone) << " , "
1144 << getPotential(type,next_zone) << G4endl;
1145 }
1146
1147 G4double qv = dv*dv + 2.0*dv*mom.e() + pr*pr;
1148 // G4double qv = dv*dv + 2.0*dv*mom.m() + pr*pr; // more correct? NAT
1149 G4double p1r = 0.;
1150
1151 // Perpendicular contribution to pr^2 after penetrating // NAT
1152 // potential, to leading order in thickness // NAT
1153 G4double qperp = 2.0*pperp2*potentialThickness/r; // NAT
1154
1155 if (verboseLevel > 3) {
1156 G4cout << " type " << type << " zone " << zone << " next " << next_zone
1157 << " qv " << qv << " dv " << dv << G4endl;
1158 }
1159
1160 G4bool adjustpperp = false; // NAT
1161 G4double smallish = 0.001; // NAT
1162
1163// if (qv <= 0.0) { // reflection
1164 if (qv <= 0.0 && qv+qperp <=0 ) { // reflection // NAT
1165 if (verboseLevel > 3) G4cout << " reflects off boundary" << G4endl;
1166 p1r = -pr;
1167 cparticle.incrementReflectionCounter();
1168// } else { // transition
1169
1170 } else if (qv > 0.0) { // standard transition // NAT
1171 if (verboseLevel > 3) G4cout << " passes thru boundary" << G4endl;
1172 p1r = std::sqrt(qv);
1173 if (pr < 0.0) p1r = -p1r;
1174 cparticle.updateZone(next_zone);
1175 cparticle.resetReflection();
1176
1177 } else { // transition via transverse kinetic energy (allowed for thick walls) // NAT
1178 if (verboseLevel > 3) G4cout << " passes thru boundary due to angular momentum" << G4endl;
1179 p1r = smallish * pr; // don't want exactly tangent momentum
1180 adjustpperp = true;
1181
1182 cparticle.updateZone(next_zone);
1183 cparticle.resetReflection();
1184 }
1185
1186 G4double prr = (p1r - pr)/r; // change to radial momentum, divided by r
1187
1188 if (verboseLevel > 3) {
1189 G4cout << " prr " << prr << " delta px " << prr*pos.x() << " py "
1190 << prr*pos.y() << " pz " << prr*pos.z() << " mag "
1191 << std::fabs(prr*r) << G4endl;
1192 }
1193
1194 if (adjustpperp) { // NAT
1195 G4ThreeVector old_pperp = mom.vect() - pos*(pr/r);
1196 G4double new_pperp_mag = std::sqrt(std::max(0.0, pperp2 + qv - p1r*p1r) );
1197 // new total momentum found by rescaling p_perp
1198 mom.setVect(old_pperp * new_pperp_mag/std::sqrt(pperp2));
1199 // add a small radial component to make sure that we propagate into new zone
1200 mom.setVect(mom.vect() + pos*p1r/r);
1201 } else {
1202 mom.setVect(mom.vect() + pos*prr);
1203 }
1204
1205 cparticle.updateParticleMomentum(mom);
1206}
1207
1208
1209// Select random point along full trajectory through nucleus
1210// NOTE: Intended for projectile photons for initial interaction
1211
1213 if (verboseLevel > 1)
1214 G4cout << " >>> G4NucleiModel::choosePointAlongTraj" << G4endl;
1215
1216 // Get trajectory through nucleus by computing exit point of line,
1217 // assuming that current position is on surface
1218
1219 // FIXME: These need to be reusable (static) buffers
1220 G4ThreeVector pos = cparticle.getPosition();
1221 G4ThreeVector rhat = pos.unit();
1222
1223 G4ThreeVector phat = cparticle.getMomentum().vect().unit();
1224 if (cparticle.getMomentum().vect().mag() < small) phat.set(0.,0.,1.);
1225
1226 if (verboseLevel > 3)
1227 G4cout << " pos " << pos << " phat " << phat << " rhat " << rhat << G4endl;
1228
1229 G4ThreeVector posout = pos;
1230 G4double prang = rhat.angle(-phat);
1231
1232 if (prang < 1e-6) posout = -pos; // Check for radial incidence
1233 else {
1234 G4double posrot = 2.*prang - pi;
1235 posout.rotate(posrot, phat.cross(rhat));
1236 if (verboseLevel > 3) G4cout << " posrot " << posrot/deg << " deg";
1237 }
1238
1239 if (verboseLevel > 3) G4cout << " posout " << posout << G4endl;
1240
1241 // Get list of zone crossings along trajectory
1242 G4ThreeVector posmid = (pos+posout)/2.; // Midpoint of trajectory
1243 G4double r2mid = posmid.mag2();
1244 G4double lenmid = (posout-pos).mag()/2.; // Half-length of trajectory
1245
1246 G4int zoneout = number_of_zones-1;
1247 G4int zonemid = getZone(posmid.mag()); // Middle zone traversed
1248
1249 // Every zone is entered then exited, so preallocate vector
1250 G4int ncross = (number_of_zones-zonemid)*2;
1251
1252 if (verboseLevel > 3) {
1253 G4cout << " posmid " << posmid << " lenmid " << lenmid
1254 << " zoneout " << zoneout << " zonemid " << zonemid
1255 << " ncross " << ncross << G4endl;
1256 }
1257
1258 // FIXME: These need to be reusable (static) buffers
1259 std::vector<G4double> wtlen(ncross,0.); // CDF from entry point
1260 std::vector<G4double> len(ncross,0.); // Distance from entry point
1261
1262 // Work from outside in, to accumulate CDF steps properly
1263 G4int i; // Loop variable, used multiple times
1264 for (i=0; i<ncross/2; i++) {
1265 G4int iz = zoneout-i;
1266 G4double ds = std::sqrt(zone_radii[iz]*zone_radii[iz]-r2mid);
1267
1268 len[i] = lenmid - ds; // Distance from entry to crossing
1269 len[ncross-1-i] = lenmid + ds; // Distance to outbound crossing
1270
1271 if (verboseLevel > 3) {
1272 G4cout << " i " << i << " iz " << iz << " ds " << ds
1273 << " len " << len[i] << G4endl;
1274 }
1275 }
1276
1277 // Compute weights for each zone along trajectory and integrate
1278 for (i=1; i<ncross; i++) {
1279 G4int iz = (i<ncross/2) ? zoneout-i+1 : zoneout-ncross+i+1;
1280
1281 G4double dlen = len[i]-len[i-1]; // Distance from previous crossing
1282
1283 // Uniform probability across entire zone
1284 //*** G4double wt = dlen*(getDensity(1,iz)+getDensity(2,iz));
1285
1286 // Probability based on interaction length through zone
1287 G4double invmfp = (inverseMeanFreePath(cparticle, neutronEP, iz)
1288 + inverseMeanFreePath(cparticle, protonEP, iz));
1289
1290 // Integral of exp(-len/mfp) from start of zone to end
1291 G4double wt = (G4Exp(-len[i-1]*invmfp)-G4Exp(-len[i]*invmfp)) / invmfp;
1292
1293 wtlen[i] = wtlen[i-1] + wt;
1294
1295 if (verboseLevel > 3) {
1296 G4cout << " i " << i << " iz " << iz << " avg.mfp " << 1./invmfp
1297 << " dlen " << dlen << " wt " << wt << " wtlen " << wtlen[i]
1298 << G4endl;
1299 }
1300 }
1301
1302 // Normalize CDF to unit integral
1303 std::transform(wtlen.begin(), wtlen.end(), wtlen.begin(),
1304 std::bind(std::divides<G4double>(), std::placeholders::_1, wtlen.back()));
1305
1306 if (verboseLevel > 3) {
1307 G4cout << " weights";
1308 for (i=0; i<ncross; i++) G4cout << " " << wtlen[i];
1309 G4cout << G4endl;
1310 }
1311
1312 // Choose random point along trajectory, weighted by density
1313 G4double rand = G4UniformRand();
1314 G4int ir = std::upper_bound(wtlen.begin(),wtlen.end(),rand) - wtlen.begin();
1315
1316 G4double frac = (rand-wtlen[ir-1]) / (wtlen[ir]-wtlen[ir-1]);
1317 G4double drand = (1.-frac)*len[ir-1] + frac*len[ir];
1318
1319 if (verboseLevel > 3) {
1320 G4cout << " rand " << rand << " ir " << ir << " frac " << frac
1321 << " drand " << drand << G4endl;
1322 }
1323
1324 pos += drand * phat; // Distance from entry point along trajectory
1325
1326 // Update input particle with new location
1327 cparticle.updatePosition(pos);
1328 cparticle.updateZone(getZone(pos.mag()));
1329
1330 if (verboseLevel > 2) {
1331 G4cout << " moved particle to zone " << cparticle.getCurrentZone()
1332 << " @ " << pos << G4endl;
1333 }
1334}
1335
1336
1337// Returns true if particle should interact immediately
1339 return (isProjectile(cparticle) &&
1340 (cparticle.getParticle().isPhoton() ||
1341 cparticle.getParticle().isMuon())
1342 );
1343}
1344
1346 return (cparticle.getGeneration() == 0); // Only initial state particles
1347}
1348
1350 if (verboseLevel > 1) {
1351 G4cout << " >>> G4NucleiModel::worthToPropagate" << G4endl;
1352 }
1353
1354 const G4double ekin_scale = 2.0;
1355
1356 G4bool worth = true;
1357
1358 if (cparticle.reflectedNow()) { // Just reflected -- keep going?
1359 G4int zone = cparticle.getCurrentZone();
1360 G4int ip = cparticle.getParticle().type();
1361
1362 // NOTE: Temporarily backing out use of potential for non-nucleons
1363 G4double ekin_cut = (cparticle.getParticle().nucleon()) ?
1364 getFermiKinetic(ip, zone) : 0.; //*** getPotential(ip, zone);
1365
1366 worth = cparticle.getParticle().getKineticEnergy()/ekin_scale > ekin_cut;
1367
1368 if (verboseLevel > 3) {
1369 G4cout << " type=" << ip
1370 << " ekin=" << cparticle.getParticle().getKineticEnergy()
1371 << " potential=" << ekin_cut
1372 << " : worth? " << worth << G4endl;
1373 }
1374 }
1375
1376 return worth;
1377}
1378
1379
1381 if (verboseLevel > 4) {
1382 G4cout << " >>> G4NucleiModel::getRatio " << ip << G4endl;
1383 }
1384
1385 switch (ip) {
1388 case diproton: return getRatio(proton)*getRatio(proton);
1389 case unboundPN: return getRatio(proton)*getRatio(neutron);
1390 case dineutron: return getRatio(neutron)*getRatio(neutron);
1391 default: return 0.;
1392 }
1393
1394 return 0.;
1395}
1396
1398 const G4double pn_spec = 1.0; // Scale factor for pn vs. pp/nn
1399 //const G4double pn_spec = 0.5;
1400
1401 G4double dens = 0.;
1402
1403 if (ip < 100) dens = getDensity(ip,izone);
1404 else { // For dibaryons, remove extra 1/volume term in density product
1405 switch (ip) {
1406 case diproton:
1407 dens = getDensity(proton,izone) * getDensity(proton,izone);
1408 break;
1409 case unboundPN:
1410 dens = getDensity(proton,izone) * getDensity(neutron,izone) * pn_spec;
1411 break;
1412 case dineutron:
1413 dens = getDensity(neutron,izone) * getDensity(neutron,izone);
1414 break;
1415 default: dens = 0.;
1416 }
1417 dens *= getVolume(izone);
1418 }
1419
1420 return getRatio(ip) * dens;
1421}
1422
1423
1426 if (verboseLevel > 1) {
1427 G4cout << " >>> G4NucleiModel::initializeCascad(particle)" << G4endl;
1428 }
1429
1430 // FIXME: Previous version generated random sin(theta), then used -cos(theta)
1431 // Using generateWithRandomAngles changes result!
1432 // G4ThreeVector pos = generateWithRandomAngles(nuclei_radius).vect();
1433 G4double costh = std::sqrt(1.0 - inuclRndm());
1435
1436 // Start particle outside nucleus, unless capture-at-rest
1437 G4int zone = number_of_zones;
1438 if (particle->getKineticEnergy() < small) zone--;
1439
1440 G4CascadParticle cpart(*particle, pos, zone, large, 0);
1441
1442 // SPECIAL CASE: Inbound photons are emplanted along through-path
1443 if (forceFirst(cpart)) choosePointAlongTraj(cpart);
1444
1445 if (verboseLevel > 2) G4cout << cpart << G4endl;
1446
1447 return cpart;
1448}
1449
1451 G4InuclNuclei* target,
1452 modelLists& output) {
1453 if (verboseLevel) {
1454 G4cout << " >>> G4NucleiModel::initializeCascad(bullet,target,output)"
1455 << G4endl;
1456 }
1457
1458 const G4double max_a_for_cascad = 5.0;
1459 const G4double ekin_cut = 2.0;
1460 const G4double small_ekin = 1.0e-6;
1461 const G4double r_large2for3 = 62.0;
1462 const G4double r0forAeq3 = 3.92;
1463 const G4double s3max = 6.5;
1464 const G4double r_large2for4 = 69.14;
1465 const G4double r0forAeq4 = 4.16;
1466 const G4double s4max = 7.0;
1467 const G4int itry_max = 100;
1468
1469 // Convenient references to input buffer contents
1470 std::vector<G4CascadParticle>& casparticles = output.first;
1471 std::vector<G4InuclElementaryParticle>& particles = output.second;
1472
1473 casparticles.clear(); // Reset input buffer to fill this event
1474 particles.clear();
1475
1476 // first decide whether it will be cascad or compound final nuclei
1477 G4int ab = bullet->getA();
1478 G4int zb = bullet->getZ();
1479 G4int at = target->getA();
1480 G4int zt = target->getZ();
1481
1482 G4double massb = bullet->getMass(); // For creating LorentzVectors below
1483
1484 if (ab < max_a_for_cascad) {
1485
1486 G4double benb = bindingEnergy(ab,zb)/GeV / G4double(ab);
1487 G4double bent = bindingEnergy(at,zt)/GeV / G4double(at);
1488 G4double ben = benb < bent ? bent : benb;
1489
1490 if (bullet->getKineticEnergy()/ab > ekin_cut*ben) {
1491 G4int itryg = 0;
1492
1493 /* Loop checking 08.06.2015 MHK */
1494 while (casparticles.size() == 0 && itryg < itry_max) {
1495 itryg++;
1496 particles.clear();
1497
1498 // nucleons coordinates and momenta in nuclei rest frame
1499 coordinates.clear();
1500 momentums.clear();
1501
1502 if (ab < 3) { // deuteron, simplest case
1503 G4double r = 2.214 - 3.4208 * G4Log(1.0 - 0.981 * inuclRndm());
1505 coordinates.push_back(coord1);
1506 coordinates.push_back(-coord1);
1507
1508 G4double p = 0.0;
1509 G4bool bad = true;
1510 G4int itry = 0;
1511
1512 while (bad && itry < itry_max) { /* Loop checking 08.06.2015 MHK */
1513 itry++;
1514 p = 456.0 * inuclRndm();
1515
1516 if (p * p / (p * p + 2079.36) / (p * p + 2079.36) > 1.2023e-4 * inuclRndm() &&
1517 p * r > 312.0) bad = false;
1518 }
1519
1520 if (itry == itry_max)
1521 if (verboseLevel > 2){
1522 G4cout << " deutron bullet generation-> itry = " << itry_max << G4endl;
1523 }
1524
1525 p = 0.0005 * p;
1526
1527 if (verboseLevel > 2){
1528 G4cout << " p nuc " << p << G4endl;
1529 }
1530
1532
1533 momentums.push_back(mom);
1534 mom.setVect(-mom.vect());
1535 momentums.push_back(-mom);
1536 } else {
1537 G4ThreeVector coord1;
1538
1539 G4bool badco = true;
1540
1541 G4int itry = 0;
1542
1543 if (ab == 3) {
1544 while (badco && itry < itry_max) {/* Loop checking 08.06.2015 MHK */
1545 if (itry > 0) coordinates.clear();
1546 itry++;
1547 G4int i(0);
1548
1549 for (i = 0; i < 2; i++) {
1550 G4int itry1 = 0;
1551 G4double ss, u, rho;
1552 G4double fmax = G4Exp(-0.5) / std::sqrt(0.5);
1553
1554 while (itry1 < itry_max) { /* Loop checking 08.06.2015 MHK */
1555 itry1++;
1556 ss = -G4Log(inuclRndm());
1557 u = fmax * inuclRndm();
1558 rho = std::sqrt(ss) * G4Exp(-ss);
1559
1560 if (rho > u && ss < s3max) {
1561 ss = r0forAeq3 * std::sqrt(ss);
1562 coord1 = generateWithRandomAngles(ss).vect();
1563 coordinates.push_back(coord1);
1564
1565 if (verboseLevel > 2){
1566 G4cout << " i " << i << " r " << coord1.mag() << G4endl;
1567 }
1568 break;
1569 }
1570 }
1571
1572 if (itry1 == itry_max) { // bad case
1573 coord1.set(10000.,10000.,10000.);
1574 coordinates.push_back(coord1);
1575 break;
1576 }
1577 }
1578
1579 coord1 = -coordinates[0] - coordinates[1];
1580 if (verboseLevel > 2) {
1581 G4cout << " 3 r " << coord1.mag() << G4endl;
1582 }
1583
1584 coordinates.push_back(coord1);
1585
1586 G4bool large_dist = false;
1587
1588 for (i = 0; i < 2; i++) {
1589 for (G4int j = i+1; j < 3; j++) {
1590 G4double r2 = (coordinates[i]-coordinates[j]).mag2();
1591
1592 if (verboseLevel > 2) {
1593 G4cout << " i " << i << " j " << j << " r2 " << r2 << G4endl;
1594 }
1595
1596 if (r2 > r_large2for3) {
1597 large_dist = true;
1598
1599 break;
1600 }
1601 }
1602
1603 if (large_dist) break;
1604 }
1605
1606 if(!large_dist) badco = false;
1607
1608 }
1609
1610 } else { // a >= 4
1611 G4double b = 3./(ab - 2.0);
1612 G4double b1 = 1.0 - b / 2.0;
1613 G4double u = b1 + std::sqrt(b1 * b1 + b);
1614 G4double fmax = (1.0 + u/b) * u * G4Exp(-u);
1615
1616 while (badco && itry < itry_max) {/* Loop checking 08.06.2015 MHK */
1617
1618 if (itry > 0) coordinates.clear();
1619 itry++;
1620 G4int i(0);
1621
1622 for (i = 0; i < ab-1; i++) {
1623 G4int itry1 = 0;
1624 G4double ss;
1625
1626 while (itry1 < itry_max) { /* Loop checking 08.06.2015 MHK */
1627 itry1++;
1628 ss = -G4Log(inuclRndm());
1629 u = fmax * inuclRndm();
1630
1631 if (std::sqrt(ss) * G4Exp(-ss) * (1.0 + ss/b) > u
1632 && ss < s4max) {
1633 ss = r0forAeq4 * std::sqrt(ss);
1634 coord1 = generateWithRandomAngles(ss).vect();
1635 coordinates.push_back(coord1);
1636
1637 if (verboseLevel > 2) {
1638 G4cout << " i " << i << " r " << coord1.mag() << G4endl;
1639 }
1640
1641 break;
1642 }
1643 }
1644
1645 if (itry1 == itry_max) { // bad case
1646 coord1.set(10000.,10000.,10000.);
1647 coordinates.push_back(coord1);
1648 break;
1649 }
1650 }
1651
1652 coord1 *= 0.0; // Cheap way to reset
1653 for(G4int j = 0; j < ab -1; j++) coord1 -= coordinates[j];
1654
1655 coordinates.push_back(coord1);
1656
1657 if (verboseLevel > 2){
1658 G4cout << " last r " << coord1.mag() << G4endl;
1659 }
1660
1661 G4bool large_dist = false;
1662
1663 for (i = 0; i < ab-1; i++) {
1664 for (G4int j = i+1; j < ab; j++) {
1665
1666 G4double r2 = (coordinates[i]-coordinates[j]).mag2();
1667
1668 if (verboseLevel > 2){
1669 G4cout << " i " << i << " j " << j << " r2 " << r2 << G4endl;
1670 }
1671
1672 if (r2 > r_large2for4) {
1673 large_dist = true;
1674
1675 break;
1676 }
1677 }
1678
1679 if (large_dist) break;
1680 }
1681
1682 if (!large_dist) badco = false;
1683 }
1684 }
1685
1686 if(badco) {
1687 G4cout << " can not generate the nucleons coordinates for a "
1688 << ab << G4endl;
1689
1690 casparticles.clear(); // Return empty buffer on error
1691 particles.clear();
1692 return;
1693
1694 } else { // momentums
1695 G4double p, u, x;
1696 G4LorentzVector mom;
1697 //G4bool badp = True;
1698
1699 for (G4int i = 0; i < ab - 1; i++) {
1700 G4int itry2 = 0;
1701
1702 while(itry2 < itry_max) { /* Loop checking 08.06.2015 MHK */
1703 itry2++;
1704 u = -G4Log(0.879853 - 0.8798502 * inuclRndm());
1705 x = u * G4Exp(-u);
1706
1707 if(x > inuclRndm()) {
1708 p = std::sqrt(0.01953 * u);
1709 mom = generateWithRandomAngles(p, massb);
1710 momentums.push_back(mom);
1711
1712 break;
1713 }
1714 } // while (itry2 < itry_max)
1715
1716 if(itry2 == itry_max) {
1717 G4cout << " can not generate proper momentum for a "
1718 << ab << G4endl;
1719
1720 casparticles.clear(); // Return empty buffer on error
1721 particles.clear();
1722 return;
1723 }
1724 } // for (i=0 ...
1725 // last momentum
1726
1727 mom *= 0.; // Cheap way to reset
1728 mom.setE(bullet->getEnergy()+target->getEnergy());
1729
1730 for(G4int j=0; j< ab-1; j++) mom -= momentums[j];
1731
1732 momentums.push_back(mom);
1733 }
1734 }
1735
1736 // Coordinates and momenta at rest are generated, now back to the lab
1737 G4double rb = 0.0;
1738 G4int i(0);
1739
1740 for(i = 0; i < G4int(coordinates.size()); i++) {
1741 G4double rp = coordinates[i].mag2();
1742
1743 if(rp > rb) rb = rp;
1744 }
1745
1746 // nuclei i.p. as a whole
1747 G4double s1 = std::sqrt(inuclRndm());
1748 G4double phi = randomPHI();
1749 G4double rz = (nuclei_radius + rb) * s1;
1750 G4ThreeVector global_pos(rz*std::cos(phi), rz*std::sin(phi),
1751 -(nuclei_radius+rb)*std::sqrt(1.0-s1*s1));
1752
1753 for (i = 0; i < G4int(coordinates.size()); i++) {
1754 coordinates[i] += global_pos;
1755 }
1756
1757 // all nucleons at rest
1758 raw_particles.clear();
1759
1760 for (G4int ipa = 0; ipa < ab; ipa++) {
1761 G4int knd = ipa < zb ? 1 : 2;
1763 }
1764
1765 G4InuclElementaryParticle dummy(small_ekin, 1);
1766 G4LorentzConvertor toTheBulletRestFrame(&dummy, bullet);
1767 toTheBulletRestFrame.toTheTargetRestFrame();
1768
1769 particleIterator ipart;
1770
1771 for (ipart = raw_particles.begin(); ipart != raw_particles.end(); ipart++) {
1772 ipart->setMomentum(toTheBulletRestFrame.backToTheLab(ipart->getMomentum()));
1773 }
1774
1775 // fill cascad particles and outgoing particles
1776
1777 for(G4int ip = 0; ip < G4int(raw_particles.size()); ip++) {
1778 G4LorentzVector mom = raw_particles[ip].getMomentum();
1779 G4double pmod = mom.rho();
1780 G4double t0 = -mom.vect().dot(coordinates[ip]) / pmod;
1781 G4double det = t0 * t0 + nuclei_radius * nuclei_radius
1782 - coordinates[ip].mag2();
1783 G4double tr = -1.0;
1784
1785 if(det > 0.0) {
1786 G4double t1 = t0 + std::sqrt(det);
1787 G4double t2 = t0 - std::sqrt(det);
1788
1789 if(std::fabs(t1) <= std::fabs(t2)) {
1790 if(t1 > 0.0) {
1791 if(coordinates[ip].z() + mom.z() * t1 / pmod <= 0.0) tr = t1;
1792 }
1793
1794 if(tr < 0.0 && t2 > 0.0) {
1795
1796 if(coordinates[ip].z() + mom.z() * t2 / pmod <= 0.0) tr = t2;
1797 }
1798
1799 } else {
1800 if(t2 > 0.0) {
1801
1802 if(coordinates[ip].z() + mom.z() * t2 / pmod <= 0.0) tr = t2;
1803 }
1804
1805 if(tr < 0.0 && t1 > 0.0) {
1806 if(coordinates[ip].z() + mom.z() * t1 / pmod <= 0.0) tr = t1;
1807 }
1808 }
1809
1810 }
1811
1812 if(tr >= 0.0) { // cascad particle
1813 coordinates[ip] += mom.vect()*tr / pmod;
1814 casparticles.push_back(G4CascadParticle(raw_particles[ip],
1815 coordinates[ip],
1816 number_of_zones, large, 0));
1817
1818 } else {
1819 particles.push_back(raw_particles[ip]);
1820 }
1821 }
1822 }
1823
1824 if(casparticles.size() == 0) {
1825 particles.clear();
1826
1827 G4cout << " can not generate proper distribution for " << itry_max
1828 << " steps " << G4endl;
1829 }
1830 }
1831 }
1832
1833 if(verboseLevel > 2){
1834 G4int ip(0);
1835
1836 G4cout << " cascad particles: " << casparticles.size() << G4endl;
1837 for(ip = 0; ip < G4int(casparticles.size()); ip++)
1838 G4cout << casparticles[ip] << G4endl;
1839
1840 G4cout << " outgoing particles: " << particles.size() << G4endl;
1841 for(ip = 0; ip < G4int(particles.size()); ip++)
1842 G4cout << particles[ip] << G4endl;
1843 }
1844
1845 return; // Buffer has been filled
1846}
1847
1848
1849// Compute mean free path for inclusive interaction of projectile and target
1850G4double
1852 const G4InuclElementaryParticle& target,
1853 G4int zone) {
1854 G4int ptype = cparticle.getParticle().type();
1855 G4int ip = target.type();
1856
1857 // Ensure that zone specified is within nucleus, for array lookups
1858 if (zone<0) zone = cparticle.getCurrentZone();
1859 if (zone>=number_of_zones) zone = number_of_zones-1;
1860
1861 // Special cases: neutrinos, and muon-on-neutron, have infinite path
1862 if (cparticle.getParticle().isNeutrino()) return 0.;
1863 if (ptype == muonMinus && ip == neutron) return 0.;
1864
1865 // Configure frame transformation to get kinetic energy for lookups
1867 dummy_convertor.setTarget(&target);
1868 dummy_convertor.toTheCenterOfMass(); // Fill internal kinematics
1870
1871 // Get cross section for interacting with target (dibaryons are absorptive)
1872 G4double csec = (ip < 100) ? totalCrossSection(ekin, ptype*ip)
1873 : absorptionCrossSection(ekin, ptype);
1874
1875 if (verboseLevel > 2) {
1876 G4cout << " ip " << ip << " zone " << zone << " ekin " << ekin
1877 << " dens " << getCurrentDensity(ip, zone)
1878 << " csec " << csec << G4endl;
1879 }
1880
1881 if (csec <= 0.) return 0.; // No interaction, avoid unnecessary work
1882
1883 // Get nuclear density and compute mean free path
1884 return csec * getCurrentDensity(ip, zone);
1885}
1886
1887// Throw random distance for interaction of particle using path and MFP
1888
1889G4double
1891 G4double path, G4double invmfp) const {
1892 // Delay interactions of newly formed secondaries (minimum int. length)
1893 const G4double young_cut = std::sqrt(10.0) * 0.25;
1894 const G4double huge_num = 50.0; // Argument to exponential
1895
1896 G4double spath = large; // Buffer for return value
1897
1898 if (invmfp < small) return spath; // No interaction, avoid unnecessary work
1899
1900 G4double pw = -path * invmfp; // Ratio of path in zone to MFP
1901 if (pw < -huge_num) pw = -huge_num;
1902 pw = 1.0 - G4Exp(pw);
1903
1904 if (verboseLevel > 2)
1905 G4cout << " mfp " << 1./invmfp << " pw " << pw << G4endl;
1906
1907 // Primary particle(s) should always interact at least once
1908 if (forceFirst(cparticle) || (inuclRndm() < pw)) {
1909 spath = -G4Log(1.0 - pw * inuclRndm()) / invmfp;
1910 if (cparticle.young(young_cut, spath)) spath = large;
1911
1912 if (verboseLevel > 2)
1913 G4cout << " spath " << spath << " path " << path << G4endl;
1914 }
1915
1916 return spath;
1917}
1918
1919
1920// Parametrized cross section for pion and photon absorption
1921
1923 if (!useQuasiDeuteron(type)) {
1924 G4cerr << "absorptionCrossSection() only valid for incident pions or gammas"
1925 << G4endl;
1926 return 0.;
1927 }
1928
1929 G4double csec = 0.;
1930
1931 // Pion absorption is parametrized for low vs. medium energy
1932 // ... use for muon capture as well
1933 if (type == pionPlus || type == pionMinus || type == pionZero ||
1934 type == muonMinus) {
1935 if (ke < 0.3) csec = (0.1106 / std::sqrt(ke) - 0.8
1936 + 0.08 / ((ke-0.123)*(ke-0.123) + 0.0056) );
1937 else if (ke < 1.0) csec = 3.6735 * (1.0-ke)*(1.0-ke);
1938 }
1939
1940 if (type == photon) {
1942 }
1943
1944 if (csec < 0.0) csec = 0.0;
1945
1946 if (verboseLevel > 2) {
1947 G4cout << " ekin " << ke << " abs. csec " << csec << " mb" << G4endl;
1948 }
1949
1950 return crossSectionUnits * csec;
1951}
1952
1953
1955{
1956 // All scattering cross-sections are available from tables
1957 const G4CascadeChannel* xsecTable = G4CascadeChannelTables::GetTable(rtype);
1958 if (!xsecTable) {
1959 G4cerr << " unknown collison type = " << rtype << G4endl;
1960 return 0.;
1961 }
1962
1963 return (crossSectionUnits * xsecTable->getCrossSection(ke));
1964}
std::vector< G4InuclElementaryParticle >::iterator particleIterator
Definition: G4BigBanger.cc:64
G4double epsilon(G4double density, G4double temperature)
static const G4double pos
static const G4double ab
static const G4double d2
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double G4Log(G4double x)
Definition: G4Log.hh:226
std::vector< G4InuclElementaryParticle >::iterator particleIterator
static constexpr double second
Definition: G4SIunits.hh:137
static constexpr double GeV
Definition: G4SIunits.hh:203
static constexpr double pi
Definition: G4SIunits.hh:55
static constexpr double deg
Definition: G4SIunits.hh:132
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
double mag2() const
Hep3Vector cross(const Hep3Vector &) const
double angle(const Hep3Vector &) const
double dot(const Hep3Vector &) const
double mag() const
void set(double x, double y, double z)
Hep3Vector & rotate(double, const Hep3Vector &)
Definition: ThreeVectorR.cc:24
Hep3Vector vect() const
void setVect(const Hep3Vector &)
G4bool reflectedNow() const
G4int getGeneration() 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)
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
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)
const G4InuclElementaryParticle protonEP
G4int getZone(G4double r) const
const G4double R_nucleon
static const G4double large
G4bool forceFirst(const G4CascadParticle &cparticle) const
G4double getRatio(G4int ip) const
const G4double skinDepth
G4double nuclei_radius
std::vector< std::vector< G4double > > zone_potentials
std::vector< G4double > binding_energies
const G4double crossSectionUnits
G4double getDensity(G4int ip, G4int izone) const
G4CascadParticle initializeCascad(G4InuclElementaryParticle *particle)
G4bool worthToPropagate(const G4CascadParticle &cparticle) const
G4double generateInteractionLength(const G4CascadParticle &cparticle, G4double path, G4double invmfp) const
G4double v[6]
const G4double radiusScale2
static const G4double small
void generateParticleFate(G4CascadParticle &cparticle, G4ElementaryParticleCollider *theEPCollider, std::vector< G4CascadParticle > &cascade)
const G4double radiusUnits
G4double inverseMeanFreePath(const G4CascadParticle &cparticle, const G4InuclElementaryParticle &target, G4int zone=-1)
void printModel() const
std::vector< G4ThreeVector > collisionPts
void fillBindingEnergies()
std::vector< G4double > acsecs
const G4double radScaleAlpha
G4CascadeInterpolator< 30 > gammaQDinterp
G4double getCurrentDensity(G4int ip, G4int izone) const
std::vector< G4InuclElementaryParticle > qdeutrons
std::pair< G4InuclElementaryParticle, G4double > partner
G4double fillZoneVolumes(G4double nuclearRadius)
static const G4double alfa3[3]
std::vector< G4LorentzVector > momentums
G4InuclElementaryParticle generateNucleon(G4int type, G4int zone) const
std::vector< G4double > pf
void boundaryTransition(G4CascadParticle &cparticle)
G4double getVolume(G4int izone) const
G4double totalCrossSection(G4double ke, G4int rtype) const
G4double getPotential(G4int ip, G4int izone) const
G4CollisionOutput EPCoutput
G4int protonNumberCurrent
virtual ~G4NucleiModel()
static const G4double piTimes4thirds
G4double getFermiMomentum(G4int ip, G4int izone) const
const G4double radiusForSmall
G4bool passTrailing(const G4ThreeVector &hit_position)
static G4bool useQuasiDeuteron(G4int ptype, G4int qdtype=0)
G4double v1[6]
std::vector< std::vector< G4double > > fermi_momenta
static const G4double alfa6[6]
G4InuclNuclei * theNucleus
void reset(G4int nHitNeutrons=0, G4int nHitProtons=0, const std::vector< G4ThreeVector > *hitPoints=0)
G4LorentzConvertor dummy_convertor
G4double nuclei_volume
G4double absorptionCrossSection(G4double e, G4int type) const
const G4double fermiMomentum
void generateModel(G4InuclNuclei *nuclei)
G4InuclElementaryParticle generateQuasiDeuteron(G4int type1, G4int type2, G4int zone) const
std::vector< G4double > vz
G4double zoneIntegralGaussian(G4double ur1, G4double ur2, G4double nuclearRadius) const
std::vector< G4double > zone_volumes
static const G4double pion_vp
static const G4double hyperon_vp
std::vector< G4InuclElementaryParticle > raw_particles
std::vector< G4ThreeVector > coordinates
G4bool passFermi(const std::vector< G4InuclElementaryParticle > &particles, G4int zone)
std::pair< std::vector< G4CascadParticle >, std::vector< G4InuclElementaryParticle > > modelLists
G4double getFermiKinetic(G4int ip, G4int izone) const
std::vector< partner > thePartners
static const G4double pion_vp_small
std::vector< G4double > zone_radii
G4double zoneIntegralWoodsSaxon(G4double ur1, G4double ur2, G4double nuclearRadius) const
void fillZoneRadii(G4double nuclearRadius)
G4LorentzVector generateNucleonMomentum(G4int type, G4int zone) const
G4bool isProjectile(const G4CascadParticle &cparticle) const
static G4bool sortPartners(const partner &p1, const partner &p2)
void fillPotentials(G4int type, G4double tot_vol)
std::vector< G4double > rod
const G4double gammaQDscale
void generateInteractionPartners(G4CascadParticle &cparticle)
G4double ur[7]
const G4double radiusScale
std::vector< std::vector< G4double > > nucleon_densities
const G4InuclElementaryParticle neutronEP
G4int neutronNumberCurrent
const G4double potentialThickness
void choosePointAlongTraj(G4CascadParticle &cparticle)
static const G4double kaon_vp
virtual void setVerboseLevel(G4int verbose=0)
void Register(T *inst)
Definition: G4AutoDelete.hh:65
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4bool nucleon(G4int ityp)
G4LorentzVector generateWithFixedTheta(G4double ct, G4double p, G4double mass=0.)
G4double bindingEnergy(G4int A, G4int Z)
G4LorentzVector generateWithRandomAngles(G4double p, G4double mass=0.)
G4bool transform(G4String &input, const G4String &type)
#define G4ThreadLocal
Definition: tls.hh:77