Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4INCL::StandardPropagationModel Class Reference

#include <G4INCLStandardPropagationModel.hh>

+ Inheritance diagram for G4INCL::StandardPropagationModel:

Public Member Functions

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

Additional Inherited Members

virtual G4double shootParticle (ParticleType const t, const G4double kineticEnergy, const G4double impactParameter, const G4double phi)=0
 
virtual G4double shootComposite (ParticleSpecies const s, const G4double kineticEnergy, const G4double impactParameter, const G4double phi)=0
 

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 70 of file G4INCLStandardPropagationModel.hh.

Constructor & Destructor Documentation

◆ StandardPropagationModel()

G4INCL::StandardPropagationModel::StandardPropagationModel ( LocalEnergyType  localEnergyType,
LocalEnergyType  localEnergyDeltaType 
)

Definition at line 64 of file G4INCLStandardPropagationModel.cc.

65 :theNucleus(0), maximumTime(70.0), currentTime(0.0), firstAvatar(true),
66 theLocalEnergyType(localEnergyType),
67 theLocalEnergyDeltaType(localEnergyDeltaType)
68 {
69 }

◆ ~StandardPropagationModel()

G4INCL::StandardPropagationModel::~StandardPropagationModel ( )
virtual

Definition at line 71 of file G4INCLStandardPropagationModel.cc.

72 {
73 delete theNucleus;
74 }

Member Function Documentation

◆ generateAllAvatars()

void G4INCL::StandardPropagationModel::generateAllAvatars ( G4bool  excludeUpdated = false)

(Re)Generate all possible avatars.

Parameters
excludeUpdatedexclude collisions between updated particles.

Definition at line 428 of file G4INCLStandardPropagationModel.cc.

428 {
429 ParticleList particles = theNucleus->getStore()->getParticles();
430 if(particles.empty()) { ERROR("No particles inside the nucleus!" << std::endl); }
431 for(ParticleIter i = particles.begin(); i != particles.end(); ++i) {
432 G4double time = this->getReflectionTime(*i);
433 if(time <= maximumTime) registerAvatar(new SurfaceAvatar(*i, time, theNucleus));
434 }
435 ParticleList except;
436 if(excludeUpdated)
437 except = theNucleus->getUpdatedParticles();
438 generateCollisions(particles,except);
439 generateDecays(particles);
440 }
#define ERROR(x)
double G4double
Definition: G4Types.hh:64
Store * getStore() const
ParticleList const & getUpdatedParticles() const
G4double getReflectionTime(G4INCL::Particle const *const aParticle)
Get the reflection time.
void registerAvatar(G4INCL::IAvatar *anAvatar)
void generateCollisions(const ParticleList &particles, const ParticleList &except)
Generate and register collisions among particles in a list, except between those in another list.
void generateDecays(const ParticleList &particles)
Generate decays for particles that can decay.
ParticleList const & getParticles() const
Definition: G4INCLStore.hh:237
std::list< G4INCL::Particle * > ParticleList
std::list< G4INCL::Particle * >::const_iterator ParticleIter

Referenced by propagate(), shootComposite(), and shootParticle().

◆ generateBinaryCollisionAvatar()

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

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) 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())
282 hasLocalEnergy = ((theLocalEnergyDeltaType == FirstCollisionLocalEnergy &&
283 theNucleus->getStore()->getBook()->getAcceptedCollisions()==0) ||
284 theLocalEnergyDeltaType == AlwaysLocalEnergy);
285 else
286 hasLocalEnergy = ((theLocalEnergyType == FirstCollisionLocalEnergy &&
287 theNucleus->getStore()->getBook()->getAcceptedCollisions()==0) ||
288 theLocalEnergyType == AlwaysLocalEnergy);
289 const G4bool p1HasLocalEnergy = (hasLocalEnergy && !p1->isPion());
290 const G4bool p2HasLocalEnergy = (hasLocalEnergy && !p2->isPion());
291
292 Particle backupParticle1 = *p1;
293 if(p1HasLocalEnergy) {
294 p1->propagate(t - currentTime);
295 if(p1->getPosition().mag() > theNucleus->getSurfaceRadius(p1)) {
296 *p1 = backupParticle1;
297 return NULL;
298 }
300 }
301 Particle backupParticle2 = *p2;
302 if(p2HasLocalEnergy) {
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);
316 const G4double squareTotalEnergyInCM = KinematicsUtils::squareTotalEnergyInCM(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)
327 if(theNucleus->getStore()->getBook()->getAcceptedCollisions()>0
328 && p1->isNucleon() && p2->isNucleon()
329 && squareTotalEnergyInCM < BinaryCollisionAvatar::cutNNSquared) return NULL;
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:67
G4int getAcceptedCollisions() const
Definition: G4INCLBook.hh:85
static G4double total(Particle const *const p1, Particle const *const p2)
static G4double squareTotalEnergyInCM(Particle const *const p1, Particle const *const p2)
static void transformToLocalEnergyFrame(Nucleus const *const n, Particle *const p)
G4double getSurfaceRadius(Particle const *const particle) const
Get the maximum allowed radius for a given particle.
G4double getTime(G4INCL::Particle const *const particleA, G4INCL::Particle const *const particleB, G4double *minDistOfApproach) const
Book * getBook()
Definition: G4INCLStore.hh:243
const G4double tenPi
@ FirstCollisionLocalEnergy

Referenced by generateCollisions(), and generateUpdatedCollisions().

◆ generateCollisions()

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 397 of file G4INCLStandardPropagationModel.cc.

397 {
398
399 G4bool haveExcept;
400 haveExcept=(except.size()!=0);
401
402 // Loop over all the particles
403 for(ParticleIter p1 = particles.begin(); p1 != particles.end(); ++p1)
404 {
405 // Loop over the rest of the particles
406 ParticleIter p2 = p1;
407 for(++p2; p2 != particles.end(); ++p2)
408 {
409 // Skip the collision if both particles must be excluded
410 if(haveExcept && (*p1)->isInList(except) && (*p2)->isInList(except)) continue;
411
413 }
414 }
415
416 }
IAvatar * generateBinaryCollisionAvatar(Particle *const p1, Particle *const p2) const
Generate a two-particle avatar.

Referenced by generateAllAvatars().

◆ 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 442 of file G4INCLStandardPropagationModel.cc.

442 {
443 for(ParticleIter i = particles.begin(); i != particles.end(); ++i) {
444 if((*i)->isDelta()) {
446 G4double time = currentTime + decayTime;
447 if(time <= maximumTime) {
448 registerAvatar(new DecayAvatar((*i), time, theNucleus));
449 }
450 }
451 }
452 }
static G4double computeDecayTime(Particle *p)

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(); updated != updatedParticles.end(); ++updated)
382 {
383 // Loop over all the particles
384 for(ParticleIter particle = particles.begin(); particle != particles.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((*particle)->isInList(updatedParticles)) continue;
391
393 }
394 }
395 }

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 }

◆ getNucleus()

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

Get the nucleus.

Implements G4INCL::IPropagationModel.

Definition at line 76 of file G4INCLStandardPropagationModel.cc.

77 {
78 return theNucleus;
79 }

◆ 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 ERROR("Imaginary reflection time for particle: " << std::endl
352 << aParticle->print());
353 time = 10000.0;
354 }
355 return time;
356 }
static Intersection getLaterTrajectoryIntersection(const ThreeVector &x0, const ThreeVector &p, const G4double r)

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 }

◆ 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

Referenced by generateBinaryCollisionAvatar().

◆ propagate()

G4INCL::IAvatar * G4INCL::StandardPropagationModel::propagate ( )
virtual

Propagate all particles and return the first avatar.

Implements G4INCL::IPropagationModel.

Definition at line 454 of file G4INCLStandardPropagationModel.cc.

455 {
456 // We update only the information related to particles that were updated
457 // by the previous avatar.
458#ifdef INCL_REGENERATE_AVATARS
459#warning "The INCL_REGENERATE_AVATARS code has not been tested in a while. Use it at your peril."
460 if(theNucleus->getUpdatedParticles().size()!=0 || theNucleus->getCreatedParticles().size()!=0) {
461 // Regenerates the entire avatar list, skipping collisions between
462 // updated particles
463 theNucleus->getStore()->clearAvatars();
465 generateAllAvatars(true);
466 }
467#else
468 // Deltas are created by transforming nucleon into a delta for
469 // efficiency reasons
470 Particle * const blockedDelta = theNucleus->getBlockedDelta();
471 ParticleList updatedParticles = theNucleus->getUpdatedParticles();
472 if(blockedDelta)
473 updatedParticles.push_back(blockedDelta);
474 generateDecays(updatedParticles);
475
476 ParticleList needNewAvatars = theNucleus->getUpdatedParticles();
477 ParticleList created = theNucleus->getCreatedParticles();
478 needNewAvatars.splice(needNewAvatars.end(), created);
479 updateAvatars(needNewAvatars);
480#endif
481
482 G4INCL::IAvatar *theAvatar = theNucleus->getStore()->findSmallestTime();
483 if(theAvatar == 0) return 0; // Avatar list is empty
484 // theAvatar->dispose();
485
486 if(theAvatar->getTime() < currentTime) {
487 ERROR("Avatar time = " << theAvatar->getTime() << ", currentTime = " << currentTime << std::endl);
488 return 0;
489 } else if(theAvatar->getTime() > currentTime) {
490 theNucleus->getStore()->timeStep(theAvatar->getTime() - currentTime);
491
492 currentTime = theAvatar->getTime();
493 theNucleus->getStore()->getBook()->setCurrentTime(currentTime);
494 }
495
496 return theAvatar;
497 }
void setCurrentTime(G4double t)
Definition: G4INCLBook.hh:82
G4double getTime() const
Particle * getBlockedDelta() const
Get the delta that could not decay.
ParticleList const & getCreatedParticles() const
void generateAllAvatars(G4bool excludeUpdated=false)
(Re)Generate all possible avatars.
void updateAvatars(const ParticleList &particles)
void clearAvatars()
Definition: G4INCLStore.cc:318
void timeStep(G4double step)
Definition: G4INCLStore.cc:263
IAvatar * findSmallestTime()
Definition: G4INCLStore.cc:212
void initialiseParticleAvatarConnections()
Initialise the particleAvatarConnections map.
Definition: G4INCLStore.cc:335

◆ 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:62

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 }

◆ 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 }

◆ shoot()

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

Implements G4INCL::IPropagationModel.

Definition at line 81 of file G4INCLStandardPropagationModel.cc.

81 {
82 if(projectileSpecies.theType==Composite)
83 return shootComposite(projectileSpecies, kineticEnergy, impactParameter, phi);
84 else
85 return shootParticle(projectileSpecies.theType, kineticEnergy, impactParameter, phi);
86 }
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)

◆ shootComposite()

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

Implements G4INCL::IPropagationModel.

Definition at line 162 of file G4INCLStandardPropagationModel.cc.

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

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 88 of file G4INCLStandardPropagationModel.cc.

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

Referenced by shoot().

◆ updateAvatars()

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

Update all avatars related to a particle.

Definition at line 418 of file G4INCLStandardPropagationModel.cc.

418 {
419
420 for(ParticleIter iter = particles.begin(); iter != particles.end(); ++iter) {
421 G4double time = this->getReflectionTime(*iter);
422 if(time <= maximumTime) registerAvatar(new SurfaceAvatar(*iter, time, theNucleus));
423 }
424 ParticleList const &p = theNucleus->getStore()->getParticles();
425 generateUpdatedCollisions(particles, p); // Predict collisions with spectators and participants
426 }
void generateUpdatedCollisions(const ParticleList &updatedParticles, const ParticleList &particles)
Generate and register collisions between a list of updated particles and all the other particles.

Referenced by propagate().


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