Geant4-11
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// Alain Boudard, CEA-Saclay, France
28// Joseph Cugnon, University of Liege, Belgium
29// Jean-Christophe David, CEA-Saclay, France
30// Pekka Kaitaniemi, CEA-Saclay, France, and Helsinki Institute of Physics, Finland
31// Sylvie Leray, CEA-Saclay, France
32// Davide Mancusi, CEA-Saclay, France
33//
34#define INCLXX_IN_GEANT4_MODE 1
35
36#include "globals.hh"
37
42#include "G4INCLCascade.hh"
43#include "G4INCLRandom.hh"
46#include "G4INCLParticle.hh"
48#include "G4INCLGlobalInfo.hh"
49
51
53
55
56#include "G4INCLLogger.hh"
57#include "G4INCLGlobals.hh"
59
61
63
64#include "G4INCLClustering.hh"
65
66#include "G4INCLIntersection.hh"
67
69
72
73#include <cstring>
74#include <cstdlib>
75#include <numeric>
76
77namespace G4INCL {
78
79 INCL::INCL(Config const * const config)
80 :propagationModel(0), theA(208), theZ(82), theS(0),
81 targetInitSuccess(false),
83 maxUniverseRadius(0.),
84 maxInteractionDistance(0.),
85 fixedImpactParameter(0.),
86 theConfig(config),
87 nucleus(NULL),
88 forceTransparent(false),
89 minRemnantSize(4)
90 {
91 // Set the logger object.
92#ifdef INCLXX_IN_GEANT4_MODE
94#else // INCLXX_IN_GEANT4_MODE
96#endif // INCLXX_IN_GEANT4_MODE
97
98 // Set the random number generator algorithm. The system can support
99 // multiple different generator algorithms in a completely
100 // transparent way.
102
103 // Select the Pauli and CDPP blocking algorithms
105
106 // Set the cross-section set
108
109 // Set the phase-space generator
111
112 // Select the Coulomb-distortion algorithm:
114
115 // Select the clustering algorithm:
117
118 // Initialize the INCL particle table:
120
121 // Initialize the value of cutNN in BinaryCollisionAvatar
123
124 // Initialize the value of strange cross section bias
126
127 // Propagation model is responsible for finding avatars and
128 // transporting the particles. In principle this step is "hidden"
129 // behind an abstract interface and the rest of the system does not
130 // care how the transportation and avatar finding is done. This
131 // should allow us to "easily" experiment with different avatar
132 // finding schemes and even to support things like curved
133 // trajectories in the future.
137 else
140
143#ifdef INCL_ROOT_USE
144 theGlobalInfo.rootSelection = theConfig->getROOTSelectionString();
145#endif
146
147#ifndef INCLXX_IN_GEANT4_MODE
148 // Fill in the global information
152 const ParticleSpecies theSpecies = theConfig->getProjectileSpecies();
153 theGlobalInfo.Ap = theSpecies.theA;
154 theGlobalInfo.Zp = theSpecies.theZ;
155 theGlobalInfo.Sp = theSpecies.theS;
158#endif
159
161 }
162
165#ifndef INCLXX_IN_GEANT4_MODE
166 NuclearMassTable::deleteTable();
167#endif
174#ifndef INCLXX_IN_GEANT4_MODE
175 Logger::deleteLoggerSlave();
176#endif
180 delete cascadeAction;
181 delete propagationModel;
182 delete theConfig;
183 }
184
185 G4bool INCL::prepareReaction(const ParticleSpecies &projectileSpecies, const G4double kineticEnergy, const G4int A, const G4int Z, const G4int S) {
186 if(A < 0 || A > 300 || Z < 1 || Z > 200) {
187 INCL_ERROR("Unsupported target: A = " << A << " Z = " << Z << " S = " << S << '\n'
188 << "Target configuration rejected." << '\n');
189 return false;
190 }
191 if(projectileSpecies.theType==Composite &&
192 (projectileSpecies.theZ==projectileSpecies.theA || projectileSpecies.theZ==0)) {
193 INCL_ERROR("Unsupported projectile: A = " << projectileSpecies.theA << " Z = " << projectileSpecies.theZ << " S = " << projectileSpecies.theS << '\n'
194 << "Projectile configuration rejected." << '\n');
195 return false;
196 }
197
198 // Reset the forced-transparent flag
199 forceTransparent = false;
200
201 // Initialise the maximum universe radius
202 initUniverseRadius(projectileSpecies, kineticEnergy, A, Z);
203
204 // Initialise the nucleus
205 theZ = Z;
206 theS = S;
209 else
210 theA = A;
212
213 // Set the maximum impact parameter
214 maxImpactParameter = CoulombDistortion::maxImpactParameter(projectileSpecies, kineticEnergy, nucleus);
215 INCL_DEBUG("Maximum impact parameter initialised: " << maxImpactParameter << '\n');
216
217 // For forced CN events
218 initMaxInteractionDistance(projectileSpecies, kineticEnergy);
219
220 // Set the geometric cross section
223
224 // Set the minimum remnant size
225 if(projectileSpecies.theA > 0)
227 else
229
230 return true;
231 }
232
234 delete nucleus;
235
239
241 return true;
242 }
243
245 ParticleSpecies const &projectileSpecies,
246 const G4double kineticEnergy,
247 const G4int targetA,
248 const G4int targetZ,
249 const G4int targetS
250 ) {
251 // ReInitialize the bias vector
253 //Particle::INCLBiasVector.Clear();
255
256 // Set the target and the projectile
257 targetInitSuccess = prepareReaction(projectileSpecies, kineticEnergy, targetA, targetZ, targetS);
258
259 if(!targetInitSuccess) {
260 INCL_WARN("Target initialisation failed for A=" << targetA << ", Z=" << targetZ << ", S=" << targetS << '\n');
262 return theEventInfo;
263 }
264
266
267 const G4bool canRunCascade = preCascade(projectileSpecies, kineticEnergy);
268 if(canRunCascade) {
269 cascade();
270 postCascade();
272 }
274 return theEventInfo;
275 }
276
277 G4bool INCL::preCascade(ParticleSpecies const &projectileSpecies, const G4double kineticEnergy) {
278 // Reset theEventInfo
280
282
283 // Fill in the event information
284 theEventInfo.projectileType = projectileSpecies.theType;
285 theEventInfo.Ap = projectileSpecies.theA;
286 theEventInfo.Zp = projectileSpecies.theZ;
287 theEventInfo.Sp = projectileSpecies.theS;
288 theEventInfo.Ep = kineticEnergy;
292
293 // Do nothing below the Coulomb barrier
294 if(maxImpactParameter<=0.) {
295 // Fill in the event information
297
298 return false;
299 }
300
301 // Randomly draw an impact parameter or use a fixed value, depending on the
302 // Config option
303 G4double impactParameter, phi;
304 if(fixedImpactParameter<0.) {
305 impactParameter = maxImpactParameter * std::sqrt(Random::shoot0());
306 phi = Random::shoot() * Math::twoPi;
307 } else {
308 impactParameter = fixedImpactParameter;
309 phi = 0.;
310 }
311 INCL_DEBUG("Selected impact parameter: " << impactParameter << '\n');
312
313 // Fill in the event information
314 theEventInfo.impactParameter = impactParameter;
315
316 const G4double effectiveImpactParameter = propagationModel->shoot(projectileSpecies, kineticEnergy, impactParameter, phi);
317 if(effectiveImpactParameter < 0.) {
318 // Fill in the event information
320
321 return false;
322 }
323
324 // Fill in the event information
326 theEventInfo.effectiveImpactParameter = effectiveImpactParameter;
327
328 return true;
329 }
330
332 FinalState *finalState = new FinalState;
333
334 unsigned long loopCounter = 0;
335 const unsigned long maxLoopCounter = 10000000;
336 do {
337 // Run book keeping actions that should take place before propagation:
339
340 // Get the avatar with the smallest time and propagate particles
341 // to that point in time.
342 IAvatar *avatar = propagationModel->propagate(finalState);
343
344 finalState->reset();
345
346 // Run book keeping actions that should take place after propagation:
348
349 if(avatar == 0) break; // No more avatars in the avatar list.
350
351 // Run book keeping actions that should take place before avatar:
353
354 // Channel is responsible for calculating the outcome of the
355 // selected avatar. There are different kinds of channels. The
356 // class IChannel is, again, an abstract interface that defines
357 // the externally observable behavior of all interaction
358 // channels.
359 // The handling of the channel is transparent to the API.
360 // Final state tells what changed...
361 avatar->fillFinalState(finalState);
362 // Run book keeping actions that should take place after avatar:
363 cascadeAction->afterAvatarAction(avatar, nucleus, finalState);
364
365 // So now we must give this information to the nucleus
366 nucleus->applyFinalState(finalState);
367 // and now we are ready to process the next avatar!
368
369 delete avatar;
370
371 ++loopCounter;
372 } while(continueCascade() && loopCounter<maxLoopCounter); /* Loop checking, 10.07.2015, D.Mancusi */
373
374 delete finalState;
375 }
376
378 // Fill in the event information
380
381 // The event bias
383
384 // Forced CN?
386 INCL_DEBUG("Trying compound nucleus" << '\n');
389 // Global checks of conservation laws
390#ifndef INCLXX_IN_GEANT4_MODE
391 if(!theEventInfo.transparent) globalConservationChecks(true);
392#endif
393 return;
394 }
395
397
399 ProjectileRemnant * const projectileRemnant = nucleus->getProjectileRemnant();
400 if(projectileRemnant) {
401 // Clear the incoming list (particles will be deleted by the ProjectileRemnant)
403 } else {
404 // Delete particles in the incoming list
406 }
407 } else {
408
409 // Check if the nucleus contains strange particles
414
415 // Capture antiKaons and Sigmas and produce Lambda instead
417
418 // Emit strange particles still inside the nucleus
421
422#ifdef INCLXX_IN_GEANT4_MODE
424#endif // INCLXX_IN_GEANT4_MODE
425
426 // Check if the nucleus contains deltas
428
429 // Take care of any remaining deltas
432
433 // Take care of any remaining etas, omegas, neutral Sigmas and/or neutral kaons
436 nucleus->decayOutgoingSigmaZero(timeThreshold);
438
439 // Apply Coulomb distortion, if appropriate
440 // Note that this will apply Coulomb distortion also on pions emitted by
441 // unphysical remnants (see decayInsideDeltas). This is at variance with
442 // what INCL4.6 does, but these events are (should be!) so rare that
443 // whatever we do doesn't (shouldn't!) make any noticeable difference.
445
446 // If the normal cascade predicted complete fusion, use the tabulated
447 // masses to compute the excitation energy, the recoil, etc.
448 if(nucleus->getStore()->getOutgoingParticles().size()==0
450 || nucleus->getProjectileRemnant()->getParticles().size()==0)) {
451
452 INCL_DEBUG("Cascade resulted in complete fusion, using realistic fusion kinematics" << '\n');
453
455
456 if(nucleus->getExcitationEnergy()<0.) {
457 // Complete fusion is energetically impossible, return a transparent
458 INCL_WARN("Complete-fusion kinematics yields negative excitation energy, returning a transparent!" << '\n');
460 return;
461 }
462
463 } else { // Normal cascade here
464
465 // Set the excitation energy
467
468 // Make a projectile pre-fragment out of the geometrical and dynamical
469 // spectators
471
472 // Compute recoil momentum, energy and spin of the nucleus
473 if(nucleus->getA()==1 && minRemnantSize>1) {
474 INCL_ERROR("Computing one-nucleon recoil kinematics. We should never be here nowadays, cascade should stop earlier than this." << '\n');
475 }
477
478#ifndef INCLXX_IN_GEANT4_MODE
479 // Global checks of conservation laws
480 globalConservationChecks(false);
481#endif
482
483 // Make room for the remnant recoil by rescaling the energies of the
484 // outgoing particles.
486
487 }
488
489 // Cluster decay
491
492#ifndef INCLXX_IN_GEANT4_MODE
493 // Global checks of conservation laws
494 globalConservationChecks(true);
495#endif
496
497 // Fill the EventInfo structure
499
500 }
501 }
502
504 // If this is not a nucleus-nucleus collision, don't attempt to make a
505 // compound nucleus.
506 //
507 // Yes, even nucleon-nucleus collisions can lead to particles entering
508 // below the Fermi level. Take e.g. 1-MeV p + He4.
510 forceTransparent = true;
511 return;
512 }
513
514 // Reset the internal Nucleus variables
520
521 // CN kinematical variables
522 // Note: the CN orbital angular momentum is neglected in what follows. We
523 // should actually take it into account!
524 ThreeVector theCNMomentum = nucleus->getIncomingMomentum();
527 G4int theCNA=theEventInfo.At, theCNZ=theEventInfo.Zt, theCNS=theEventInfo.St;
528 Cluster * const theProjectileRemnant = nucleus->getProjectileRemnant();
529 G4double theCNEnergy = theTargetMass + theProjectileRemnant->getEnergy();
530
531 // Loop over the potential participants
532 ParticleList const &initialProjectileComponents = theProjectileRemnant->getParticles();
533 std::vector<Particle *> shuffledComponents(initialProjectileComponents.begin(), initialProjectileComponents.end());
534 // Shuffle the list of potential participants
535 std::shuffle(shuffledComponents.begin(), shuffledComponents.end(), Random::getAdapter());
536
537 G4bool success = true;
538 G4bool atLeastOneNucleonEntering = false;
539 for(std::vector<Particle*>::const_iterator p=shuffledComponents.begin(), e=shuffledComponents.end(); p!=e; ++p) {
540 // Skip particles that miss the interaction distance
542 (*p)->getPosition(),
543 (*p)->getPropagationVelocity(),
545 if(!intersectionInteractionDistance.exists)
546 continue;
547
548 // Build an entry avatar for this nucleon
549 atLeastOneNucleonEntering = true;
550 ParticleEntryAvatar *theAvatar = new ParticleEntryAvatar(0.0, nucleus, *p);
552 FinalState *fs = theAvatar->getFinalState();
554 FinalStateValidity validity = fs->getValidity();
555 delete fs;
556 switch(validity) {
557 case ValidFS:
560 // Add the particle to the CN
561 theCNA++;
562 theCNZ += (*p)->getZ();
563 theCNS += (*p)->getS();
564 break;
565 case PauliBlockedFS:
567 default:
568 success = false;
569 break;
570 }
571 }
572
573 if(!success || !atLeastOneNucleonEntering) {
574 INCL_DEBUG("No nucleon entering in forced CN, forcing a transparent" << '\n');
575 forceTransparent = true;
576 return;
577 }
578
579// assert(theCNA==nucleus->getA());
580// assert(theCNA<=theEventInfo.At+theEventInfo.Ap);
581// assert(theCNZ<=theEventInfo.Zt+theEventInfo.Zp);
582// assert(theCNS>=theEventInfo.St+theEventInfo.Sp);
583
584 // Update the kinematics of the CN
585 theCNEnergy -= theProjectileRemnant->getEnergy();
586 theCNMomentum -= theProjectileRemnant->getMomentum();
587
588 // Deal with the projectile remnant
590
591 // Subtract the angular momentum of the projectile remnant
592// assert(nucleus->getStore()->getOutgoingParticles().empty());
593 theCNSpin -= theProjectileRemnant->getAngularMomentum();
594
595 // Compute the excitation energy of the CN
596 const G4double theCNMass = ParticleTable::getTableMass(theCNA,theCNZ,theCNS);
597 const G4double theCNInvariantMassSquared = theCNEnergy*theCNEnergy-theCNMomentum.mag2();
598 if(theCNInvariantMassSquared<0.) {
599 // Negative invariant mass squared, return a transparent
600 forceTransparent = true;
601 return;
602 }
603 const G4double theCNExcitationEnergy = std::sqrt(theCNInvariantMassSquared) - theCNMass;
604 if(theCNExcitationEnergy<0.) {
605 // Negative excitation energy, return a transparent
606 INCL_DEBUG("CN excitation energy is negative, forcing a transparent" << '\n'
607 << " theCNA = " << theCNA << '\n'
608 << " theCNZ = " << theCNZ << '\n'
609 << " theCNS = " << theCNS << '\n'
610 << " theCNEnergy = " << theCNEnergy << '\n'
611 << " theCNMomentum = (" << theCNMomentum.getX() << ", "<< theCNMomentum.getY() << ", " << theCNMomentum.getZ() << ")" << '\n'
612 << " theCNExcitationEnergy = " << theCNExcitationEnergy << '\n'
613 << " theCNSpin = (" << theCNSpin.getX() << ", "<< theCNSpin.getY() << ", " << theCNSpin.getZ() << ")" << '\n'
614 );
615 forceTransparent = true;
616 return;
617 } else {
618 // Positive excitation energy, can make a CN
619 INCL_DEBUG("CN excitation energy is positive, forcing a CN" << '\n'
620 << " theCNA = " << theCNA << '\n'
621 << " theCNZ = " << theCNZ << '\n'
622 << " theCNS = " << theCNS << '\n'
623 << " theCNEnergy = " << theCNEnergy << '\n'
624 << " theCNMomentum = (" << theCNMomentum.getX() << ", "<< theCNMomentum.getY() << ", " << theCNMomentum.getZ() << ")" << '\n'
625 << " theCNExcitationEnergy = " << theCNExcitationEnergy << '\n'
626 << " theCNSpin = (" << theCNSpin.getX() << ", "<< theCNSpin.getY() << ", " << theCNSpin.getZ() << ")" << '\n'
627 );
628 nucleus->setA(theCNA);
629 nucleus->setZ(theCNZ);
630 nucleus->setS(theCNS);
631 nucleus->setMomentum(theCNMomentum);
632 nucleus->setEnergy(theCNEnergy);
633 nucleus->setExcitationEnergy(theCNExcitationEnergy);
634 nucleus->setMass(theCNMass+theCNExcitationEnergy);
635 nucleus->setSpin(theCNSpin); // neglects any orbital angular momentum of the CN
636
637 // Take care of any remaining deltas
639
640 // Take care of any remaining etas and/or omegas
643
644 // Take care of any remaining Kaons
646
647 // Cluster decay
649
650 // Fill the EventInfo structure
652 }
653 }
654
656 RecoilCMFunctor theRecoilFunctor(nucleus, theEventInfo);
657
658 // Apply the root-finding algorithm
659 const RootFinder::Solution theSolution = RootFinder::solve(&theRecoilFunctor, 1.0);
660 if(theSolution.success) {
661 theRecoilFunctor(theSolution.x); // Apply the solution
662 } else {
663 INCL_WARN("Couldn't accommodate remnant recoil while satisfying energy conservation, root-finding algorithm failed." << '\n');
664 }
665
666 }
667
668#ifndef INCLXX_IN_GEANT4_MODE
669 void INCL::globalConservationChecks(G4bool afterRecoil) {
671
672 // Global conservation checks
673 const G4double pLongBalance = theBalance.momentum.getZ();
674 const G4double pTransBalance = theBalance.momentum.perp();
675 if(theBalance.Z != 0) {
676 INCL_ERROR("Violation of charge conservation! ZBalance = " << theBalance.Z << " eventNumber=" << theEventInfo.eventNumber << '\n');
677 }
678 if(theBalance.A != 0) {
679 INCL_ERROR("Violation of baryon-number conservation! ABalance = " << theBalance.A << " Emit Lambda=" << theEventInfo.emitLambda << " eventNumber=" << theEventInfo.eventNumber << '\n');
680 }
681 if(theBalance.S != 0) {
682 INCL_ERROR("Violation of strange-number conservation! SBalance = " << theBalance.S << " eventNumber=" << theEventInfo.eventNumber << '\n');
683 }
684 G4double EThreshold, pLongThreshold, pTransThreshold;
685 if(afterRecoil) {
686 // Less stringent checks after accommodating recoil
687 EThreshold = 10.; // MeV
688 pLongThreshold = 1.; // MeV/c
689 pTransThreshold = 1.; // MeV/c
690 } else {
691 // More stringent checks before accommodating recoil
692 EThreshold = 0.1; // MeV
693 pLongThreshold = 0.1; // MeV/c
694 pTransThreshold = 0.1; // MeV/c
695 }
696 if(std::abs(theBalance.energy)>EThreshold) {
697 INCL_WARN("Violation of energy conservation > " << EThreshold << " MeV. EBalance = " << theBalance.energy << " Emit Lambda=" << theEventInfo.emitLambda << " afterRecoil = " << afterRecoil << " eventNumber=" << theEventInfo.eventNumber << '\n');
698 }
699 if(std::abs(pLongBalance)>pLongThreshold) {
700 INCL_WARN("Violation of longitudinal momentum conservation > " << pLongThreshold << " MeV/c. pLongBalance = " << pLongBalance << " afterRecoil = " << afterRecoil << " eventNumber=" << theEventInfo.eventNumber << '\n');
701 }
702 if(std::abs(pTransBalance)>pTransThreshold) {
703 INCL_WARN("Violation of transverse momentum conservation > " << pTransThreshold << " MeV/c. pTransBalance = " << pTransBalance << " afterRecoil = " << afterRecoil << " eventNumber=" << theEventInfo.eventNumber << '\n');
704 }
705
706 // Feed the EventInfo variables
707 theEventInfo.EBalance = theBalance.energy;
708 theEventInfo.pLongBalance = pLongBalance;
709 theEventInfo.pTransBalance = pTransBalance;
710 }
711#endif
712
714 // Stop if we have passed the stopping time
716 INCL_DEBUG("Cascade time (" << propagationModel->getCurrentTime()
717 << ") exceeded stopping time (" << propagationModel->getStoppingTime()
718 << "), stopping cascade" << '\n');
719 return false;
720 }
721 // Stop if there are no participants and no pions inside the nucleus
722 if(nucleus->getStore()->getBook().getCascading()==0 &&
723 nucleus->getStore()->getIncomingParticles().empty()) {
724 INCL_DEBUG("No participants in the nucleus and no incoming particles left, stopping cascade" << '\n');
725 return false;
726 }
727 // Stop if the remnant is smaller than minRemnantSize
728 if(nucleus->getA() <= minRemnantSize) {
729 INCL_DEBUG("Remnant size (" << nucleus->getA()
730 << ") smaller than or equal to minimum (" << minRemnantSize
731 << "), stopping cascade" << '\n');
732 return false;
733 }
734 // Stop if we have to try and make a compound nucleus or if we have to
735 // force a transparent
737 INCL_DEBUG("Trying to make a compound nucleus, stopping cascade" << '\n');
738 return false;
739 }
740
741 return true;
742 }
743
745 const G4double normalisationFactor = theGlobalInfo.geometricCrossSection /
747 theGlobalInfo.nucleonAbsorptionCrossSection = normalisationFactor *
749 theGlobalInfo.pionAbsorptionCrossSection = normalisationFactor *
751 theGlobalInfo.reactionCrossSection = normalisationFactor *
753 theGlobalInfo.errorReactionCrossSection = normalisationFactor *
755 theGlobalInfo.forcedCNCrossSection = normalisationFactor *
757 theGlobalInfo.errorForcedCNCrossSection = normalisationFactor *
759 theGlobalInfo.completeFusionCrossSection = normalisationFactor *
761 theGlobalInfo.errorCompleteFusionCrossSection = normalisationFactor *
765
766 theGlobalInfo.initialRandomSeeds.assign(initialSeeds.begin(), initialSeeds.end());
767
769 theGlobalInfo.finalRandomSeeds.assign(theSeeds.begin(), theSeeds.end());
770 }
771
773 // Do nothing if this is not a nucleus-nucleus reaction
775 return 0;
776
777 // Get the spectators (geometrical+dynamical) from the Store
780
781 G4int nUnmergedSpectators = 0;
782
783 // If there are no spectators, do nothing
784 if(dynSpectators.empty() && geomSpectators.empty()) {
785 return 0;
786 } else if(dynSpectators.size()==1 && geomSpectators.empty()) {
787 // No geometrical spectators, one dynamical spectator
788 // Just put it back in the outgoing list
789 nucleus->getStore()->addToOutgoing(dynSpectators.front());
790 } else {
791 // Make a cluster out of the geometrical spectators
792 ProjectileRemnant *theProjectileRemnant = nucleus->getProjectileRemnant();
793
794 // Add the dynamical spectators to the bunch
795 ParticleList rejected = theProjectileRemnant->addAllDynamicalSpectators(dynSpectators);
796 // Put back the rejected spectators into the outgoing list
797 nUnmergedSpectators = rejected.size();
798 nucleus->getStore()->addToOutgoing(rejected);
799
800 // Deal with the projectile remnant
802
803 }
804
805 return nUnmergedSpectators;
806 }
807
808 void INCL::initMaxInteractionDistance(ParticleSpecies const &projectileSpecies, const G4double kineticEnergy) {
809 if(projectileSpecies.theType != Composite) {
811 return;
812 }
813
816
817 const G4double theNNDistance = CrossSections::interactionDistanceNN(projectileSpecies, kineticEnergy);
818 maxInteractionDistance = r0 + theNNDistance;
819 INCL_DEBUG("Initialised interaction distance: r0 = " << r0 << '\n'
820 << " theNNDistance = " << theNNDistance << '\n'
821 << " maxInteractionDistance = " << maxInteractionDistance << '\n');
822 }
823
824 void INCL::initUniverseRadius(ParticleSpecies const &p, const G4double kineticEnergy, const G4int A, const G4int Z) {
825 G4double rMax = 0.0;
826 if(A==0) {
827 IsotopicDistribution const &anIsotopicDistribution =
829 IsotopeVector theIsotopes = anIsotopicDistribution.getIsotopes();
830 for(IsotopeIter i=theIsotopes.begin(), e=theIsotopes.end(); i!=e; ++i) {
831 const G4double pMaximumRadius = ParticleTable::getMaximumNuclearRadius(Proton, i->theA, Z);
832 const G4double nMaximumRadius = ParticleTable::getMaximumNuclearRadius(Neutron, i->theA, Z);
833 const G4double maximumRadius = std::max(pMaximumRadius, nMaximumRadius);
834 rMax = std::max(maximumRadius, rMax);
835 }
836 } else {
837 const G4double pMaximumRadius = ParticleTable::getMaximumNuclearRadius(Proton, A, Z);
839 const G4double maximumRadius = std::max(pMaximumRadius, nMaximumRadius);
840 rMax = std::max(maximumRadius, rMax);
841 }
842 if(p.theType==Composite || p.theType==Proton || p.theType==Neutron) {
845 } else if(p.theType==PiPlus
846 || p.theType==PiZero
847 || p.theType==PiMinus) {
850 } else if(p.theType==KPlus
851 || p.theType==KZero) {
854 } else if(p.theType==KZeroBar
855 || p.theType==KMinus) {
858 } else if(p.theType==Lambda
859 ||p.theType==SigmaPlus
860 || p.theType==SigmaZero
861 || p.theType==SigmaMinus) {
864 }
865 INCL_DEBUG("Initialised universe radius: " << maxUniverseRadius << '\n');
866 }
867
869 // Increment the global counter for the number of shots
871
873 // Increment the global counter for the number of transparents
875 // Increment the global counter for the number of forced transparents
878 return;
879 }
880
881 // Check if we have an absorption:
884
885 // Count complete-fusion events
887
890
891 // Counters for the number of violations of energy conservation in
892 // collisions
894 }
895
896}
G4double S(G4double temp)
Alternative CascadeAction for dumping avatars.
Class containing default actions to be performed at intermediate cascade steps.
Static class for cluster formation.
Static class for selecting Coulomb distortion.
Simple container for output of calculation-wide results.
Abstract interface to the nuclear potential.
Simple class for computing intersections between a straight line and a sphere.
#define INCL_ERROR(x)
#define INCL_WARN(x)
#define INCL_DEBUG(x)
Functions that encapsulate a mass table.
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]
static void setBias(const G4double b)
Set the global bias factor.
static void setCutNN(const G4double c)
G4int getCascading() const
Definition: G4INCLBook.hh:105
void reset()
Definition: G4INCLBook.hh:52
void beforeCascadeAction(IPropagationModel *)
void afterAvatarAction(IAvatar *a, Nucleus *n, FinalState *fs)
void afterPropagationAction(IPropagationModel *pm, IAvatar *avatar)
void afterCascadeAction(Nucleus *)
void beforeRunAction(Config const *config)
void beforeAvatarAction(IAvatar *a, Nucleus *n)
void beforePropagationAction(IPropagationModel *pm)
ParticleList const & getParticles() const
void setExcitationEnergy(const G4double e)
Set the excitation energy of the cluster.
void setS(const G4int S)
Set the strangess number of the cluster.
void setSpin(const ThreeVector &j)
Set the spin of the nucleus.
G4INCL::ThreeVector getAngularMomentum() const
Get the total angular momentum (orbital + spin)
void setA(const G4int A)
Set the mass number of the cluster.
void setZ(const G4int Z)
Set the charge number of the cluster.
static std::string const getVersionString()
Get the INCL version string.
G4double getImpactParameter() const
CascadeActionType getCascadeActionType() const
Get the cascade-action type.
G4int getTargetA() const
Get the target mass number.
Definition: G4INCLConfig.hh:94
G4double getDecayTimeThreshold() const
Get the decay time threshold time.
G4double getProjectileKineticEnergy() const
Get the projectile kinetic energy.
G4double getHadronizationTime() const
Get the hadronization time.
ParticleSpecies getProjectileSpecies() const
Get the projectile species.
LocalEnergyType getLocalEnergyPiType() const
Get the type of local energy for pi-N and decay avatars.
G4bool isNaturalTarget() const
Natural targets.
Definition: G4INCLConfig.hh:87
G4int getTargetS() const
Get the target strangess number.
G4double getCutNN() const
LocalEnergyType getLocalEnergyBBType() const
Get the type of local energy for N-N avatars.
G4double getBias() const
Get the bias.
G4int getTargetZ() const
Get the target charge number.
Definition: G4INCLConfig.hh:97
std::string getDeExcitationString() const
Get the de-excitation string.
FinalStateValidity getValidity() const
void fillFinalState(FinalState *fs)
FinalState * getFinalState()
Class to adjust remnant recoil in the reaction CM system.
void initUniverseRadius(ParticleSpecies const &p, const G4double kineticEnergy, const G4int A, const G4int Z)
Initialize the universe radius.
void finalizeGlobalInfo(Random::SeedVector const &initialSeeds)
void initMaxInteractionDistance(ParticleSpecies const &p, const G4double kineticEnergy)
Initialise the maximum interaction distance.
G4bool prepareReaction(const ParticleSpecies &projectileSpecies, const G4double kineticEnergy, const G4int A, const G4int Z, const G4int S)
G4double maxImpactParameter
G4bool initializeTarget(const G4int A, const G4int Z, const G4int S)
INCL(Config const *const config)
EventInfo theEventInfo
const EventInfo & processEvent()
G4bool continueCascade()
Stopping criterion for the cascade.
Config const *const theConfig
GlobalInfo theGlobalInfo
CascadeAction * cascadeAction
G4int minRemnantSize
Remnant size below which cascade stops.
G4double maxUniverseRadius
G4bool targetInitSuccess
void makeCompoundNucleus()
Make a compound nucleus.
Nucleus * nucleus
G4double maxInteractionDistance
void cascade()
The actual cascade loop.
G4double fixedImpactParameter
IPropagationModel * propagationModel
G4int makeProjectileRemnant()
Make a projectile pre-fragment out of geometrical spectators.
void updateGlobalInfo()
Update global counters and other members of theGlobalInfo object.
void postCascade()
Finalise the cascade and clean up.
void rescaleOutgoingForRecoil()
Rescale the energies of the outgoing particles.
G4bool preCascade(ParticleSpecies const &projectileSpecies, const G4double kineticEnergy)
Initialise the cascade.
G4bool forceTransparent
virtual G4double shoot(ParticleSpecies const &projectileSpecies, const G4double kineticEnergy, const G4double impactParameter, const G4double phi)=0
virtual G4INCL::IAvatar * propagate(FinalState const *const fs)=0
virtual G4double getCurrentTime()=0
virtual G4double getStoppingTime()=0
virtual void setNucleus(G4INCL::Nucleus *nucleus)=0
static void deleteBackupParticles()
Release the memory allocated for the backup particles.
Class that stores isotopic abundances for a given element.
IsotopeVector const & getIsotopes() const
Get the isotope vector.
G4bool decayOutgoingNeutralKaon()
Force the transformation of outgoing Neutral Kaon into propation eigenstate.
G4bool containsSigma()
Returns true if the nucleus contains any Sigma.
const ThreeVector & getIncomingAngularMomentum() const
Get the incoming angular-momentum vector.
Store * getStore() const
ProjectileRemnant * getProjectileRemnant() const
Get the projectile remnant.
void fillEventInfo(EventInfo *eventInfo)
G4bool decayMe()
Force the phase-space decay of the Nucleus.
G4bool decayInsideDeltas()
Force the decay of deltas inside the nucleus.
void computeRecoilKinematics()
Compute the recoil momentum and spin of the nucleus.
void finalizeProjectileRemnant(const G4double emissionTime)
Finalise the projectile remnant.
G4bool containsDeltas()
Returns true if the nucleus contains any deltas.
G4bool getTryCompoundNucleus()
G4bool emitInsideKaon()
Force emission of all Kaon inside the nucleus.
G4bool isEventTransparent() const
Is the event transparent?
void applyFinalState(FinalState *)
G4bool isNucleusNucleusCollision() const
Is it a nucleus-nucleus collision?
G4bool containsKaon()
Returns true if the nucleus contains any Kaons.
void emitInsideStrangeParticles()
Force emission of all strange particles inside the nucleus.
G4bool decayOutgoingPionResonances(G4double timeThreshold)
Force the decay of outgoing PionResonances (eta/omega).
void useFusionKinematics()
Adjust the kinematics for complete-fusion events.
G4bool containsLambda()
Returns true if the nucleus contains any Lambda.
const ThreeVector & getIncomingMomentum() const
Get the incoming momentum vector.
G4bool decayOutgoingClusters()
Force the decay of unstable outgoing clusters.
G4double getExcitationEnergy() const
Get the excitation energy of the nucleus.
G4int emitInsideLambda()
Force emission of all Lambda (desexitation code with strangeness not implanted yet)
G4bool decayInsideStrangeParticles()
Force the transformation of strange particles into a Lambda;.
ConservationBalance getConservationBalance(EventInfo const &theEventInfo, const G4bool afterRecoil) const
Compute charge, mass, energy and momentum balance.
void initializeParticles()
G4bool hasRemnant() const
Does the nucleus give a cascade remnant?
G4double computeExcitationEnergy() const
Compute the current excitation energy.
G4bool containsAntiKaon()
Returns true if the nucleus contains any anti Kaons.
G4bool decayOutgoingDeltas()
Force the decay of outgoing deltas.
G4bool decayOutgoingSigmaZero(G4double timeThreshold)
Force the decay of outgoing Neutral Sigma.
G4int getS() const
Returns the strangeness number.
G4double getEnergy() const
static std::vector< G4double > INCLBiasVector
Time ordered vector of all bias applied.
void setMass(G4double mass)
G4int getZ() const
Returns the charge number.
static G4double getTotalBias()
General bias vector function.
const G4INCL::ThreeVector & getMomentum() const
virtual void setMomentum(const G4INCL::ThreeVector &momentum)
void setEnergy(G4double energy)
static G4ThreadLocal G4int nextBiasedCollisionID
G4int getA() const
Returns the baryon number.
ParticleList addAllDynamicalSpectators(ParticleList const &pL)
Add back all dynamical spectators to the projectile remnant.
void reset()
Reset the projectile remnant to the state at the beginning of the cascade.
ParticleList const & getIncomingParticles() const
Definition: G4INCLStore.hh:217
void deleteIncoming()
Clear the incoming list and delete the particles.
Definition: G4INCLStore.hh:150
ParticleList const & getOutgoingParticles() const
Definition: G4INCLStore.hh:223
void addToOutgoing(Particle *p)
add the particle to the outgoing particle list.
Definition: G4INCLStore.hh:190
void clearOutgoing()
Definition: G4INCLStore.cc:223
Book & getBook()
Definition: G4INCLStore.hh:259
void clearIncoming()
Clear the incoming list.
Definition: G4INCLStore.hh:145
void addParticleEntryAvatar(IAvatar *a)
Add one ParticleEntry avatar.
Definition: G4INCLStore.cc:66
ParticleList extractDynamicalSpectators()
Returns a list of dynamical spectators.
Definition: G4INCLStore.hh:232
G4double getY() const
G4double getZ() const
G4double perp() const
G4double mag2() const
G4double getX() const
void initialize(Config const *const theConfig)
Initialize the clustering model based on the Config object.
void deleteClusteringModel()
Delete the clustering model.
void initialize(Config const *const theConfig)
Initialize the Coulomb-distortion algorithm.
void deleteCoulomb()
Delete the Coulomb-distortion object.
void distortOut(ParticleList const &pL, Nucleus const *const n)
Modify the momentum of an outgoing particle.
G4double maxImpactParameter(ParticleSpecies const &p, const G4double kinE, Nucleus const *const n)
Return the maximum impact parameter for Coulomb-distorted trajectories.
G4double interactionDistanceKbarN(const G4double projectileKineticEnergy)
Compute the "interaction distance".
G4double interactionDistancePiN(const G4double projectileKineticEnergy)
Compute the "interaction distance".
G4double interactionDistanceKN(const G4double projectileKineticEnergy)
Compute the "interaction distance".
G4double interactionDistanceYN(const G4double projectileKineticEnergy)
Compute the "interaction distance".
void initialize(Config const *const theConfig)
G4double interactionDistanceNN(const ParticleSpecies &aSpecies, const G4double kineticEnergy)
Compute the "interaction distance".
Intersection getEarlierTrajectoryIntersection(const ThreeVector &x0, const ThreeVector &p, const G4double r)
Compute the first intersection of a straight particle trajectory with a sphere.
void initVerbosityLevelFromEnvvar()
T max(const T t1, const T t2)
brief Return the largest of the two arguments
const G4double twoPi
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
const G4double tenPi
void clearCache()
Clear the INuclearPotential cache.
G4ThreadLocal NuclearMassFn getTableMass
Static pointer to the mass function for nuclei.
void initialize(Config const *const theConfig=0)
Initialize the particle table.
G4int drawRandomNaturalIsotope(const G4int Z)
G4double getMaximumNuclearRadius(const ParticleType t, const G4int A, const G4int Z)
IsotopicDistribution const & getNaturalIsotopicDistribution(const G4int Z)
G4double getNuclearRadius(const ParticleType t, const G4int A, const G4int Z)
void initialize(Config const *const aConfig)
Initialise blockers according to a Config object.
void deleteBlockers()
Delete blockers.
void initialize(Config const *const theConfig)
SeedVector getSeeds()
Definition: G4INCLRandom.cc:89
G4double shoot0()
Adapter const & getAdapter()
G4double shoot()
Definition: G4INCLRandom.cc:93
void deleteGenerator()
void initialize(Config const *const)
Initialize generator according to a Config object.
Solution solve(RootFunctor const *const f, const G4double x0)
Numerically solve a one-dimensional equation.
IsotopeVector::iterator IsotopeIter
@ ParticleBelowZeroFS
@ NoEnergyConservationFS
@ ParticleBelowFermiFS
G4double Double_t
std::vector< Isotope > IsotopeVector
Bool_t pionAbsorption
True if the event is a pion absorption.
void reset()
Reset the EventInfo members.
Bool_t lambdasInside
Event involved lambdas in the nucleus at the end of the cascade.
Bool_t forcedDeltasOutside
Event involved forced delta decays outside the nucleus.
Short_t At
Mass number of the target nucleus.
Short_t Zp
Charge number of the projectile nucleus.
Int_t projectileType
Projectile particle type.
Float_t pTransBalance
Transverse momentum-conservation balance [MeV/c].
Float_t Ep
Projectile kinetic energy given as input.
Short_t nCascadeParticles
Number of cascade particles.
Float_t eventBias
Event bias.
Bool_t transparent
True if the event is transparent.
Bool_t sigmasInside
Event involved sigmas in the nucleus at the end of the cascade.
Bool_t deltasInside
Event involved deltas in the nucleus at the end of the cascade.
Float_t EBalance
Energy-conservation balance [MeV].
Bool_t forcedPionResonancesOutside
Event involved forced eta/omega decays outside the nucleus.
Bool_t nucleonAbsorption
True if the event is a nucleon absorption.
Int_t emitLambda
Number of forced Lambda emit out of the nucleus.
Short_t St
Strangeness number of the target nucleus.
Float_t effectiveImpactParameter
Effective (Coulomb-distorted) impact parameter [fm].
Bool_t kaonsInside
Event involved kaons in the nucleus at the end of the cascade.
Float_t stoppingTime
Cascade stopping time [fm/c].
Int_t nEnergyViolationInteraction
Number of attempted collisions/decays for which the energy-conservation algorithm failed to find a so...
Short_t Ap
Mass number of the projectile nucleus.
Short_t Sp
Strangeness number of the projectile nucleus.
Bool_t emitKaon
Event involved forced Kaon emission.
Bool_t antikaonsInside
Event involved antikaons in the nucleus at the end of the cascade.
Bool_t clusterDecay
Event involved cluster decay.
Bool_t forcedDeltasInside
Event involved forced delta decays inside the nucleus.
Short_t Zt
Charge number of the target nucleus.
static G4ThreadLocal Int_t eventNumber
Number of the event.
Float_t impactParameter
Impact parameter [fm].
Bool_t absorbedStrangeParticle
Event involved forced strange absorption inside the nucleus.
Int_t nUnmergedSpectators
Number of dynamical spectators that were merged back into the projectile remnant.
Float_t pLongBalance
Longitudinal momentum-conservation balance [MeV/c].
Int_t nNucleonAbsorptions
Number of nucleon absorptions (no outcoming particles)
std::string deexcitationModel
Name of the de-excitation model.
Float_t forcedCNCrossSection
Calculated forced-compound-nucleus cross section.
Float_t pionAbsorptionCrossSection
Pion absorption cross section.
Float_t errorCompleteFusionCrossSection
Error on the calculated complete-fusion cross section (nParticles==0)
Short_t St
Target strangeness number given as input.
Float_t Ep
Projectile kinetic energy given as input.
Int_t nEnergyViolationInteraction
Number of attempted collisions/decays for which the energy-conservation algorithm failed to find a so...
Float_t biasFactor
Bias factor.
std::vector< Int_t > initialRandomSeeds
Initial seeds for the pseudo-random-number generator.
Short_t At
Target mass number given as input.
Int_t nShots
Number of shots.
Float_t completeFusionCrossSection
Calculated complete-fusion cross section (nParticles==0)
Float_t nucleonAbsorptionCrossSection
Nucleon absorption cross section.
Float_t errorReactionCrossSection
Error on the calculated reaction cross section.
Short_t Zt
Target charge number given as input.
std::vector< Int_t > finalRandomSeeds
Final seeds for the pseudo-random-number generator.
Int_t nPionAbsorptions
Number of nucleon absorptions (no outcoming pions)
Short_t Ap
Projectile mass number given as input.
Float_t energyViolationInteractionCrossSection
Cross section for attempted collisions/decays for which the energy-conservation algorithm failed to f...
Int_t nCompleteFusion
Number of complete-fusion events (nParticles==0)
Int_t nTransparents
Number of transparent shots.
Short_t Sp
Projectile strangeness number given as input.
Int_t nForcedCompoundNucleus
Number of forced compound-nucleus events.
Short_t Zp
Projectile charge number given as input.
Int_t nForcedTransparents
Number of forced transparents.
std::string cascadeModel
Name of the cascade model.
Float_t reactionCrossSection
Calculated reaction cross section.
Float_t errorForcedCNCrossSection
Error on the calculated forced-compound-nucleus cross section.
Float_t geometricCrossSection
Geometric cross section.
Intersection-point structure.
Struct for conservation laws.