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