Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4INCLStandardPropagationModel.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // INCL++ intra-nuclear cascade model
27 // Pekka Kaitaniemi, CEA and Helsinki Institute of Physics
28 // Davide Mancusi, CEA
29 // Alain Boudard, CEA
30 // Sylvie Leray, CEA
31 // Joseph Cugnon, University of Liege
32 //
33 #define INCLXX_IN_GEANT4_MODE 1
34 
35 #include "globals.hh"
36 
37 /*
38  * StandardPropagationModel.cpp
39  *
40  * \date 4 juin 2009
41  * \author Pekka Kaitaniemi
42  */
43 
45 #include "G4INCLSurfaceAvatar.hh"
47 #include "G4INCLDecayAvatar.hh"
48 #include "G4INCLCrossSections.hh"
49 #include "G4INCLRandom.hh"
50 #include <iostream>
51 #include <algorithm>
52 #include "G4INCLLogger.hh"
53 #include "G4INCLGlobals.hh"
54 #include "G4INCLKinematicsUtils.hh"
58 #include "G4INCLIntersection.hh"
59 
60 namespace G4INCL {
61 
63  :theNucleus(0), maximumTime(70.0), currentTime(0.0), firstAvatar(true),
64  theLocalEnergyType(localEnergyType),
65  theLocalEnergyDeltaType(localEnergyDeltaType)
66  {
67  }
68 
70  {
71  delete theNucleus;
72  }
73 
75  {
76  return theNucleus;
77  }
78 
79  G4double StandardPropagationModel::shoot(ParticleSpecies const projectileSpecies, const G4double kineticEnergy, const G4double impactParameter, const G4double phi) {
80  if(projectileSpecies.theType==Composite)
81  return shootComposite(projectileSpecies, kineticEnergy, impactParameter, phi);
82  else
83  return shootParticle(projectileSpecies.theType, kineticEnergy, impactParameter, phi);
84  }
85 
86  G4double StandardPropagationModel::shootParticle(ParticleType const type, const G4double kineticEnergy, const G4double impactParameter, const G4double phi) {
87  theNucleus->setParticleNucleusCollision();
88  currentTime = 0.0;
89 
90  // Create the projectile particle
91  const G4double projectileMass = ParticleTable::getTableParticleMass(type);
92  G4double energy = kineticEnergy + projectileMass;
93  G4double momentumZ = std::sqrt(energy*energy - projectileMass*projectileMass);
94  ThreeVector momentum(0.0, 0.0, momentumZ);
95  Particle *p= new G4INCL::Particle(type, energy, momentum, ThreeVector());
96 
97  G4double temfin = 0.0;
98  if( p->isNucleon() )
99  temfin = 29.8 * std::pow(theNucleus->getA(), 0.16);
100  else {
101  const G4double tlab = p->getEnergy() - p->getMass();
102  temfin = 30.18 * std::pow(theNucleus->getA(), 0.17*(1.0 - 5.7E-5*tlab));
103  }
104 
105  maximumTime = temfin;
106 
107  // If the incoming particle is slow, use a larger stopping time
108  const G4double rMax = theNucleus->getUniverseRadius();
109  const G4double distance = 2.*rMax;
110  const G4double projectileVelocity = p->boostVector().mag();
111  const G4double traversalTime = distance / projectileVelocity;
112  if(maximumTime < traversalTime)
113  maximumTime = traversalTime;
114  INCL_DEBUG("Cascade stopping time is " << maximumTime << std::endl);
115 
116  // If Coulomb is activated, do not process events with impact
117  // parameter larger than the maximum impact parameter, taking into
118  // account Coulomb distortion.
119  if(impactParameter>CoulombDistortion::maxImpactParameter(p->getSpecies(), kineticEnergy, theNucleus)) {
120  INCL_DEBUG("impactParameter>CoulombDistortion::maxImpactParameter" << std::endl);
121  delete p;
122  return -1.;
123  }
124 
125  ThreeVector position(impactParameter * std::cos(phi),
126  impactParameter * std::sin(phi),
127  0.);
128  p->setPosition(position);
129 
130  // Fill in the relevant kinematic variables
132  theNucleus->setIncomingMomentum(p->getMomentum());
133  theNucleus->setInitialEnergy(p->getEnergy()
134  + ParticleTable::getTableMass(theNucleus->getA(),theNucleus->getZ()));
135 
136  // Reset the particle kinematics to the INCL values
137  p->setINCLMass();
138  p->setEnergy(p->getMass() + kineticEnergy);
140 
143  firstAvatar = false;
144 
145  // Get the entry avatars from Coulomb and put them in the Store
146  ParticleEntryAvatar *theEntryAvatar = CoulombDistortion::bringToSurface(p, theNucleus);
147  if(theEntryAvatar) {
148  theNucleus->getStore()->addParticleEntryAvatar(theEntryAvatar);
149 
150  theNucleus->setProjectileChargeNumber(p->getZ());
151  theNucleus->setProjectileMassNumber(p->getA());
152 
153  return p->getTransversePosition().mag();
154  } else {
155  delete p;
156  return -1.;
157  }
158  }
159 
160  G4double StandardPropagationModel::shootComposite(ParticleSpecies const species, const G4double kineticEnergy, const G4double impactParameter, const G4double phi) {
161  theNucleus->setNucleusNucleusCollision();
162  currentTime = 0.0;
163 
164  // Create the ProjectileRemnant object
165  ProjectileRemnant *pr = new ProjectileRemnant(species, kineticEnergy);
166 
167  // Same stopping time as for nucleon-nucleus
168  maximumTime = 29.8 * std::pow(theNucleus->getA(), 0.16);
169 
170  // If the incoming cluster is slow, use a larger stopping time
171  const G4double rms = ParticleTable::getLargestNuclearRadius(pr->getA(), pr->getZ());
172  const G4double rMax = theNucleus->getUniverseRadius();
173  const G4double distance = 2.*rMax + 2.725*rms;
174  const G4double projectileVelocity = pr->boostVector().mag();
175  const G4double traversalTime = distance / projectileVelocity;
176  if(maximumTime < traversalTime)
177  maximumTime = traversalTime;
178  INCL_DEBUG("Cascade stopping time is " << maximumTime << std::endl);
179 
180  // If Coulomb is activated, do not process events with impact
181  // parameter larger than the maximum impact parameter, taking into
182  // account Coulomb distortion.
183  if(impactParameter>CoulombDistortion::maxImpactParameter(pr,theNucleus)) {
184  INCL_DEBUG("impactParameter>CoulombDistortion::maxImpactParameter" << std::endl);
185  delete pr;
186  return -1.;
187  }
188 
189  // Position the cluster at the right impact parameter
190  ThreeVector position(impactParameter * std::cos(phi),
191  impactParameter * std::sin(phi),
192  0.);
193  pr->setPosition(position);
194 
195  // Fill in the relevant kinematic variables
197  theNucleus->setIncomingMomentum(pr->getMomentum());
198  theNucleus->setInitialEnergy(pr->getEnergy()
199  + ParticleTable::getTableMass(theNucleus->getA(),theNucleus->getZ()));
200 
202  firstAvatar = false;
203 
204  // Get the entry avatars from Coulomb
205  IAvatarList theAvatarList
206  = CoulombDistortion::bringToSurface(pr, theNucleus);
207 
208  if(theAvatarList.empty()) {
209  INCL_DEBUG("No ParticleEntryAvatar found, transparent event" << std::endl);
210  delete pr;
211  return -1.;
212  }
213 
214  /* Store the internal kinematics of the projectile remnant.
215  *
216  * Note that this is at variance with the Fortran version, which stores
217  * the initial kinematics of the particles *after* putting the spectators
218  * on mass shell, but *before* removing the same energy from the
219  * participants. Due to the different code flow, doing so in the C++
220  * version leads to wrong excitation energies for the forced compound
221  * nucleus.
222  */
223  pr->storeComponents();
224 
225  // Tell the Nucleus about the ProjectileRemnant
226  theNucleus->setProjectileRemnant(pr);
227 
228  // Set the number of projectile particles
229  theNucleus->setProjectileChargeNumber(pr->getZ());
230  theNucleus->setProjectileMassNumber(pr->getA());
231 
232  // Register the ParticleEntryAvatars
233  theNucleus->getStore()->addParticleEntryAvatars(theAvatarList);
234 
235  return pr->getTransversePosition().mag();
236  }
237 
239  return maximumTime;
240  }
241 
243 // assert(time>0.0);
244  maximumTime = time;
245  }
246 
248  return currentTime;
249  }
250 
252  {
253  theNucleus = nucleus;
254  }
255 
257  {
258  if(anAvatar) theNucleus->getStore()->add(anAvatar);
259  }
260 
262  // Is either particle a participant?
263  if(!p1->isParticipant() && !p2->isParticipant() && p1->getParticipantType()==p2->getParticipantType()) return NULL;
264 
265  // Is it a pi-resonance collision (we don't treat them)?
266  if((p1->isResonance() && p2->isPion()) || (p1->isPion() && p2->isResonance()))
267  return NULL;
268 
269  // Will the avatar take place between now and the end of the cascade?
270  G4double minDistOfApproachSquared = 0.0;
271  G4double t = getTime(p1, p2, &minDistOfApproachSquared);
272  if(t>maximumTime || t<currentTime) return NULL;
273 
274  // Local energy. Jump through some hoops to calculate the cross section
275  // at the collision point, and clean up after yourself afterwards.
276  G4bool hasLocalEnergy;
277  if(p1->isPion() || p2->isPion())
278  hasLocalEnergy = ((theLocalEnergyDeltaType == FirstCollisionLocalEnergy &&
279  theNucleus->getStore()->getBook().getAcceptedCollisions()==0) ||
280  theLocalEnergyDeltaType == AlwaysLocalEnergy);
281  else
282  hasLocalEnergy = ((theLocalEnergyType == FirstCollisionLocalEnergy &&
283  theNucleus->getStore()->getBook().getAcceptedCollisions()==0) ||
284  theLocalEnergyType == AlwaysLocalEnergy);
285  const G4bool p1HasLocalEnergy = (hasLocalEnergy && !p1->isPion());
286  const G4bool p2HasLocalEnergy = (hasLocalEnergy && !p2->isPion());
287 
288  if(p1HasLocalEnergy) {
289  backupParticle1 = *p1;
290  p1->propagate(t - currentTime);
291  if(p1->getPosition().mag() > theNucleus->getSurfaceRadius(p1)) {
292  *p1 = backupParticle1;
293  return NULL;
294  }
296  }
297  if(p2HasLocalEnergy) {
298  backupParticle2 = *p2;
299  p2->propagate(t - currentTime);
300  if(p2->getPosition().mag() > theNucleus->getSurfaceRadius(p2)) {
301  *p2 = backupParticle2;
302  if(p1HasLocalEnergy) {
303  *p1 = backupParticle1;
304  }
305  return NULL;
306  }
308  }
309 
310  // Compute the total cross section
311  const G4double totalCrossSection = CrossSections::total(p1, p2);
313 
314  // Restore particles to their state before the local-energy tweak
315  if(p1HasLocalEnergy) {
316  *p1 = backupParticle1;
317  }
318  if(p2HasLocalEnergy) {
319  *p2 = backupParticle2;
320  }
321 
322  // Is the CM energy > cutNN? (no cutNN on the first collision)
323  if(theNucleus->getStore()->getBook().getAcceptedCollisions()>0
324  && p1->isNucleon() && p2->isNucleon()
325  && squareTotalEnergyInCM < BinaryCollisionAvatar::getCutNNSquared()) return NULL;
326 
327  // Do the particles come close enough to each other?
328  if(Math::tenPi*minDistOfApproachSquared > totalCrossSection) return NULL;
329 
330  // Bomb out if the two collision partners are the same particle
331 // assert(p1->getID() != p2->getID());
332 
333  // Return a new avatar, then!
334  return new G4INCL::BinaryCollisionAvatar(t, totalCrossSection, theNucleus, p1, p2);
335  }
336 
338  Intersection theIntersection(
340  aParticle->getPosition(),
341  aParticle->getPropagationVelocity(),
342  theNucleus->getSurfaceRadius(aParticle)));
343  G4double time;
344  if(theIntersection.exists) {
345  time = currentTime + theIntersection.time;
346  } else {
347  INCL_ERROR("Imaginary reflection time for particle: " << std::endl
348  << aParticle->print());
349  time = 10000.0;
350  }
351  return time;
352  }
353 
355  G4INCL::Particle const * const particleB, G4double *minDistOfApproach) const
356  {
357  G4double time;
358  G4INCL::ThreeVector t13 = particleA->getPropagationVelocity();
359  t13 -= particleB->getPropagationVelocity();
360  G4INCL::ThreeVector distance = particleA->getPosition();
361  distance -= particleB->getPosition();
362  const G4double t7 = t13.dot(distance);
363  const G4double dt = t13.mag2();
364  if(dt <= 1.0e-10) {
365  (*minDistOfApproach) = 100000.0;
366  return currentTime + 100000.0;
367  } else {
368  time = -t7/dt;
369  }
370  (*minDistOfApproach) = distance.mag2() + time * t7;
371  return currentTime + time;
372  }
373 
374  void StandardPropagationModel::generateUpdatedCollisions(const ParticleList &updatedParticles, const ParticleList &particles) {
375 
376  // Loop over all the updated particles
377  for(ParticleIter updated=updatedParticles.begin(), e=updatedParticles.end(); updated!=e; ++updated)
378  {
379  // Loop over all the particles
380  for(ParticleIter particle=particles.begin(), end=particles.end(); particle!=end; ++particle)
381  {
382  /* Consider the generation of a collision avatar only if (*particle)
383  * is not one of the updated particles.
384  * The criterion makes sure that you don't generate avatars between
385  * updated particles. */
386  if((*particle)->isInList(updatedParticles)) continue;
387 
388  registerAvatar(generateBinaryCollisionAvatar(*particle,*updated));
389  }
390  }
391  }
392 
394  // Loop over all the particles
395  for(ParticleIter p1=particles.begin(), e=particles.end(); p1!=e; ++p1) {
396  // Loop over the rest of the particles
397  for(ParticleIter p2 = p1 + 1; p2 != particles.end(); ++p2) {
399  }
400  }
401  }
402 
404 
405  const G4bool haveExcept = !except.empty();
406 
407  // Loop over all the particles
408  for(ParticleIter p1=particles.begin(), e=particles.end(); p1!=e; ++p1)
409  {
410  // Loop over the rest of the particles
411  ParticleIter p2 = p1;
412  for(++p2; p2 != particles.end(); ++p2)
413  {
414  // Skip the collision if both particles must be excluded
415  if(haveExcept && (*p1)->isInList(except) && (*p2)->isInList(except)) continue;
416 
418  }
419  }
420 
421  }
422 
424 
425  for(ParticleIter iter=particles.begin(), e=particles.end(); iter!=e; ++iter) {
426  G4double time = this->getReflectionTime(*iter);
427  if(time <= maximumTime) registerAvatar(new SurfaceAvatar(*iter, time, theNucleus));
428  }
429  ParticleList const &p = theNucleus->getStore()->getParticles();
430  generateUpdatedCollisions(particles, p); // Predict collisions with spectators and participants
431  }
432 
434  ParticleList const &particles = theNucleus->getStore()->getParticles();
435 // assert(!particles.empty());
436  for(ParticleIter i=particles.begin(), e=particles.end(); i!=e; ++i) {
437  G4double time = this->getReflectionTime(*i);
438  if(time <= maximumTime) registerAvatar(new SurfaceAvatar(*i, time, theNucleus));
439  }
440  generateCollisions(particles);
441  generateDecays(particles);
442  }
443 
444 #ifdef INCL_REGENERATE_AVATARS
445  void StandardPropagationModel::generateAllAvatarsExceptUpdated() {
446  ParticleList const &particles = theNucleus->getStore()->getParticles();
447 // assert(!particles.empty());
448  for(ParticleIter i=particles.begin(), e=particles.end(); i!=e; ++i) {
449  G4double time = this->getReflectionTime(*i);
450  if(time <= maximumTime) registerAvatar(new SurfaceAvatar(*i, time, theNucleus));
451  }
452  ParticleList except = theNucleus->getUpdatedParticles();
453  generateCollisions(particles,except);
454  generateDecays(particles);
455  }
456 #endif
457 
459  for(ParticleIter i=particles.begin(), e=particles.end(); i!=e; ++i) {
460  if((*i)->isDelta()) {
462  G4double time = currentTime + decayTime;
463  if(time <= maximumTime) {
464  registerAvatar(new DecayAvatar((*i), time, theNucleus));
465  }
466  }
467  }
468  }
469 
471  {
472  // We update only the information related to particles that were updated
473  // by the previous avatar.
474 #ifdef INCL_REGENERATE_AVATARS
475 #warning "The INCL_REGENERATE_AVATARS code has not been tested in a while. Use it at your peril."
476  if(theNucleus->getUpdatedParticles().size()!=0 || theNucleus->getCreatedParticles().size()!=0) {
477  // Regenerates the entire avatar list, skipping collisions between
478  // updated particles
479  theNucleus->getStore()->clearAvatars();
481  generateAllAvatarsExceptUpdated();
482  }
483 #else
484  // Deltas are created by transforming nucleon into a delta for
485  // efficiency reasons
486  Particle * const blockedDelta = theNucleus->getBlockedDelta();
487  ParticleList updatedParticles = theNucleus->getUpdatedParticles();
488  if(blockedDelta)
489  updatedParticles.push_back(blockedDelta);
490  generateDecays(updatedParticles);
491 
492  ParticleList needNewAvatars = theNucleus->getUpdatedParticles();
493  ParticleList const &created = theNucleus->getCreatedParticles();
494  needNewAvatars.insert(needNewAvatars.end(), created.begin(), created.end());
495  updateAvatars(needNewAvatars);
496 #endif
497 
498  G4INCL::IAvatar *theAvatar = theNucleus->getStore()->findSmallestTime();
499  if(theAvatar == 0) return 0; // Avatar list is empty
500  // theAvatar->dispose();
501 
502  if(theAvatar->getTime() < currentTime) {
503  INCL_ERROR("Avatar time = " << theAvatar->getTime() << ", currentTime = " << currentTime << std::endl);
504  return 0;
505  } else if(theAvatar->getTime() > currentTime) {
506  theNucleus->getStore()->timeStep(theAvatar->getTime() - currentTime);
507 
508  currentTime = theAvatar->getTime();
509  theNucleus->getStore()->getBook().setCurrentTime(currentTime);
510  }
511 
512  return theAvatar;
513  }
514 
515  void StandardPropagationModel::putSpectatorsOnShell(IAvatarList const &entryAvatars, ParticleList const &spectators) {
516  G4double deltaE = 0.0;
517  for(ParticleIter p=spectators.begin(), e=spectators.end(); p!=e; ++p) {
518  // put the spectators on shell (conserving their momentum)
519  const G4double oldEnergy = (*p)->getEnergy();
520  (*p)->setTableMass();
521  (*p)->adjustEnergyFromMomentum();
522  deltaE += (*p)->getEnergy() - oldEnergy;
523  }
524 
525  deltaE /= entryAvatars.size(); // energy to remove from each participant
526 
527  for(IAvatarIter a=entryAvatars.begin(), e=entryAvatars.end(); a!=e; ++a) {
528  // remove the energy from the participant
529  Particle *p = (*a)->getParticles().front();
530  ParticleType const t = p->getType();
531  // we also need to slightly correct the participant mass
532  const G4double energy = p->getEnergy() - deltaE
534  p->setEnergy(energy);
535  const G4double newMass = std::sqrt(energy*energy - p->getMomentum().mag2());
536  p->setMass(newMass);
537  }
538  }
539 }
G4int getA() const
Returns the baryon number.
void clearAvatars()
Definition: G4INCLStore.cc:252
void registerAvatar(G4INCL::IAvatar *anAvatar)
G4bool isResonance() const
Is it a resonance?
G4double shootComposite(ParticleSpecies const s, const G4double kineticEnergy, const G4double impactParameter, const G4double phi)
void generateAllAvatars()
(Re)Generate all possible avatars.
G4double getMass() const
Get the cached particle mass.
void setIncomingAngularMomentum(const ThreeVector &j)
Set the incoming angular-momentum vector.
void setIncomingMomentum(const ThreeVector &p)
Set the incoming momentum vector.
G4ThreadLocal ParticleMassFn getTableParticleMass
Static pointer to the mass function for particles.
G4double dot(const ThreeVector &v) const
ParticleList const & getParticles() const
Definition: G4INCLStore.hh:231
void transformToLocalEnergyFrame(Nucleus const *const n, Particle *const p)
G4int getAcceptedCollisions() const
Definition: G4INCLBook.hh:99
ThreeVector getPropagationVelocity() const
Get the propagation velocity of the particle.
void setInitialEnergy(const G4double e)
Set the initial energy.
G4double squareTotalEnergyInCM(Particle const *const p1, Particle const *const p2)
StandardPropagationModel(LocalEnergyType localEnergyType, LocalEnergyType localEnergyDeltaType)
void generateCollisions(const ParticleList &particles)
Generate and register collisions among particles in a list, except between those in another list...
const char * p
Definition: xmltok.h:285
const G4double tenPi
UnorderedVector< IAvatar * >::const_iterator IAvatarIter
#define INCL_ERROR(x)
IAvatar * findSmallestTime()
Definition: G4INCLStore.cc:173
const G4INCL::ThreeVector & getMomentum() const
Store * getStore() const
void setParticleNucleusCollision()
Set a particle-nucleus collision.
std::string print() const
void add(Particle *p)
Definition: G4INCLStore.cc:57
G4double getReflectionTime(G4INCL::Particle const *const aParticle)
Get the reflection time.
void updateAvatars(const ParticleList &particles)
void storeComponents()
Store the projectile components.
G4double getEnergy() const
void setINCLMass()
Set the mass of the Particle to its table mass.
virtual void makeProjectileSpectator()
void propagate(G4double step)
G4double getSurfaceRadius(Particle const *const particle) const
Get the maximum allowed radius for a given particle.
G4bool isParticipant() const
ThreeVector boostVector() const
static ParticleEntryAvatar * bringToSurface(Particle *p, Nucleus *const n)
Modify the momentum of an incoming particle and position it on the surface of the nucleus...
G4double mag2() const
G4double getTime(G4INCL::Particle const *const particleA, G4INCL::Particle const *const particleB, G4double *minDistOfApproach) const
void setEnergy(G4double energy)
Intersection getLaterTrajectoryIntersection(const ThreeVector &x0, const ThreeVector &p, const G4double r)
Compute the second intersection of a straight particle trajectory with a sphere.
G4INCL::ThreeVector getAngularMomentum() const
Get the total angular momentum (orbital + spin)
static G4double computeDecayTime(Particle *p)
void addParticleEntryAvatars(IAvatarList const &al)
Add one ParticleEntry avatar.
Definition: G4INCLStore.cc:78
G4double getTime() const
G4double shoot(ParticleSpecies const projectileSpecies, const G4double kineticEnergy, const G4double impactParameter, const G4double phi)
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
Static class for selecting Coulomb distortion.
UnorderedVector< Particle * > ParticleList
static G4double maxImpactParameter(ParticleSpecies const &p, const G4double kinE, Nucleus const *const n)
Return the maximum impact parameter for Coulomb-distorted trajectories.
Book & getBook()
Definition: G4INCLStore.hh:237
void initialiseParticleAvatarConnections()
Initialise the particleAvatarConnections map.
Definition: G4INCLStore.cc:267
ParticipantType getParticipantType() const
G4double getLargestNuclearRadius(const G4int A, const G4int Z)
bool G4bool
Definition: G4Types.hh:79
virtual G4INCL::ThreeVector getAngularMomentum() const
virtual void setPosition(const G4INCL::ThreeVector &position)
G4int getZ() const
Returns the charge number.
ParticleList const & getUpdatedParticles() const
ParticleList const & getCreatedParticles() const
ThreeVector getTransversePosition() const
Transverse component of the position w.r.t. the momentum.
IAvatar * generateBinaryCollisionAvatar(Particle *const p1, Particle *const p2)
Generate a two-particle avatar.
void generateDecays(const ParticleList &particles)
Generate decays for particles that can decay.
void addParticleEntryAvatar(IAvatar *a)
Add one ParticleEntry avatar.
Definition: G4INCLStore.cc:66
void setProjectileChargeNumber(G4int n)
Set the charge number of the projectile.
G4double total(Particle const *const p1, Particle const *const p2)
G4double getINCLMass(const G4int A, const G4int Z)
Get INCL nuclear mass (in MeV/c^2)
G4ThreadLocal NuclearMassFn getTableMass
Static pointer to the mass function for nuclei.
void setProjectileRemnant(ProjectileRemnant *const c)
Set the projectile remnant.
int position
Definition: filter.cc:7
const G4INCL::ThreeVector & getPosition() const
G4bool isNucleon() const
void setProjectileMassNumber(G4int n)
Set the mass number of the projectile.
void setCurrentTime(G4double t)
Definition: G4INCLBook.hh:96
virtual G4INCL::ParticleSpecies getSpecies() const
Get the particle species.
G4double getUniverseRadius() const
Getter for theUniverseRadius.
Intersection-point structure.
G4double mag() const
double G4double
Definition: G4Types.hh:76
#define INCL_DEBUG(x)
G4double shootParticle(ParticleType const t, const G4double kineticEnergy, const G4double impactParameter, const G4double phi)
void generateUpdatedCollisions(const ParticleList &updatedParticles, const ParticleList &particles)
Generate and register collisions between a list of updated particles and all the other particles...
Particle * getBlockedDelta() const
Get the delta that could not decay.
G4bool isPion() const
Is this a pion?
const ThreeVector & adjustMomentumFromEnergy()
Rescale the momentum to match the total energy.
void setNucleusNucleusCollision()
Set a nucleus-nucleus collision.
void timeStep(G4double step)
Definition: G4INCLStore.cc:224
void setPosition(const ThreeVector &position)
Set the position of the cluster.
ParticleList::const_iterator ParticleIter
Simple class for computing intersections between a straight line and a sphere.