Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4INCLCascade.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 // INCL++ intra-nuclear cascade model
27 // Pekka Kaitaniemi, CEA and Helsinki Institute of Physics
28 // Davide Mancusi, CEA
29 // Alain Boudard, CEA
30 // Sylvie Leray, CEA
31 // Joseph Cugnon, University of Liege
32 //
33 #define INCLXX_IN_GEANT4_MODE 1
34 
35 #include "globals.hh"
36 
41 #include "G4INCLCascade.hh"
42 #include "G4INCLRandom.hh"
43 #include "G4INCLRanecu.hh"
44 #include "G4INCLGeant4Random.hh"
46 #include "G4INCLParticleTable.hh"
47 #include "G4INCLGlobalInfo.hh"
48 
49 #include "G4INCLPauliBlocking.hh"
50 #include "G4INCLIPauli.hh"
51 #include "G4INCLPauliStrict.hh"
52 #include "G4INCLPauliStandard.hh"
54 #include "G4INCLPauliGlobal.hh"
55 #include "G4INCLCDPP.hh"
56 
57 #include "G4INCLLogger.hh"
58 #include "G4INCLGlobals.hh"
60 
62 #include "G4INCLICoulomb.hh"
63 #include "G4INCLCoulombNone.hh"
65 
66 #include "G4INCLClustering.hh"
69 
70 #include "G4INCLIntersection.hh"
71 
72 #include "G4INCLCrossSections.hh"
73 
74 #include <cstring>
75 #include <cstdlib>
76 #include <numeric>
77 
78 namespace G4INCL {
79 
81  :propagationModel(0), theA(208), theZ(82),
82  targetInitSuccess(false),
83  maxImpactParameter(0.),
84  maxUniverseRadius(0.),
85  maxInteractionDistance(0.),
86  fixedImpactParameter(0.),
87  theConfig(config),
88  nucleus(NULL),
89  minRemnantSize(4)
90  {
91  // Set the logger object.
94 
95  // Set the random number generator algorithm. The system can support
96  // multiple different generator algorithms in a completely
97  // transparent way.
98 #ifdef INCLXX_IN_GEANT4_MODE
100 #else
102 #endif // INCLXX_IN_GEANT4_MODE
103 
104  // Select the Pauli blocking algorithm:
105  G4INCL::PauliType pauli = theConfig->getPauliType();
106  if(pauli == G4INCL::StrictStatisticalPauli)
108  else if(pauli == G4INCL::StatisticalPauli)
110  else if(pauli == G4INCL::StrictPauli)
112  else if(pauli == G4INCL::GlobalPauli)
114  else if(pauli == G4INCL::NoPauli)
116 
117  if(theConfig->getCDPP())
119  else
121 
122  // Select the Coulomb-distortion algorithm:
123  G4INCL::CoulombType coulombType = theConfig->getCoulombType();
124  if(coulombType == G4INCL::NonRelativisticCoulomb)
126  else // if(coulombType == G4INCL::NoCoulomb)
128 
129  // Select the clustering algorithm:
130  G4INCL::ClusterAlgorithmType clusterAlgorithm = theConfig->getClusterAlgorithm();
131  if(clusterAlgorithm == G4INCL::IntercomparisonClusterAlgorithm)
133  else // if(clusterAlgorithm == G4INCL::NoClusterAlgorithm)
135 
136  // Initialize the INCL particle table:
138 
139  // Propagation model is responsible for finding avatars and
140  // transporting the particles. In principle this step is "hidden"
141  // behind an abstract interface and the rest of the system does not
142  // care how the transportation and avatar finding is done. This
143  // should allow us to "easily" experiment with different avatar
144  // finding schemes and even to support things like curved
145  // trajectories in the future.
146  propagationModel = new G4INCL::StandardPropagationModel(theConfig->getLocalEnergyBBType(),theConfig->getLocalEnergyPiType());
147  eventAction = new EventAction();
148  propagationAction = new PropagationAction();
149  avatarAction = new AvatarAction();
150 
151  theGlobalInfo.cascadeModel = theConfig->getVersionID().c_str();
152  theGlobalInfo.deexcitationModel = "none";
153 
154 #ifndef INCLXX_IN_GEANT4_MODE
155  // Fill in the global information
156  theGlobalInfo.At = theConfig->getTargetA();
157  theGlobalInfo.Zt = theConfig->getTargetZ();
158  const ParticleSpecies theSpecies = theConfig->getProjectileSpecies();
159  theGlobalInfo.Ap = theSpecies.theA;
160  theGlobalInfo.Zp = theSpecies.theZ;
161  theGlobalInfo.Ep = theConfig->getProjectileKineticEnergy();
162  // Echo the input parameters to the log file
163  INFO(theConfig->echo() << std::endl);
164 #endif
165 
166  fixedImpactParameter = theConfig->getImpactParameter();
167  }
168 
176  delete avatarAction;
177  delete propagationAction;
178  delete eventAction;
179  delete propagationModel;
180  delete theConfig;
181  }
182 
183  G4bool INCL::prepareReaction(const ParticleSpecies &projectileSpecies, const G4double kineticEnergy, const G4int A, const G4int Z) {
184  if(A < 0 || A > 300 || Z < 1 || Z > 200) {
185  ERROR("Unsupported target: A = " << A << " Z = " << Z << std::endl);
186  ERROR("Target configuration rejected." << std::endl);
187  return false;
188  }
189 
190  // Initialise the maximum universe radius
191  initUniverseRadius(projectileSpecies, kineticEnergy, A, Z);
192 
193  // Initialise the nucleus
194  theZ = Z;
195  if(theConfig->isNaturalTarget())
197  else
198  theA = A;
199  initializeTarget(theA, theZ);
200 
201  // Set the maximum impact parameter
202  maxImpactParameter = CoulombDistortion::maxImpactParameter(projectileSpecies, kineticEnergy, nucleus);
203  initMaxInteractionDistance(projectileSpecies, kineticEnergy); // for forced CN events
204 
205  // Set the geometric cross section
206  theGlobalInfo.geometricCrossSection =
207  Math::tenPi*std::pow(maxImpactParameter,2);
208 
209  // Set the minimum remnant size
210  if(projectileSpecies.theA > 0)
211  minRemnantSize = std::min(theA, 4);
212  else
213  minRemnantSize = std::min(theA-1, 4);
214 
215  return true;
216  }
217 
219  delete nucleus;
220 
221  nucleus = new Nucleus(A, Z, theConfig, maxUniverseRadius);
222  nucleus->getStore()->getBook()->reset();
223  nucleus->initializeParticles();
224 
225  propagationModel->setNucleus(nucleus);
226  return true;
227  }
228 
230  ParticleSpecies const &projectileSpecies,
231  const G4double kineticEnergy,
232  const G4int targetA,
233  const G4int targetZ
234  ) {
235  // Set the target and the projectile
236  targetInitSuccess = prepareReaction(projectileSpecies, kineticEnergy, targetA, targetZ);
237 
238  if(!targetInitSuccess) {
239  WARN("Target initialisation failed for A=" << targetA << ", Z=" << targetZ << std::endl);
240  theEventInfo.transparent=true;
241  return theEventInfo;
242  }
243 
244  const G4bool canRunCascade = preCascade(projectileSpecies, kineticEnergy);
245  if(canRunCascade) {
246  cascade();
247  postCascade();
248  }
249  return theEventInfo;
250  }
251 
252  G4bool INCL::preCascade(ParticleSpecies const projectileSpecies, const G4double kineticEnergy) {
253  // Reset theEventInfo
254  theEventInfo.reset();
255 
257 
258  // Increment the global counter for the number of shots
259  theGlobalInfo.nShots++;
260 
261  // Fill in the event information
262  theEventInfo.projectileType = projectileSpecies.theType;
263  theEventInfo.Ap = projectileSpecies.theA;
264  theEventInfo.Zp = projectileSpecies.theZ;
265  theEventInfo.Ep = kineticEnergy;
266  theEventInfo.At = nucleus->getA();
267  theEventInfo.Zt = nucleus->getZ();
268 
269  // Do nothing below the Coulomb barrier
270  if(maxImpactParameter<=0.) {
271  // Increment the global counter for the number of transparents
272  theGlobalInfo.nTransparents++;
273 
274  // Fill in the event information
275  theEventInfo.transparent = true;
276 
277  return false;
278  }
279 
280  // Randomly draw an impact parameter or use a fixed value, depending on the
281  // Config option
282  G4double impactParameter, phi;
283  if(fixedImpactParameter<0.) {
284  impactParameter = maxImpactParameter * std::sqrt(Random::shoot0());
285  phi = Random::shoot() * Math::twoPi;
286  } else {
287  impactParameter = fixedImpactParameter;
288  phi = 0.;
289  }
290 
291  // Fill in the event information
292  theEventInfo.impactParameter = impactParameter;
293 
294  const G4double effectiveImpactParameter = propagationModel->shoot(projectileSpecies, kineticEnergy, impactParameter, phi);
295  if(effectiveImpactParameter < 0.) {
296  // Increment the global counter for the number of transparents
297  theGlobalInfo.nTransparents++;
298 
299  // Fill in the event information
300  theEventInfo.transparent = true;
301 
302  return false;
303  }
304 
305  // Fill in the event information
306  theEventInfo.transparent = false;
307  theEventInfo.effectiveImpactParameter = effectiveImpactParameter;
308 
309  return true;
310  }
311 
312  void INCL::cascade() {
313  do {
314  // Run book keeping actions that should take place before propagation:
315  propagationAction->beforePropagationAction(propagationModel);
316 
317  // Get the avatar with the smallest time and propagate particles
318  // to that point in time.
319  G4INCL::IAvatar *avatar = propagationModel->propagate();
320 
321  // Run book keeping actions that should take place after propagation:
322  propagationAction->afterPropagationAction(propagationModel, avatar);
323 
324  if(avatar == 0) break; // No more avatars in the avatar list.
325 
326  // Run book keeping actions that should take place before avatar:
327  avatarAction->beforeAvatarAction(avatar, nucleus);
328 
329  // Channel is responsible for calculating the outcome of the
330  // selected avatar. There are different kinds of channels. The
331  // class IChannel is, again, an abstract interface that defines
332  // the externally observable behavior of all interaction
333  // channels.
334  // The handling of the channel is transparent to the API.
335  // Final state tells what changed...
336  G4INCL::FinalState *finalState = avatar->getFinalState();
337 
338  // Run book keeping actions that should take place after avatar:
339  avatarAction->afterAvatarAction(avatar, nucleus, finalState);
340 
341  // So now we must give this information to the nucleus
342  nucleus->applyFinalState(finalState);
343  // and now we are ready to process the next avatar!
344 
345  delete avatar;
346  delete finalState;
347  } while(continueCascade());
348 
349  }
350 
351  void INCL::postCascade() {
352  // Fill in the event information
353  theEventInfo.stoppingTime = propagationModel->getCurrentTime();
354 
355  // Forced CN?
356  if(nucleus->getTryCompoundNucleus()) {
357  DEBUG("Trying compound nucleus" << std::endl);
358  makeCompoundNucleus();
359  // Global checks of conservation laws
360 #ifndef INCLXX_IN_GEANT4_MODE
361  if(!theEventInfo.transparent) globalConservationChecks(true);
362 #endif
363  return;
364  }
365 
366  theEventInfo.transparent = nucleus->isEventTransparent();
367 
368  if(theEventInfo.transparent) {
369  // Increment the global counter for the number of transparents
370  theGlobalInfo.nTransparents++;
371  if(nucleus->isForcedTransparent())
372  theGlobalInfo.nForcedTransparents++;
373  ProjectileRemnant * const projectileRemnant = nucleus->getProjectileRemnant();
374  if(projectileRemnant) {
375  // Delete the projectile remnant and the particles it contains
376  projectileRemnant->deleteParticles();
377  nucleus->deleteProjectileRemnant();
378  nucleus->getStore()->clearIncoming();
379  } else {
380  // Clean up the incoming list and force a transparent gracefully
381  nucleus->getStore()->deleteIncoming();
382  }
383  } else {
384  // Check if the nucleus contains deltas
385  theEventInfo.deltasInside = nucleus->containsDeltas();
386 
387  // Take care of any remaining deltas
388  theEventInfo.forcedDeltasOutside = nucleus->decayOutgoingDeltas();
389  theEventInfo.forcedDeltasInside = nucleus->decayInsideDeltas();
390 
391  // Apply Coulomb distortion, if appropriate
392  // Note that this will apply Coulomb distortion also on pions emitted by
393  // unphysical remnants (see decayInsideDeltas). This is at variance with
394  // what INCL4.6 does, but these events are (should be!) so rare that
395  // whatever we do doesn't (shouldn't!) make any noticeable difference.
397 
398  // If the normal cascade predicted complete fusion, use the tabulated
399  // masses to compute the excitation energy, the recoil, etc.
400  if(nucleus->getStore()->getOutgoingParticles().size()==0
401  && nucleus->getProjectileRemnant()
402  && nucleus->getProjectileRemnant()->getParticles().size()==0) {
403 
404  DEBUG("Cascade resulted in complete fusion, using realistic fusion kinematics" << std::endl);
405 
406  nucleus->useFusionKinematics();
407  nucleus->deleteProjectileRemnant();
408 
409  if(nucleus->getExcitationEnergy()<0.) {
410  // Complete fusion is energetically impossible, return a transparent
411  WARN("Complete-fusion kinematics yields negative excitation energy, returning a transparent!" << std::endl);
412  theEventInfo.transparent = true;
413  return;
414  }
415 
416  } else { // Normal cascade here
417 
418  // Set the excitation energy
419  nucleus->setExcitationEnergy(nucleus->computeExcitationEnergy());
420 
421  // Make a projectile pre-fragment out of the geometrical and dynamical
422  // spectators
423  theEventInfo.nUnmergedSpectators = makeProjectileRemnant();
424 
425  // Compute recoil momentum, energy and spin of the nucleus
426  nucleus->computeRecoilKinematics();
427 
428 #ifndef INCLXX_IN_GEANT4_MODE
429  // Global checks of conservation laws
430  globalConservationChecks(false);
431 #endif
432 
433  // Make room for the remnant recoil by rescaling the energies of the
434  // outgoing particles.
435  if(nucleus->hasRemnant()) rescaleOutgoingForRecoil();
436 
437  }
438 
439  // Cluster decay
440  theEventInfo.clusterDecay = nucleus->decayOutgoingClusters() | nucleus->decayMe();
441 
442 #ifndef INCLXX_IN_GEANT4_MODE
443  // Global checks of conservation laws
444  globalConservationChecks(true);
445 #endif
446 
447  // Fill the EventInfo structure
448  nucleus->fillEventInfo(&theEventInfo);
449 
450  // Check if we have an absorption:
451  if(theEventInfo.nucleonAbsorption) theGlobalInfo.nNucleonAbsorptions++;
452  if(theEventInfo.pionAbsorption) theGlobalInfo.nPionAbsorptions++;
453  }
454  }
455 
456  void INCL::makeCompoundNucleus() {
457  // Reset the internal Nucleus variables
458  nucleus->getStore()->clearIncoming();
459  nucleus->getStore()->clearOutgoing();
460  nucleus->getProjectileRemnant()->reset();
461  nucleus->setA(theEventInfo.At);
462  nucleus->setZ(theEventInfo.Zt);
463 
464  // CN kinematical variables
465  // Note: the CN orbital angular momentum is neglected in what follows. We
466  // should actually take it into account!
467  ThreeVector theCNMomentum = nucleus->getIncomingMomentum();
468  ThreeVector theCNSpin = nucleus->getIncomingAngularMomentum();
469  const G4double theTargetMass = ParticleTable::getTableMass(theEventInfo.At, theEventInfo.Zt);
470  G4int theCNA=theEventInfo.At, theCNZ=theEventInfo.Zt;
471  Cluster * const theProjectileRemnant = nucleus->getProjectileRemnant();
472  G4double theCNEnergy = theTargetMass + theProjectileRemnant->getEnergy();
473 
474  // Loop over the potential participants
475  ParticleList initialProjectileComponents = theProjectileRemnant->getParticles();
476  std::vector<Particle *> shuffledComponents(initialProjectileComponents.begin(), initialProjectileComponents.end());
477  // Shuffle the list of potential participants
478  std::random_shuffle(shuffledComponents.begin(), shuffledComponents.end(), shuffleComponentsHelper);
479 
480  G4bool success = true;
481  G4bool atLeastOneNucleonEntering = false;
482  for(std::vector<Particle*>::const_iterator p=shuffledComponents.begin(); p!=shuffledComponents.end(); ++p) {
483  // Skip geometrical spectators
485  (*p)->getPosition(),
486  (*p)->getPropagationVelocity(),
487  maxUniverseRadius));
488  if(!intersectionUniverse.exists)
489  continue;
490 
491  // At least one nucleon must enter the interaction distance
493  (*p)->getPosition(),
494  (*p)->getPropagationVelocity(),
495  maxInteractionDistance));
496  if(intersectionInteraction.exists)
497  atLeastOneNucleonEntering = true;
498 
499  // Build an entry avatar for this nucleon
500  ParticleEntryAvatar theAvatar(0.0, nucleus, *p);
501  FinalState *fs = theAvatar.getFinalState();
502  nucleus->applyFinalState(fs);
503  FinalStateValidity validity = fs->getValidity();
504  delete fs;
505  switch(validity) {
506  case ValidFS:
508  // Add the particle to the CN
509  theCNA++;
510  theCNZ += (*p)->getZ();
511  break;
512  case ParticleBelowZeroFS:
513  case PauliBlockedFS:
515  default:
516  success = false;
517  break;
518  }
519  }
520 
521  if(!success || !atLeastOneNucleonEntering) {
522  DEBUG("No nucleon entering in forced CN, or some nucleons entering below zero, forcing a transparent" << std::endl);
523  theEventInfo.transparent = true;
524  theGlobalInfo.nTransparents++;
525  theGlobalInfo.nForcedTransparents++;
527  nucleus->deleteProjectileRemnant();
528  return;
529  }
530 
531 // assert(theCNA==nucleus->getA());
532 // assert(theCNA>theEventInfo.At);
533 
534  // Update the kinematics of the CN
535  theCNEnergy -= theProjectileRemnant->getEnergy();
536  theCNMomentum -= theProjectileRemnant->getMomentum();
537 
538  // Deal with the projectile remnant
539  nucleus->finalizeProjectileRemnant(propagationModel->getCurrentTime());
540 
541  // Subtract the angular momentum of the projectile remnant
542  ParticleList const &outgoing = nucleus->getStore()->getOutgoingParticles();
543 // assert(outgoing.size()==0 || outgoing.size()==1);
544  for(ParticleIter i=outgoing.begin(); i!=outgoing.end(); ++i) {
545  theCNSpin -= (*i)->getAngularMomentum();
546  }
547 
548  // Compute the excitation energy of the CN
549  const G4double theCNMass = ParticleTable::getTableMass(theCNA,theCNZ);
550  const G4double theCNInvariantMassSquared = theCNEnergy*theCNEnergy-theCNMomentum.mag2();
551  if(theCNInvariantMassSquared<0.) {
552  // Negative invariant mass squared, return a transparent
553  theGlobalInfo.nTransparents++;
554  theGlobalInfo.nForcedTransparents++;
555  theEventInfo.transparent = true;
556  return;
557  }
558  const G4double theCNExcitationEnergy = std::sqrt(theCNInvariantMassSquared) - theCNMass;
559  if(theCNExcitationEnergy<0.) {
560  // Negative excitation energy, return a transparent
561  theGlobalInfo.nTransparents++;
562  theGlobalInfo.nForcedTransparents++;
563  theEventInfo.transparent = true;
564  return;
565  } else {
566  // Positive excitation energy, can make a CN
567  nucleus->setA(theCNA);
568  nucleus->setZ(theCNZ);
569  nucleus->setMomentum(theCNMomentum);
570  nucleus->setEnergy(theCNEnergy);
571  nucleus->setExcitationEnergy(theCNExcitationEnergy);
572  nucleus->setMass(theCNMass+theCNExcitationEnergy);
573  nucleus->setSpin(theCNSpin); // neglects any orbital angular momentum of the CN
574 
575  // Take care of any remaining deltas
576  theEventInfo.forcedDeltasOutside = nucleus->decayOutgoingDeltas();
577 
578  // Cluster decay
579  theEventInfo.clusterDecay = nucleus->decayOutgoingClusters() | nucleus->decayMe();
580 
581  // Fill the EventInfo structure
582  nucleus->fillEventInfo(&theEventInfo);
583  theGlobalInfo.nForcedCompoundNucleus++;
584  }
585  }
586 
587  void INCL::rescaleOutgoingForRecoil() {
588  RecoilCMFunctor theRecoilFunctor(nucleus, theEventInfo);
589 
590  // Apply the root-finding algorithm
591  const G4bool success = RootFinder::solve(&theRecoilFunctor, 1.0);
592  if(success) {
593  std::pair<G4double,G4double> theSolution = RootFinder::getSolution();
594  theRecoilFunctor(theSolution.first); // Apply the solution
595  } else {
596  WARN("Couldn't accommodate remnant recoil while satisfying energy conservation, root-finding algorithm failed." << std::endl);
597  }
598 
599  }
600 
601 #ifndef INCLXX_IN_GEANT4_MODE
602  void INCL::globalConservationChecks(G4bool afterRecoil) {
603  Nucleus::ConservationBalance theBalance = nucleus->getConservationBalance(theEventInfo,afterRecoil);
604 
605  // Global conservation checks
606  const G4double pLongBalance = theBalance.momentum.getZ();
607  const G4double pTransBalance = theBalance.momentum.perp();
608  if(theBalance.Z != 0) {
609  ERROR("Violation of charge conservation! ZBalance = " << theBalance.Z << std::endl);
610  }
611  if(theBalance.A != 0) {
612  ERROR("Violation of baryon-number conservation! ABalance = " << theBalance.A << std::endl);
613  }
614  G4double EThreshold, pLongThreshold, pTransThreshold;
615  if(afterRecoil) {
616  // Less stringent checks after accommodating recoil
617  EThreshold = 10.; // MeV
618  pLongThreshold = 1.; // MeV/c
619  pTransThreshold = 1.; // MeV/c
620  } else {
621  // More stringent checks before accommodating recoil
622  EThreshold = 0.1; // MeV
623  pLongThreshold = 0.1; // MeV/c
624  pTransThreshold = 0.1; // MeV/c
625  }
626  if(std::abs(theBalance.energy)>EThreshold) {
627  WARN("Violation of energy conservation > " << EThreshold << " MeV. EBalance = " << theBalance.energy << " afterRecoil = " << afterRecoil << " eventNumber=" << theEventInfo.eventNumber << std::endl);
628  }
629  if(std::abs(pLongBalance)>pLongThreshold) {
630  WARN("Violation of longitudinal momentum conservation > " << pLongThreshold << " MeV/c. pLongBalance = " << pLongBalance << " afterRecoil = " << afterRecoil << " eventNumber=" << theEventInfo.eventNumber << std::endl);
631  }
632  if(std::abs(pTransBalance)>pTransThreshold) {
633  WARN("Violation of transverse momentum conservation > " << pTransThreshold << " MeV/c. pTransBalance = " << pTransBalance << " afterRecoil = " << afterRecoil << " eventNumber=" << theEventInfo.eventNumber << std::endl);
634  }
635 
636  // Feed the EventInfo variables
637  theEventInfo.EBalance = theBalance.energy;
638  theEventInfo.pLongBalance = pLongBalance;
639  theEventInfo.pTransBalance = pTransBalance;
640  }
641 #endif
642 
643  G4bool INCL::continueCascade() {
644  // Stop if we have passed the stopping time
645  if(propagationModel->getCurrentTime() > propagationModel->getStoppingTime()) {
646  DEBUG("Cascade time (" << propagationModel->getCurrentTime()
647  << ") exceeded stopping time (" << propagationModel->getStoppingTime()
648  << "), stopping cascade" << std::endl);
649  return false;
650  }
651  // Stop if there are no participants and no pions inside the nucleus
652  if(nucleus->getStore()->getBook()->getCascading()==0 &&
653  nucleus->getStore()->getIncomingParticles().empty()) {
654  DEBUG("No participants in the nucleus and no incoming particles left, stopping cascade" << std::endl);
655  return false;
656  }
657  // Stop if the remnant is smaller than minRemnantSize
658  if(nucleus->getA() <= minRemnantSize) {
659  DEBUG("Remnant size (" << nucleus->getA()
660  << ") smaller than or equal to minimum (" << minRemnantSize
661  << "), stopping cascade" << std::endl);
662  return false;
663  }
664  // Stop if we have to try and make a compound nucleus or if we have to
665  // force a transparent
666  if(nucleus->getTryCompoundNucleus()) {
667  DEBUG("Trying to make a compound nucleus, stopping cascade" << std::endl);
668  return false;
669  }
670  if(nucleus->isForcedTransparent()) {
671  DEBUG("Forcing a transparent, stopping cascade" << std::endl);
672  return false;
673  }
674 
675  return true;
676  }
677 
679  theGlobalInfo.nucleonAbsorptionCrossSection = theGlobalInfo.geometricCrossSection *
680  ((G4double) theGlobalInfo.nNucleonAbsorptions) / ((G4double) theGlobalInfo.nShots);
681  theGlobalInfo.pionAbsorptionCrossSection = theGlobalInfo.geometricCrossSection *
682  ((G4double) theGlobalInfo.nPionAbsorptions) / ((G4double) theGlobalInfo.nShots);
683  theGlobalInfo.reactionCrossSection = theGlobalInfo.geometricCrossSection *
684  ((G4double) (theGlobalInfo.nShots - theGlobalInfo.nTransparents)) /
685  ((G4double) theGlobalInfo.nShots);
686  theGlobalInfo.errorReactionCrossSection = theGlobalInfo.geometricCrossSection *
687  std::sqrt((G4double) (theGlobalInfo.nShots - theGlobalInfo.nTransparents)) /
688  ((G4double) theGlobalInfo.nShots);
689  }
690 
691  G4int INCL::makeProjectileRemnant() {
692  G4int nUnmergedSpectators = 0;
693 
694  // Do nothing if this is not a nucleus-nucleus reaction
695  if(!nucleus->getProjectileRemnant())
696  return 0;
697 
698  // Get the spectators (geometrical+dynamical) from the Store
699  ParticleList geomSpectators(nucleus->getProjectileRemnant()->getParticles());
700  ParticleList dynSpectators(nucleus->getStore()->extractDynamicalSpectators());
701 
702  // If there are no spectators, do nothing
703  if(dynSpectators.empty() && geomSpectators.empty()) {
704  nucleus->deleteProjectileRemnant();
705  return 0;
706  } else if(geomSpectators.size()+dynSpectators.size()==1) {
707  if(dynSpectators.empty()) {
708  // No dynamical spectators, one geometrical spectator
709  // It should already be on shell.
710 #if !defined(NDEBUG) && !defined(INCLXX_IN_GEANT4_MODE)
711  Particle *theSpectator = geomSpectators.front();
712 #endif
713 // assert(std::abs(theSpectator->getTableMass()-theSpectator->getInvariantMass())<1.e-3);
715  } else {
716  // No geometrical spectators, one dynamical spectator
717  // Just put it back in the outgoing list
718  nucleus->getStore()->addToOutgoing(dynSpectators.front());
719  }
720  nucleus->deleteProjectileRemnant();
721  } else {
722  // Make a cluster out of the geometrical spectators
723  ProjectileRemnant *theProjectileRemnant = nucleus->getProjectileRemnant();
724 
725  // Add the dynamical spectators to the bunch
726  ParticleList rejected = theProjectileRemnant->addMostDynamicalSpectators(dynSpectators);
727  // Put back the rejected spectators into the outgoing list
728  nUnmergedSpectators = rejected.size();
729  nucleus->getStore()->addToOutgoing(rejected);
730 
731  // Deal with the projectile remnant
732  nucleus->finalizeProjectileRemnant(propagationModel->getCurrentTime());
733 
734  }
735 
736  return nUnmergedSpectators;
737  }
738 
739  void INCL::initMaxInteractionDistance(ParticleSpecies const &projectileSpecies, const G4double kineticEnergy) {
740  if(projectileSpecies.theType != Composite) {
741  maxInteractionDistance = 0.;
742  return;
743  }
744 
745  const G4double projectileKineticEnergyPerNucleon = kineticEnergy/projectileSpecies.theA;
747 
748  maxInteractionDistance = r0 + CrossSections::interactionDistanceNN(projectileKineticEnergyPerNucleon);
749  }
750 
751  void INCL::initUniverseRadius(ParticleSpecies const &p, const G4double kineticEnergy, const G4int A, const G4int Z) {
752  G4double rMax = 0.0;
753  if(A==0) {
754  IsotopicDistribution const &anIsotopicDistribution =
756  IsotopeVector theIsotopes = anIsotopicDistribution.getIsotopes();
757  for(IsotopeIter i=theIsotopes.begin(); i!=theIsotopes.end(); ++i) {
758  NuclearDensity *theDensity = NuclearDensityFactory::createDensity(i->theA,Z);
759  if(!theDensity) {
760  FATAL("NULL density in initUniverseRadius. "
761  << "Projectile type=" << p.theType
762  << ", A=" << p.theA
763  << ", Z=" << p.theZ
764  << ", kinE=" << kineticEnergy
765  << ", target A=" << A
766  << ", Z=" << Z
767  << std::endl);
768  }
769  rMax = std::max(theDensity->getMaximumRadius(), rMax);
770  }
771  } else {
772  NuclearDensity *theDensity = NuclearDensityFactory::createDensity(A,Z);
773  if(!theDensity) {
774  FATAL("NULL density in initUniverseRadius. "
775  << "Projectile type=" << p.theType
776  << ", A=" << p.theA
777  << ", Z=" << p.theZ
778  << ", kinE=" << kineticEnergy
779  << ", target A=" << A
780  << ", Z=" << Z
781  << std::endl);
782  }
783  rMax = theDensity->getMaximumRadius();
784  }
785  if(p.theType==Composite) {
786  maxUniverseRadius = rMax;
787  } else if(p.theType==Proton || p.theType==Neutron) {
788  const G4double interactionDistanceNN = CrossSections::interactionDistanceNN(kineticEnergy);
789  if(interactionDistanceNN>CrossSections::interactionDistanceNN1GeV()) {
790  maxUniverseRadius = rMax
792  + interactionDistanceNN;
793  } else
794  maxUniverseRadius = rMax;
795  } else if(p.theType==PiPlus
796  || p.theType==PiZero
797  || p.theType==PiMinus) {
798  const G4double interactionDistancePiN = CrossSections::interactionDistancePiN(kineticEnergy);
799  if(interactionDistancePiN>CrossSections::interactionDistancePiN1GeV()) {
800  maxUniverseRadius = rMax
802  + interactionDistancePiN;
803  } else
804  maxUniverseRadius = rMax;
805  }
806  }
807 
808 }