Geant4-11
Public Member Functions | Private Attributes
G4INCL::StandardPropagationModel Class Reference

#include <G4INCLStandardPropagationModel.hh>

Inheritance diagram for G4INCL::StandardPropagationModel:
G4INCL::IPropagationModel

Public Member Functions

void generateAllAvatars ()
 (Re)Generate all possible avatars. More...
 
IAvatargenerateBinaryCollisionAvatar (Particle *const p1, Particle *const p2)
 Generate a two-particle avatar. More...
 
void generateCollisions (const ParticleList &particles)
 Generate and register collisions among particles in a list, except between those in another list. More...
 
void generateCollisions (const ParticleList &particles, const ParticleList &except)
 Generate and register collisions among particles in a list, except between those in another list. More...
 
void generateDecays (const ParticleList &particles)
 Generate decays for particles that can decay. More...
 
void generateUpdatedCollisions (const ParticleList &updatedParticles, const ParticleList &particles)
 Generate and register collisions between a list of updated particles and all the other particles. More...
 
G4double getCurrentTime ()
 
G4INCL::NucleusgetNucleus ()
 
G4double getReflectionTime (G4INCL::Particle const *const aParticle)
 Get the reflection time. More...
 
G4double getStoppingTime ()
 
G4double getTime (G4INCL::Particle const *const particleA, G4INCL::Particle const *const particleB, G4double *minDistOfApproach) const
 
G4INCL::IAvatarpropagate (FinalState const *const fs)
 
void registerAvatar (G4INCL::IAvatar *anAvatar)
 
void setNucleus (G4INCL::Nucleus *nucleus)
 
void setStoppingTime (G4double)
 
G4double shoot (ParticleSpecies const &projectileSpecies, const G4double kineticEnergy, const G4double impactParameter, const G4double phi)
 
G4double shootComposite (ParticleSpecies const &s, const G4double kineticEnergy, const G4double impactParameter, const G4double phi)
 
G4double shootParticle (ParticleType const t, const G4double kineticEnergy, const G4double impactParameter, const G4double phi)
 
 StandardPropagationModel (LocalEnergyType localEnergyType, LocalEnergyType localEnergyDeltaType, const G4double hTime=0.0)
 
void updateAvatars (const ParticleList &particles)
 
virtual ~StandardPropagationModel ()
 

Private Attributes

Particle backupParticle1
 
Particle backupParticle2
 
G4double currentTime
 
G4bool firstAvatar
 
G4double hadronizationTime
 
G4double maximumTime
 
LocalEnergyType theLocalEnergyDeltaType
 
LocalEnergyType theLocalEnergyType
 
G4INCL::NucleustheNucleus
 

Detailed Description

Standard INCL4 particle propagation and avatar prediction

This class implements the standard INCL4 avatar prediction and particle propagation logic. The main idea is to predict all collisions between particles and their reflections from the potential wall. After this we select the avatar with the smallest time, propagate all particles to their positions at that time and return the avatar to the INCL kernel

See also
G4INCL::Kernel.

The particle trajectories in this propagation model are straight lines and all particles are assumed to move with constant velocity.

Definition at line 69 of file G4INCLStandardPropagationModel.hh.

Constructor & Destructor Documentation

◆ StandardPropagationModel()

G4INCL::StandardPropagationModel::StandardPropagationModel ( LocalEnergyType  localEnergyType,
LocalEnergyType  localEnergyDeltaType,
const G4double  hTime = 0.0 
)

◆ ~StandardPropagationModel()

G4INCL::StandardPropagationModel::~StandardPropagationModel ( )
virtual

Definition at line 74 of file G4INCLStandardPropagationModel.cc.

75 {
76 delete theNucleus;
77 }

References theNucleus.

Member Function Documentation

◆ generateAllAvatars()

void G4INCL::StandardPropagationModel::generateAllAvatars ( )

(Re)Generate all possible avatars.

Definition at line 437 of file G4INCLStandardPropagationModel.cc.

437 {
438 ParticleList const &particles = theNucleus->getStore()->getParticles();
439// assert(!particles.empty());
440 G4double time;
441 for(ParticleIter i=particles.begin(), e=particles.end(); i!=e; ++i) {
442 time = this->getReflectionTime(*i);
443 if(time <= maximumTime) registerAvatar(new SurfaceAvatar(*i, time, theNucleus));
444 }
445 generateCollisions(particles);
446 generateDecays(particles);
447 }
double G4double
Definition: G4Types.hh:83
Store * getStore() const
G4double getReflectionTime(G4INCL::Particle const *const aParticle)
Get the reflection time.
void generateCollisions(const ParticleList &particles)
Generate and register collisions among particles in a list, except between those in another list.
void registerAvatar(G4INCL::IAvatar *anAvatar)
void generateDecays(const ParticleList &particles)
Generate decays for particles that can decay.
ParticleList const & getParticles() const
Definition: G4INCLStore.hh:253
ParticleList::const_iterator ParticleIter

References generateCollisions(), generateDecays(), G4INCL::Store::getParticles(), getReflectionTime(), G4INCL::Nucleus::getStore(), maximumTime, registerAvatar(), and theNucleus.

Referenced by shootComposite(), and shootParticle().

◆ generateBinaryCollisionAvatar()

IAvatar * G4INCL::StandardPropagationModel::generateBinaryCollisionAvatar ( Particle *const  p1,
Particle *const  p2 
)

Generate a two-particle avatar.

Generate a two-particle avatar, if all the appropriate conditions are met.

Definition at line 265 of file G4INCLStandardPropagationModel.cc.

265 {
266 // Is either particle a participant?
267 if(!p1->isParticipant() && !p2->isParticipant() && p1->getParticipantType()==p2->getParticipantType()) return NULL;
268
269 // Is it a pi-resonance collision (we don't treat them)?
270 if((p1->isResonance() && p2->isPion()) || (p1->isPion() && p2->isResonance()))
271 return NULL;
272
273 // Will the avatar take place between now and the end of the cascade?
274 G4double minDistOfApproachSquared = 0.0;
275 G4double t = getTime(p1, p2, &minDistOfApproachSquared);
276 if(t>maximumTime || t<currentTime+hadronizationTime) return NULL;
277
278 // Local energy. Jump through some hoops to calculate the cross section
279 // at the collision point, and clean up after yourself afterwards.
280 G4bool hasLocalEnergy;
281 if(p1->isPion() || p2->isPion())
285 else
286 hasLocalEnergy = ((theLocalEnergyType == FirstCollisionLocalEnergy &&
289 const G4bool p1HasLocalEnergy = (hasLocalEnergy && !p1->isMeson());
290 const G4bool p2HasLocalEnergy = (hasLocalEnergy && !p2->isMeson());
291
292 if(p1HasLocalEnergy) {
293 backupParticle1 = *p1;
294 p1->propagate(t - currentTime);
295 if(p1->getPosition().mag() > theNucleus->getSurfaceRadius(p1)) {
296 *p1 = backupParticle1;
297 return NULL;
298 }
300 }
301 if(p2HasLocalEnergy) {
302 backupParticle2 = *p2;
303 p2->propagate(t - currentTime);
304 if(p2->getPosition().mag() > theNucleus->getSurfaceRadius(p2)) {
305 *p2 = backupParticle2;
306 if(p1HasLocalEnergy) {
307 *p1 = backupParticle1;
308 }
309 return NULL;
310 }
312 }
313
314 // Compute the total cross section
315 const G4double totalCrossSection = CrossSections::total(p1, p2);
317
318 // Restore particles to their state before the local-energy tweak
319 if(p1HasLocalEnergy) {
320 *p1 = backupParticle1;
321 }
322 if(p2HasLocalEnergy) {
323 *p2 = backupParticle2;
324 }
325
326 // Is the CM energy > cutNN? (no cutNN on the first collision)
328 && p1->isNucleon() && p2->isNucleon()
330
331 // Do the particles come close enough to each other?
332 if(Math::tenPi*minDistOfApproachSquared > totalCrossSection) return NULL;
333
334 // Bomb out if the two collision partners are the same particle
335// assert(p1->getID() != p2->getID());
336
337 // Return a new avatar, then!
338 return new G4INCL::BinaryCollisionAvatar(t, totalCrossSection, theNucleus, p1, p2);
339 }
bool G4bool
Definition: G4Types.hh:86
G4int getAcceptedCollisions() const
Definition: G4INCLBook.hh:100
G4double getSurfaceRadius(Particle const *const particle) const
Get the maximum allowed radius for a given particle.
void propagate(G4double step)
G4double getTime(G4INCL::Particle const *const particleA, G4INCL::Particle const *const particleB, G4double *minDistOfApproach) const
Book & getBook()
Definition: G4INCLStore.hh:259
G4double total(Particle const *const p1, Particle const *const p2)
G4double squareTotalEnergyInCM(Particle const *const p1, Particle const *const p2)
void transformToLocalEnergyFrame(Nucleus const *const n, Particle *const p)
const G4double tenPi
@ FirstCollisionLocalEnergy

References G4INCL::AlwaysLocalEnergy, backupParticle1, backupParticle2, currentTime, G4INCL::FirstCollisionLocalEnergy, G4INCL::Book::getAcceptedCollisions(), G4INCL::Store::getBook(), G4INCL::BinaryCollisionAvatar::getCutNNSquared(), G4INCL::Particle::getParticipantType(), G4INCL::Particle::getPosition(), G4INCL::Nucleus::getStore(), G4INCL::Nucleus::getSurfaceRadius(), getTime(), hadronizationTime, G4INCL::Particle::isMeson(), G4INCL::Particle::isNucleon(), G4INCL::Particle::isParticipant(), G4INCL::Particle::isPion(), G4INCL::Particle::isResonance(), G4INCL::ThreeVector::mag(), maximumTime, G4INCL::Particle::propagate(), G4INCL::KinematicsUtils::squareTotalEnergyInCM(), G4INCL::Math::tenPi, theLocalEnergyDeltaType, theLocalEnergyType, theNucleus, G4INCL::CrossSections::total(), and G4INCL::KinematicsUtils::transformToLocalEnergyFrame().

Referenced by generateCollisions(), and generateUpdatedCollisions().

◆ generateCollisions() [1/2]

void G4INCL::StandardPropagationModel::generateCollisions ( const ParticleList particles)

Generate and register collisions among particles in a list, except between those in another list.

This method generates all possible collisions among the particles. Each collision is generated only once.

Parameters
particleslist of particles

Definition at line 397 of file G4INCLStandardPropagationModel.cc.

397 {
398 // Loop over all the particles
399 for(ParticleIter p1=particles.begin(), e=particles.end(); p1!=e; ++p1) {
400 // Loop over the rest of the particles
401 for(ParticleIter p2 = p1 + 1; p2 != particles.end(); ++p2) {
403 }
404 }
405 }
IAvatar * generateBinaryCollisionAvatar(Particle *const p1, Particle *const p2)
Generate a two-particle avatar.

References generateBinaryCollisionAvatar(), and registerAvatar().

Referenced by generateAllAvatars().

◆ generateCollisions() [2/2]

void G4INCL::StandardPropagationModel::generateCollisions ( const ParticleList particles,
const ParticleList except 
)

Generate and register collisions among particles in a list, except between those in another list.

This method generates all possible collisions among the particles. Each collision is generated only once. The collision is NOT generated if BOTH collision partners belong to the except list.

You should pass an empty list as the except parameter if you want to generate all possible collisions among particles.

Parameters
particleslist of particles
exceptlist of excluded particles

Definition at line 407 of file G4INCLStandardPropagationModel.cc.

407 {
408
409 const G4bool haveExcept = !except.empty();
410
411 // Loop over all the particles
412 for(ParticleIter p1=particles.begin(), e=particles.end(); p1!=e; ++p1)
413 {
414 // Loop over the rest of the particles
415 ParticleIter p2 = p1;
416 for(++p2; p2 != particles.end(); ++p2)
417 {
418 // Skip the collision if both particles must be excluded
419 if(haveExcept && except.contains(*p1) && except.contains(*p2)) continue;
420
422 }
423 }
424
425 }

References G4INCL::UnorderedVector< T >::contains(), generateBinaryCollisionAvatar(), and registerAvatar().

◆ generateDecays()

void G4INCL::StandardPropagationModel::generateDecays ( const ParticleList particles)

Generate decays for particles that can decay.

The list of particles given as an argument is allowed to contain also stable particles.

Parameters
particleslist of particles to (possibly) generate decays for

Definition at line 465 of file G4INCLStandardPropagationModel.cc.

465 {
466 for(ParticleIter i=particles.begin(), e=particles.end(); i!=e; ++i) {
467 if((*i)->isDelta()) {
468 G4double decayTime = DeltaDecayChannel::computeDecayTime((*i)); // time in fm/c
469 G4double time = currentTime + decayTime;
470 if(time <= maximumTime) {
471 registerAvatar(new DecayAvatar((*i), time, theNucleus));
472 }
473 }
474 else if((*i)->getType() == SigmaZero) {
475 G4double decayTime = SigmaZeroDecayChannel::computeDecayTime((*i)); // time in fm/c
476 G4double time = currentTime + decayTime;
477 if(time <= maximumTime) {
478 registerAvatar(new DecayAvatar((*i), time, theNucleus));
479 }
480 }
481 if((*i)->isOmega()) {
483 G4double timeOmega = currentTime + decayTimeOmega;
484 if(timeOmega <= maximumTime) {
485 registerAvatar(new DecayAvatar((*i), timeOmega, theNucleus));
486 }
487 }
488 }
489 }
static G4double computeDecayTime(Particle *p)
static G4double computeDecayTime(Particle *p)

References G4INCL::DeltaDecayChannel::computeDecayTime(), G4INCL::PionResonanceDecayChannel::computeDecayTime(), G4INCL::SigmaZeroDecayChannel::computeDecayTime(), currentTime, maximumTime, registerAvatar(), G4INCL::SigmaZero, and theNucleus.

Referenced by generateAllAvatars(), and propagate().

◆ generateUpdatedCollisions()

void G4INCL::StandardPropagationModel::generateUpdatedCollisions ( const ParticleList updatedParticles,
const ParticleList particles 
)

Generate and register collisions between a list of updated particles and all the other particles.

This method does not generate collisions among the particles in updatedParticles; in other words, it generates a collision between one of the updatedParticles and one of the particles ONLY IF the latter does not belong to updatedParticles.

If you intend to generate all possible collisions among particles in a list, use generateCollisions().

Parameters
updatedParticleslist of updated particles
particleslist of particles

Definition at line 378 of file G4INCLStandardPropagationModel.cc.

378 {
379
380 // Loop over all the updated particles
381 for(ParticleIter updated=updatedParticles.begin(), e=updatedParticles.end(); updated!=e; ++updated)
382 {
383 // Loop over all the particles
384 for(ParticleIter particle=particles.begin(), end=particles.end(); particle!=end; ++particle)
385 {
386 /* Consider the generation of a collision avatar only if (*particle)
387 * is not one of the updated particles.
388 * The criterion makes sure that you don't generate avatars between
389 * updated particles. */
390 if(updatedParticles.contains(*particle)) continue;
391
393 }
394 }
395 }

References G4INCL::UnorderedVector< T >::contains(), generateBinaryCollisionAvatar(), and registerAvatar().

Referenced by updateAvatars().

◆ getCurrentTime()

G4double G4INCL::StandardPropagationModel::getCurrentTime ( )
virtual

Returns the current global time of the system.

Implements G4INCL::IPropagationModel.

Definition at line 251 of file G4INCLStandardPropagationModel.cc.

251 {
252 return currentTime;
253 }

References currentTime.

◆ getNucleus()

G4INCL::Nucleus * G4INCL::StandardPropagationModel::getNucleus ( )
virtual

Get the nucleus.

Implements G4INCL::IPropagationModel.

Definition at line 79 of file G4INCLStandardPropagationModel.cc.

80 {
81 return theNucleus;
82 }

References theNucleus.

◆ getReflectionTime()

G4double G4INCL::StandardPropagationModel::getReflectionTime ( G4INCL::Particle const *const  aParticle)

Get the reflection time.

Returns the reflection time of a particle on the potential wall.

Parameters
aParticlepointer to the particle

Definition at line 341 of file G4INCLStandardPropagationModel.cc.

341 {
342 Intersection theIntersection(
344 aParticle->getPosition(),
345 aParticle->getPropagationVelocity(),
346 theNucleus->getSurfaceRadius(aParticle)));
347 G4double time;
348 if(theIntersection.exists) {
349 time = currentTime + theIntersection.time;
350 } else {
351 INCL_ERROR("Imaginary reflection time for particle: " << '\n'
352 << aParticle->print());
353 time = 10000.0;
354 }
355 return time;
356 }
#define INCL_ERROR(x)
Intersection getLaterTrajectoryIntersection(const ThreeVector &x0, const ThreeVector &p, const G4double r)
Compute the second intersection of a straight particle trajectory with a sphere.

References currentTime, G4INCL::Intersection::exists, G4INCL::IntersectionFactory::getLaterTrajectoryIntersection(), G4INCL::Particle::getPosition(), G4INCL::Particle::getPropagationVelocity(), G4INCL::Nucleus::getSurfaceRadius(), INCL_ERROR, G4INCL::Particle::print(), theNucleus, and G4INCL::Intersection::time.

Referenced by generateAllAvatars(), and updateAvatars().

◆ getStoppingTime()

G4double G4INCL::StandardPropagationModel::getStoppingTime ( )
virtual

Get the current stopping time.

Implements G4INCL::IPropagationModel.

Definition at line 242 of file G4INCLStandardPropagationModel.cc.

242 {
243 return maximumTime;
244 }

References maximumTime.

◆ getTime()

G4double G4INCL::StandardPropagationModel::getTime ( G4INCL::Particle const *const  particleA,
G4INCL::Particle const *const  particleB,
G4double minDistOfApproach 
) const

Get the predicted time of the collision between two particles.

Definition at line 358 of file G4INCLStandardPropagationModel.cc.

360 {
361 G4double time;
362 G4INCL::ThreeVector t13 = particleA->getPropagationVelocity();
363 t13 -= particleB->getPropagationVelocity();
364 G4INCL::ThreeVector distance = particleA->getPosition();
365 distance -= particleB->getPosition();
366 const G4double t7 = t13.dot(distance);
367 const G4double dt = t13.mag2();
368 if(dt <= 1.0e-10) {
369 (*minDistOfApproach) = 100000.0;
370 return currentTime + 100000.0;
371 } else {
372 time = -t7/dt;
373 }
374 (*minDistOfApproach) = distance.mag2() + time * t7;
375 return currentTime + time;
376 }
G4double dot(const ThreeVector &v) const
G4double mag2() const

References currentTime, G4INCL::ThreeVector::dot(), G4INCL::Particle::getPosition(), G4INCL::Particle::getPropagationVelocity(), and G4INCL::ThreeVector::mag2().

Referenced by generateBinaryCollisionAvatar().

◆ propagate()

G4INCL::IAvatar * G4INCL::StandardPropagationModel::propagate ( FinalState const *const  fs)
virtual

Propagate all particles and return the first avatar.

Implements G4INCL::IPropagationModel.

Definition at line 491 of file G4INCLStandardPropagationModel.cc.

492 {
493 if(fs) {
494 // We update only the information related to particles that were updated
495 // by the previous avatar.
496#ifdef INCL_REGENERATE_AVATARS
497#warning "The INCL_REGENERATE_AVATARS code has not been tested in a while. Use it at your peril."
498 if(!fs->getModifiedParticles().empty() || !fs->getEnteringParticles().empty() || !fs->getCreatedParticles().empty()) {
499 // Regenerates the entire avatar list, skipping collisions between
500 // updated particles
502 theNucleus->getStore()->initialiseParticleAvatarConnections();
503 generateAllAvatarsExceptUpdated(fs);
504 }
505#else
506 ParticleList const &updatedParticles = fs->getModifiedParticles();
507 if(fs->getValidity()==PauliBlockedFS) {
508 // This final state might represents the outcome of a Pauli-blocked delta
509 // decay
510// assert(updatedParticles.empty() || (updatedParticles.size()==1 && updatedParticles.front()->isResonance()));
511// assert(fs->getEnteringParticles().empty() && fs->getCreatedParticles().empty() && fs->getOutgoingParticles().empty() && fs->getDestroyedParticles().empty());
512 generateDecays(updatedParticles);
513 } else {
514 ParticleList const &entering = fs->getEnteringParticles();
515 generateDecays(updatedParticles);
516 generateDecays(entering);
517
518 ParticleList const &created = fs->getCreatedParticles();
519 if(created.empty() && entering.empty())
520 updateAvatars(updatedParticles);
521 else {
522 ParticleList updatedParticlesCopy = updatedParticles;
523 updatedParticlesCopy.insert(updatedParticlesCopy.end(), entering.begin(), entering.end());
524 updatedParticlesCopy.insert(updatedParticlesCopy.end(), created.begin(), created.end());
525 updateAvatars(updatedParticlesCopy);
526 }
527 }
528#endif
529 }
530
532 if(theAvatar == 0) return 0; // Avatar list is empty
533 // theAvatar->dispose();
534
535 if(theAvatar->getTime() < currentTime) {
536 INCL_ERROR("Avatar time = " << theAvatar->getTime() << ", currentTime = " << currentTime << '\n');
537 return 0;
538 } else if(theAvatar->getTime() > currentTime) {
539 theNucleus->getStore()->timeStep(theAvatar->getTime() - currentTime);
540
541 currentTime = theAvatar->getTime();
543 }
544
545 return theAvatar;
546 }
void setCurrentTime(G4double t)
Definition: G4INCLBook.hh:97
G4double getTime() const
void updateAvatars(const ParticleList &particles)
void clearAvatars()
Definition: G4INCLStore.cc:193
void timeStep(G4double step)
Definition: G4INCLStore.cc:168
IAvatar * findSmallestTime()
Definition: G4INCLStore.cc:142

References G4INCL::Store::clearAvatars(), currentTime, G4INCL::Store::findSmallestTime(), generateDecays(), G4INCL::Store::getBook(), G4INCL::FinalState::getCreatedParticles(), G4INCL::FinalState::getEnteringParticles(), G4INCL::FinalState::getModifiedParticles(), G4INCL::Nucleus::getStore(), G4INCL::IAvatar::getTime(), G4INCL::FinalState::getValidity(), INCL_ERROR, G4INCL::PauliBlockedFS, G4INCL::Book::setCurrentTime(), theNucleus, G4INCL::Store::timeStep(), and updateAvatars().

◆ registerAvatar()

void G4INCL::StandardPropagationModel::registerAvatar ( G4INCL::IAvatar anAvatar)

Add an avatar to the storage.

Definition at line 260 of file G4INCLStandardPropagationModel.cc.

261 {
262 if(anAvatar) theNucleus->getStore()->add(anAvatar);
263 }
void add(Particle *p)
Definition: G4INCLStore.cc:58

References G4INCL::Store::add(), G4INCL::Nucleus::getStore(), and theNucleus.

Referenced by generateAllAvatars(), generateCollisions(), generateDecays(), generateUpdatedCollisions(), and updateAvatars().

◆ setNucleus()

void G4INCL::StandardPropagationModel::setNucleus ( G4INCL::Nucleus nucleus)
virtual

Set the nucleus for this propagation model.

Implements G4INCL::IPropagationModel.

Definition at line 255 of file G4INCLStandardPropagationModel.cc.

256 {
257 theNucleus = nucleus;
258 }

References theNucleus.

◆ setStoppingTime()

void G4INCL::StandardPropagationModel::setStoppingTime ( G4double  time)
virtual

Set the stopping time of the simulation.

Implements G4INCL::IPropagationModel.

Definition at line 246 of file G4INCLStandardPropagationModel.cc.

246 {
247// assert(time>0.0);
248 maximumTime = time;
249 }

References maximumTime.

◆ shoot()

G4double G4INCL::StandardPropagationModel::shoot ( ParticleSpecies const &  projectileSpecies,
const G4double  kineticEnergy,
const G4double  impactParameter,
const G4double  phi 
)
virtual

Implements G4INCL::IPropagationModel.

Definition at line 84 of file G4INCLStandardPropagationModel.cc.

84 {
85 if(projectileSpecies.theType==Composite)
86 return shootComposite(projectileSpecies, kineticEnergy, impactParameter, phi);
87 else
88 return shootParticle(projectileSpecies.theType, kineticEnergy, impactParameter, phi);
89 }
G4double shootComposite(ParticleSpecies const &s, const G4double kineticEnergy, const G4double impactParameter, const G4double phi)
G4double shootParticle(ParticleType const t, const G4double kineticEnergy, const G4double impactParameter, const G4double phi)

References G4INCL::Composite, shootComposite(), shootParticle(), and G4INCL::ParticleSpecies::theType.

◆ shootComposite()

G4double G4INCL::StandardPropagationModel::shootComposite ( ParticleSpecies const &  s,
const G4double  kineticEnergy,
const G4double  impactParameter,
const G4double  phi 
)
virtual

Implements G4INCL::IPropagationModel.

Definition at line 168 of file G4INCLStandardPropagationModel.cc.

168 {
170 currentTime = 0.0;
171
172 // Create the ProjectileRemnant object
173 ProjectileRemnant *pr = new ProjectileRemnant(species, kineticEnergy);
174
175 // Same stopping time as for nucleon-nucleus
176 maximumTime = 29.8 * std::pow(theNucleus->getA(), 0.16);
177
178 // If the incoming cluster is slow, use a larger stopping time
179 const G4double rms = ParticleTable::getLargestNuclearRadius(pr->getA(), pr->getZ());
180 const G4double rMax = theNucleus->getUniverseRadius();
181 const G4double distance = 2.*rMax + 2.725*rms;
182 const G4double projectileVelocity = pr->boostVector().mag();
183 const G4double traversalTime = distance / projectileVelocity;
184 if(maximumTime < traversalTime)
185 maximumTime = traversalTime;
186 INCL_DEBUG("Cascade stopping time is " << maximumTime << '\n');
187
188 // If Coulomb is activated, do not process events with impact
189 // parameter larger than the maximum impact parameter, taking into
190 // account Coulomb distortion.
191 if(impactParameter>CoulombDistortion::maxImpactParameter(pr,theNucleus)) {
192 INCL_DEBUG("impactParameter>CoulombDistortion::maxImpactParameter" << '\n');
193 delete pr;
194 return -1.;
195 }
196
197 // Position the cluster at the right impact parameter
198 ThreeVector position(impactParameter * std::cos(phi),
199 impactParameter * std::sin(phi),
200 0.);
201 pr->setPosition(position);
202
203 // Fill in the relevant kinematic variables
204 theNucleus->setIncomingAngularMomentum(pr->getAngularMomentum());
205 theNucleus->setIncomingMomentum(pr->getMomentum());
206 theNucleus->setInitialEnergy(pr->getEnergy()
208
210 firstAvatar = false;
211
212 // Get the entry avatars from Coulomb
213 IAvatarList theAvatarList
215
216 if(theAvatarList.empty()) {
217 INCL_DEBUG("No ParticleEntryAvatar found, transparent event" << '\n');
218 delete pr;
219 return -1.;
220 }
221
222 /* Store the internal kinematics of the projectile remnant.
223 *
224 * Note that this is at variance with the Fortran version, which stores
225 * the initial kinematics of the particles *after* putting the spectators
226 * on mass shell, but *before* removing the same energy from the
227 * participants. Due to the different code flow, doing so in the C++
228 * version leads to wrong excitation energies for the forced compound
229 * nucleus.
230 */
231 pr->storeComponents();
232
233 // Tell the Nucleus about the ProjectileRemnant
235
236 // Register the ParticleEntryAvatars
238
239 return pr->getTransversePosition().mag();
240 }
#define INCL_DEBUG(x)
void setIncomingAngularMomentum(const ThreeVector &j)
Set the incoming angular-momentum vector.
void setInitialEnergy(const G4double e)
Set the initial energy.
void setIncomingMomentum(const ThreeVector &p)
Set the incoming momentum vector.
void setNucleusNucleusCollision()
Set a nucleus-nucleus collision.
void setProjectileRemnant(ProjectileRemnant *const c)
Set the projectile remnant.
G4double getUniverseRadius() const
Getter for theUniverseRadius.
G4int getS() const
Returns the strangeness number.
G4int getZ() const
Returns the charge number.
G4int getA() const
Returns the baryon number.
void generateAllAvatars()
(Re)Generate all possible avatars.
void addParticleEntryAvatars(IAvatarList const &al)
Add one ParticleEntry avatar.
Definition: G4INCLStore.cc:78
ParticleEntryAvatar * bringToSurface(Particle *p, Nucleus *const n)
Modify the momentum of an incoming particle and position it on the surface of the nucleus.
G4double maxImpactParameter(ParticleSpecies const &p, const G4double kinE, Nucleus const *const n)
Return the maximum impact parameter for Coulomb-distorted trajectories.
G4ThreadLocal NuclearMassFn getTableMass
Static pointer to the mass function for nuclei.
G4double getLargestNuclearRadius(const G4int A, const G4int Z)
UnorderedVector< IAvatar * > IAvatarList
#define position
Definition: xmlparse.cc:622

References G4INCL::Store::addParticleEntryAvatars(), G4INCL::Particle::boostVector(), G4INCL::CoulombDistortion::bringToSurface(), currentTime, firstAvatar, generateAllAvatars(), G4INCL::Particle::getA(), G4INCL::Cluster::getAngularMomentum(), G4INCL::Particle::getEnergy(), G4INCL::ParticleTable::getLargestNuclearRadius(), G4INCL::Particle::getMomentum(), G4INCL::Particle::getS(), G4INCL::Nucleus::getStore(), G4INCL::ParticleTable::getTableMass, G4INCL::Particle::getTransversePosition(), G4INCL::Nucleus::getUniverseRadius(), G4INCL::Particle::getZ(), INCL_DEBUG, G4INCL::ThreeVector::mag(), G4INCL::CoulombDistortion::maxImpactParameter(), maximumTime, position, G4INCL::Nucleus::setIncomingAngularMomentum(), G4INCL::Nucleus::setIncomingMomentum(), G4INCL::Nucleus::setInitialEnergy(), G4INCL::Nucleus::setNucleusNucleusCollision(), G4INCL::Cluster::setPosition(), G4INCL::Nucleus::setProjectileRemnant(), G4INCL::ProjectileRemnant::storeComponents(), and theNucleus.

Referenced by shoot().

◆ shootParticle()

G4double G4INCL::StandardPropagationModel::shootParticle ( ParticleType const  t,
const G4double  kineticEnergy,
const G4double  impactParameter,
const G4double  phi 
)
virtual

Implements G4INCL::IPropagationModel.

Definition at line 91 of file G4INCLStandardPropagationModel.cc.

91 {
93 currentTime = 0.0;
94
95 // Create the projectile particle
96 const G4double projectileMass = ParticleTable::getTableParticleMass(type);
97 G4double energy = kineticEnergy + projectileMass;
98 G4double momentumZ = std::sqrt(energy*energy - projectileMass*projectileMass);
99 ThreeVector momentum(0.0, 0.0, momentumZ);
100 Particle *p= new G4INCL::Particle(type, energy, momentum, ThreeVector());
101
102 G4double temfin;
103 G4double TLab;
104 if( p->isMeson()) {
105 temfin = 30.18 * std::pow(theNucleus->getA(), 0.17);
106 TLab = p->getKineticEnergy();
107 } else {
108 temfin = 29.8 * std::pow(theNucleus->getA(), 0.16);
109 TLab = p->getKineticEnergy()/p->getA();
110 }
111
112 // energy-dependent stopping time above 2 AGeV
113 if(TLab>2000.)
114 temfin *= (5.8E4-TLab)/5.6E4;
115
116 maximumTime = temfin;
117
118 // If the incoming particle is slow, use a larger stopping time
119 const G4double rMax = theNucleus->getUniverseRadius();
120 const G4double distance = 2.*rMax;
121 const G4double projectileVelocity = p->boostVector().mag();
122 const G4double traversalTime = distance / projectileVelocity;
123 if(maximumTime < traversalTime)
124 maximumTime = traversalTime;
125 INCL_DEBUG("Cascade stopping time is " << maximumTime << '\n');
126
127 // If Coulomb is activated, do not process events with impact
128 // parameter larger than the maximum impact parameter, taking into
129 // account Coulomb distortion.
130 if(impactParameter>CoulombDistortion::maxImpactParameter(p->getSpecies(), kineticEnergy, theNucleus)) {
131 INCL_DEBUG("impactParameter>CoulombDistortion::maxImpactParameter" << '\n');
132 delete p;
133 return -1.;
134 }
135
136 ThreeVector position(impactParameter * std::cos(phi),
137 impactParameter * std::sin(phi),
138 0.);
139 p->setPosition(position);
140
141 // Fill in the relevant kinematic variables
142 theNucleus->setIncomingAngularMomentum(p->getAngularMomentum());
143 theNucleus->setIncomingMomentum(p->getMomentum());
144 theNucleus->setInitialEnergy(p->getEnergy()
146
147 // Reset the particle kinematics to the INCL values
148 p->setINCLMass();
149 p->setEnergy(p->getMass() + kineticEnergy);
150 p->adjustMomentumFromEnergy();
151
152 p->makeProjectileSpectator();
154 firstAvatar = false;
155
156 // Get the entry avatars from Coulomb and put them in the Store
157 ParticleEntryAvatar *theEntryAvatar = CoulombDistortion::bringToSurface(p, theNucleus);
158 if(theEntryAvatar) {
159 theNucleus->getStore()->addParticleEntryAvatar(theEntryAvatar);
160
161 return p->getTransversePosition().mag();
162 } else {
163 delete p;
164 return -1.;
165 }
166 }
void setParticleNucleusCollision()
Set a particle-nucleus collision.
void addParticleEntryAvatar(IAvatar *a)
Add one ParticleEntry avatar.
Definition: G4INCLStore.cc:66
G4double energy(const ThreeVector &p, const G4double m)
G4ThreadLocal ParticleMassFn getTableParticleMass
Static pointer to the mass function for particles.

References G4INCL::Store::addParticleEntryAvatar(), G4INCL::Particle::adjustMomentumFromEnergy(), G4INCL::Particle::boostVector(), G4INCL::CoulombDistortion::bringToSurface(), currentTime, G4INCL::KinematicsUtils::energy(), firstAvatar, generateAllAvatars(), G4INCL::Particle::getA(), G4INCL::Particle::getAngularMomentum(), G4INCL::Particle::getEnergy(), G4INCL::Particle::getKineticEnergy(), G4INCL::Particle::getMass(), G4INCL::Particle::getMomentum(), G4INCL::Particle::getS(), G4INCL::Particle::getSpecies(), G4INCL::Nucleus::getStore(), G4INCL::ParticleTable::getTableMass, G4INCL::ParticleTable::getTableParticleMass, G4INCL::Particle::getTransversePosition(), G4INCL::Nucleus::getUniverseRadius(), G4INCL::Particle::getZ(), INCL_DEBUG, G4INCL::Particle::isMeson(), G4INCL::ThreeVector::mag(), G4INCL::Particle::makeProjectileSpectator(), G4INCL::CoulombDistortion::maxImpactParameter(), maximumTime, position, G4INCL::Particle::setEnergy(), G4INCL::Particle::setINCLMass(), G4INCL::Nucleus::setIncomingAngularMomentum(), G4INCL::Nucleus::setIncomingMomentum(), G4INCL::Nucleus::setInitialEnergy(), G4INCL::Nucleus::setParticleNucleusCollision(), G4INCL::Particle::setPosition(), and theNucleus.

Referenced by shoot().

◆ updateAvatars()

void G4INCL::StandardPropagationModel::updateAvatars ( const ParticleList particles)

Update all avatars related to a particle.

Definition at line 427 of file G4INCLStandardPropagationModel.cc.

427 {
428
429 for(ParticleIter iter=particles.begin(), e=particles.end(); iter!=e; ++iter) {
430 G4double time = this->getReflectionTime(*iter);
431 if(time <= maximumTime) registerAvatar(new SurfaceAvatar(*iter, time, theNucleus));
432 }
433 ParticleList const &p = theNucleus->getStore()->getParticles();
434 generateUpdatedCollisions(particles, p); // Predict collisions with spectators and participants
435 }
void generateUpdatedCollisions(const ParticleList &updatedParticles, const ParticleList &particles)
Generate and register collisions between a list of updated particles and all the other particles.

References generateUpdatedCollisions(), G4INCL::Store::getParticles(), getReflectionTime(), G4INCL::Nucleus::getStore(), maximumTime, registerAvatar(), and theNucleus.

Referenced by propagate().

Field Documentation

◆ backupParticle1

Particle G4INCL::StandardPropagationModel::backupParticle1
private

Definition at line 200 of file G4INCLStandardPropagationModel.hh.

Referenced by generateBinaryCollisionAvatar().

◆ backupParticle2

Particle G4INCL::StandardPropagationModel::backupParticle2
private

Definition at line 200 of file G4INCLStandardPropagationModel.hh.

Referenced by generateBinaryCollisionAvatar().

◆ currentTime

G4double G4INCL::StandardPropagationModel::currentTime
private

◆ firstAvatar

G4bool G4INCL::StandardPropagationModel::firstAvatar
private

Definition at line 198 of file G4INCLStandardPropagationModel.hh.

Referenced by shootComposite(), and shootParticle().

◆ hadronizationTime

G4double G4INCL::StandardPropagationModel::hadronizationTime
private

Definition at line 197 of file G4INCLStandardPropagationModel.hh.

Referenced by generateBinaryCollisionAvatar().

◆ maximumTime

G4double G4INCL::StandardPropagationModel::maximumTime
private

◆ theLocalEnergyDeltaType

LocalEnergyType G4INCL::StandardPropagationModel::theLocalEnergyDeltaType
private

Definition at line 199 of file G4INCLStandardPropagationModel.hh.

Referenced by generateBinaryCollisionAvatar().

◆ theLocalEnergyType

LocalEnergyType G4INCL::StandardPropagationModel::theLocalEnergyType
private

Definition at line 199 of file G4INCLStandardPropagationModel.hh.

Referenced by generateBinaryCollisionAvatar().

◆ theNucleus

G4INCL::Nucleus* G4INCL::StandardPropagationModel::theNucleus
private

The documentation for this class was generated from the following files: