Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4INCL::ProjectileRemnant Class Reference

#include <G4INCLProjectileRemnant.hh>

Inheritance diagram for G4INCL::ProjectileRemnant:
Collaboration diagram for G4INCL::ProjectileRemnant:

Public Types

typedef std::vector< G4doubleEnergyLevels
 
typedef std::map< long, G4doubleEnergyLevelMap
 

Public Member Functions

 ProjectileRemnant (ParticleSpecies const &species, const G4double kineticEnergy)
 
 ~ProjectileRemnant ()
 
void reset ()
 Reset the projectile remnant to the state at the beginning of the cascade. More...
 
void removeParticle (Particle *const p, const G4double theProjectileCorrection)
 Remove a nucleon from the projectile remnant. More...
 
ParticleList addDynamicalSpectators (ParticleList pL)
 Add back dynamical spectators to the projectile remnant. More...
 
ParticleList addMostDynamicalSpectators (ParticleList pL)
 Add back dynamical spectators to the projectile remnant. More...
 
ParticleList addAllDynamicalSpectators (ParticleList const &pL)
 Add back all dynamical spectators to the projectile remnant. More...
 
void deleteStoredComponents ()
 Clear the stored projectile components and delete the particles. More...
 
void clearStoredComponents ()
 Clear the stored projectile components. More...
 
void clearEnergyLevels ()
 Clear the stored energy levels. More...
 
G4double computeExcitationEnergyExcept (const long exceptID) const
 Compute the excitation energy when a nucleon is removed. More...
 
G4double computeExcitationEnergyWith (const ParticleList &pL) const
 Compute the excitation energy if some nucleons are put back. More...
 
void storeComponents ()
 Store the projectile components. More...
 
G4int getNumberStoredComponents () const
 Get the number of the stored components. More...
 
void storeEnergyLevels ()
 Store the energy levels. More...
 
EnergyLevels const & getGroundStateEnergies () const
 
- Public Member Functions inherited from G4INCL::Cluster
 Cluster (const G4int Z, const G4int A, const G4bool createParticleSampler=true)
 Standard Cluster constructor. More...
 
template<class Iterator >
 Cluster (Iterator begin, Iterator end)
 
virtual ~Cluster ()
 
 Cluster (const Cluster &rhs)
 Copy constructor. More...
 
Clusteroperator= (const Cluster &rhs)
 Assignment operator. More...
 
void swap (Cluster &rhs)
 Helper method for the assignment operator. More...
 
ParticleSpecies getSpecies () const
 Get the particle species. More...
 
void deleteParticles ()
 
void clearParticles ()
 
void setZ (const G4int Z)
 Set the charge number of the cluster. More...
 
void setA (const G4int A)
 Set the mass number of the cluster. More...
 
G4double getExcitationEnergy () const
 Get the excitation energy of the cluster. More...
 
void setExcitationEnergy (const G4double e)
 Set the excitation energy of the cluster. More...
 
virtual G4double getTableMass () const
 Get the real particle mass. More...
 
ParticleList const & getParticles () const
 
void removeParticle (Particle *const p)
 Remove a particle from the cluster components. More...
 
void addParticle (Particle *const p)
 
void updateClusterParameters ()
 Set total cluster mass, energy, size, etc. from the particles. More...
 
void addParticles (ParticleList const &pL)
 Add a list of particles to the cluster. More...
 
ParticleList getParticleList () const
 Returns the list of particles that make up the cluster. More...
 
std::string print () const
 
virtual void initializeParticles ()
 Initialise the NuclearDensity pointer and sample the particles. More...
 
void internalBoostToCM ()
 Boost to the CM of the component particles. More...
 
void putParticlesOffShell ()
 Put the cluster components off shell. More...
 
void setPosition (const ThreeVector &position)
 Set the position of the cluster. More...
 
void boost (const ThreeVector &aBoostVector)
 Boost the cluster with the indicated velocity. More...
 
void freezeInternalMotion ()
 Freeze the internal motion of the particles. More...
 
virtual void rotatePosition (const G4double angle, const ThreeVector &axis)
 Rotate position of all the particles. More...
 
virtual void rotateMomentum (const G4double angle, const ThreeVector &axis)
 Rotate momentum of all the particles. More...
 
virtual void makeProjectileSpectator ()
 Make all the components projectile spectators, too. More...
 
virtual void makeTargetSpectator ()
 Make all the components target spectators, too. More...
 
virtual void makeParticipant ()
 Make all the components participants, too. More...
 
ThreeVector const & getSpin () const
 Get the spin of the nucleus. More...
 
void setSpin (const ThreeVector &j)
 Set the spin of the nucleus. More...
 
G4INCL::ThreeVector getAngularMomentum () const
 Get the total angular momentum (orbital + spin) More...
 
- Public Member Functions inherited from G4INCL::Particle
 Particle ()
 
 Particle (ParticleType t, G4double energy, ThreeVector const &momentum, ThreeVector const &position)
 
 Particle (ParticleType t, ThreeVector const &momentum, ThreeVector const &position)
 
virtual ~Particle ()
 
 Particle (const Particle &rhs)
 Copy constructor. More...
 
Particleoperator= (const Particle &rhs)
 Assignment operator. More...
 
G4INCL::ParticleType getType () const
 
void setType (ParticleType t)
 
G4bool isNucleon () const
 
ParticipantType getParticipantType () const
 
void setParticipantType (ParticipantType const p)
 
G4bool isParticipant () const
 
G4bool isTargetSpectator () const
 
G4bool isProjectileSpectator () const
 
G4bool isPion () const
 Is this a pion? More...
 
G4bool isEta () const
 Is this a eta? More...
 
G4bool isOmega () const
 Is this a omega? More...
 
G4bool isEtaPrime () const
 Is this a etaprime? More...
 
G4bool isResonance () const
 Is it a resonance? More...
 
G4bool isDelta () const
 Is it a Delta? More...
 
G4int getA () const
 Returns the baryon number. More...
 
G4int getZ () const
 Returns the charge number. More...
 
G4double getBeta () const
 
ThreeVector boostVector () const
 
void boost (const ThreeVector &aBoostVector)
 
void lorentzContract (const ThreeVector &aBoostVector, const ThreeVector &refPos)
 Lorentz-contract the particle position around some center. More...
 
G4double getMass () const
 Get the cached particle mass. More...
 
G4double getINCLMass () const
 Get the INCL particle mass. More...
 
G4double getRealMass () const
 Get the real particle mass. More...
 
void setRealMass ()
 Set the mass of the Particle to its real mass. More...
 
void setTableMass ()
 Set the mass of the Particle to its table mass. More...
 
void setINCLMass ()
 Set the mass of the Particle to its table mass. More...
 
G4double getEmissionQValueCorrection (const G4int AParent, const G4int ZParent) const
 Computes correction on the emission Q-value. More...
 
G4double getTransferQValueCorrection (const G4int AFrom, const G4int ZFrom, const G4int ATo, const G4int ZTo) const
 Computes correction on the transfer Q-value. More...
 
G4double getInvariantMass () const
 Get the the particle invariant mass. More...
 
G4double getKineticEnergy () const
 Get the particle kinetic energy. More...
 
G4double getPotentialEnergy () const
 Get the particle potential energy. More...
 
void setPotentialEnergy (G4double v)
 Set the particle potential energy. More...
 
G4double getEnergy () const
 
void setMass (G4double mass)
 
void setEnergy (G4double energy)
 
const G4INCL::ThreeVectorgetMomentum () const
 
virtual void setMomentum (const G4INCL::ThreeVector &momentum)
 
const G4INCL::ThreeVectorgetPosition () const
 
G4double getHelicity ()
 
void setHelicity (G4double h)
 
void propagate (G4double step)
 
G4int getNumberOfCollisions () const
 Return the number of collisions undergone by the particle. More...
 
void setNumberOfCollisions (G4int n)
 Set the number of collisions undergone by the particle. More...
 
void incrementNumberOfCollisions ()
 Increment the number of collisions undergone by the particle. More...
 
G4int getNumberOfDecays () const
 Return the number of decays undergone by the particle. More...
 
void setNumberOfDecays (G4int n)
 Set the number of decays undergone by the particle. More...
 
void incrementNumberOfDecays ()
 Increment the number of decays undergone by the particle. More...
 
void setOutOfWell ()
 Mark the particle as out of its potential well. More...
 
G4bool isOutOfWell () const
 Check if the particle is out of its potential well. More...
 
void setEmissionTime (G4double t)
 
G4double getEmissionTime ()
 
ThreeVector getTransversePosition () const
 Transverse component of the position w.r.t. the momentum. More...
 
ThreeVector getLongitudinalPosition () const
 Longitudinal component of the position w.r.t. the momentum. More...
 
const ThreeVectoradjustMomentumFromEnergy ()
 Rescale the momentum to match the total energy. More...
 
G4double adjustEnergyFromMomentum ()
 Recompute the energy to match the momentum. More...
 
G4bool isCluster () const
 
void setFrozenMomentum (const ThreeVector &momentum)
 Set the frozen particle momentum. More...
 
void setFrozenEnergy (const G4double energy)
 Set the frozen particle momentum. More...
 
ThreeVector getFrozenMomentum () const
 Get the frozen particle momentum. More...
 
G4double getFrozenEnergy () const
 Get the frozen particle momentum. More...
 
ThreeVector getPropagationVelocity () const
 Get the propagation velocity of the particle. More...
 
void freezePropagation ()
 Freeze particle propagation. More...
 
void thawPropagation ()
 Unfreeze particle propagation. More...
 
virtual void rotatePositionAndMomentum (const G4double angle, const ThreeVector &axis)
 Rotate the particle position and momentum. More...
 
std::string print () const
 
std::string dump () const
 
long getID () const
 
ParticleList const * getParticles () const
 
G4double getReflectionMomentum () const
 Return the reflection momentum. More...
 
void setUncorrelatedMomentum (const G4double p)
 Set the uncorrelated momentum. More...
 
void rpCorrelate ()
 Make the particle follow a strict r-p correlation. More...
 
void rpDecorrelate ()
 Make the particle not follow a strict r-p correlation. More...
 
G4double getCosRPAngle () const
 Get the cosine of the angle between position and momentum. More...
 

Additional Inherited Members

- Protected Member Functions inherited from G4INCL::Particle
void swap (Particle &rhs)
 Helper method for the assignment operator. More...
 
- Protected Attributes inherited from G4INCL::Cluster
ParticleList particles
 
G4double theExcitationEnergy
 
ThreeVector theSpin
 
ParticleSamplertheParticleSampler
 
- Protected Attributes inherited from G4INCL::Particle
G4int theZ
 
G4int theA
 
ParticipantType theParticipantType
 
G4INCL::ParticleType theType
 
G4double theEnergy
 
G4doublethePropagationEnergy
 
G4double theFrozenEnergy
 
G4INCL::ThreeVector theMomentum
 
G4INCL::ThreeVectorthePropagationMomentum
 
G4INCL::ThreeVector theFrozenMomentum
 
G4INCL::ThreeVector thePosition
 
G4int nCollisions
 
G4int nDecays
 
G4double thePotentialEnergy
 
long ID
 
G4bool rpCorrelated
 
G4double uncorrelatedMomentum
 

Detailed Description

Definition at line 58 of file G4INCLProjectileRemnant.hh.

Member Typedef Documentation

Definition at line 63 of file G4INCLProjectileRemnant.hh.

Definition at line 62 of file G4INCLProjectileRemnant.hh.

Constructor & Destructor Documentation

G4INCL::ProjectileRemnant::ProjectileRemnant ( ParticleSpecies const &  species,
const G4double  kineticEnergy 
)
inline

Definition at line 65 of file G4INCLProjectileRemnant.hh.

66  : Cluster(species.theZ, species.theA) {
67 
68  // Use the table mass
69  setTableMass();
70 
71  // Set the kinematics
72  const G4double projectileMass = getMass();
73  const G4double energy = kineticEnergy + projectileMass;
74  const G4double momentumZ = std::sqrt(energy*energy - projectileMass*projectileMass);
75 
76  // Initialise the particles
80 
81  // Store the energy levels of the ProjectileRemnant (used to compute its
82  // excitation energy)
84 
85  // Boost the whole thing
86  const ThreeVector aBoostVector = ThreeVector(0.0, 0.0, momentumZ / energy);
87  boost(-aBoostVector);
88 
89  // Freeze the internal motion of the particles
91 
92  // Set as projectile spectator
94  }
G4double getMass() const
Get the cached particle mass.
void storeEnergyLevels()
Store the energy levels.
void putParticlesOffShell()
Put the cluster components off shell.
void boost(const ThreeVector &aBoostVector)
Boost the cluster with the indicated velocity.
virtual void initializeParticles()
Initialise the NuclearDensity pointer and sample the particles.
void freezeInternalMotion()
Freeze the internal motion of the particles.
G4double energy(const ThreeVector &p, const G4double m)
void setTableMass()
Set the mass of the Particle to its table mass.
void internalBoostToCM()
Boost to the CM of the component particles.
virtual void makeProjectileSpectator()
Make all the components projectile spectators, too.
double G4double
Definition: G4Types.hh:76
Cluster(const G4int Z, const G4int A, const G4bool createParticleSampler=true)
Standard Cluster constructor.

Here is the call graph for this function:

G4INCL::ProjectileRemnant::~ProjectileRemnant ( )
inline

Definition at line 96 of file G4INCLProjectileRemnant.hh.

96  {
98  // The ProjectileRemnant owns its particles
101  }
void deleteStoredComponents()
Clear the stored projectile components and delete the particles.
void clearEnergyLevels()
Clear the stored energy levels.

Here is the call graph for this function:

Member Function Documentation

ParticleList G4INCL::ProjectileRemnant::addAllDynamicalSpectators ( ParticleList const &  pL)

Add back all dynamical spectators to the projectile remnant.

Return a list of rejected dynamical spectators.

Definition at line 144 of file G4INCLProjectileRemnant.cc.

144  {
145  // Put all the spectators in the projectile
146  ThreeVector theNewMomentum = theMomentum;
147  G4double theNewEnergy = theEnergy;
148  G4int theNewA = theA;
149  G4int theNewZ = theZ;
150  for(ParticleIter p=pL.begin(), e=pL.end(); p!=e; ++p) {
151 // assert((*p)->isNucleon());
152  // Add the initial (off-shell) momentum and energy to the projectile remnant
153  theNewMomentum += getStoredMomentum(*p);
154  theNewEnergy += (*p)->getEnergy();
155  theNewA += (*p)->getA();
156  theNewZ += (*p)->getZ();
157  }
158 
159  // Check that the excitation energy of the new projectile remnant is non-negative
160  const G4double theNewMass = ParticleTable::getTableMass(theNewA,theNewZ);
161  const G4double theNewExcitationEnergy = computeExcitationEnergyWith(pL);
162  const G4double theNewEffectiveMass = theNewMass + theNewExcitationEnergy;
163 
164  // If this condition is satisfied, there is no solution. Fall back on the
165  // "most" method
166  if(theNewEnergy<theNewEffectiveMass) {
167  INCL_WARN("Could not add all the dynamical spectators back into the projectile remnant."
168  << " Falling back to the \"most\" method." << '\n');
169  return addMostDynamicalSpectators(pL);
170  }
171 
172  // Add all the participants to the projectile remnant
173  for(ParticleIter p=pL.begin(), e=pL.end(); p!=e; ++p) {
174  particles.push_back(*p);
175  }
176 
177  // Rescale the momentum of the projectile remnant so that sqrt(s) has the
178  // correct value
179  const G4double scalingFactorSquared = (theNewEnergy*theNewEnergy-theNewEffectiveMass*theNewEffectiveMass)/theNewMomentum.mag2();
180  const G4double scalingFactor = std::sqrt(scalingFactorSquared);
181  INCL_DEBUG("Scaling factor for the projectile-remnant momentum = " << scalingFactor << '\n');
182 
183  theA = theNewA;
184  theZ = theNewZ;
185  theMomentum = theNewMomentum * scalingFactor;
186  theEnergy = theNewEnergy;
187 
188  return ParticleList();
189  }
const char * p
Definition: xmltok.h:285
#define INCL_WARN(x)
int G4int
Definition: G4Types.hh:78
G4double computeExcitationEnergyWith(const ParticleList &pL) const
Compute the excitation energy if some nucleons are put back.
ParticleList particles
G4ThreadLocal NuclearMassFn getTableMass
Static pointer to the mass function for nuclei.
ParticleList addMostDynamicalSpectators(ParticleList pL)
Add back dynamical spectators to the projectile remnant.
G4INCL::ThreeVector theMomentum
double G4double
Definition: G4Types.hh:76
#define INCL_DEBUG(x)
ParticleList::const_iterator ParticleIter

Here is the call graph for this function:

ParticleList G4INCL::ProjectileRemnant::addDynamicalSpectators ( ParticleList  pL)

Add back dynamical spectators to the projectile remnant.

Try to add the dynamical spectators back to the projectile remnant. Refuse to do so if this leads to a negative projectile excitation energy.

Return a list of rejected dynamical spectators.

Definition at line 121 of file G4INCLProjectileRemnant.cc.

121  {
122  // Try as hard as possible to add back all the dynamical spectators.
123  // Don't add spectators that lead to negative excitation energies, but
124  // iterate over the spectators as many times as possible, until
125  // absolutely sure that all of them were rejected.
126  unsigned int accepted;
127  unsigned long loopCounter = 0;
128  const unsigned long maxLoopCounter = 10000000;
129  do {
130  accepted = 0;
131  ParticleList toBeAdded = pL;
132  for(ParticleIter p=toBeAdded.begin(), e=toBeAdded.end(); p!=e; ++p) {
133  G4bool isAccepted = addDynamicalSpectator(*p);
134  if(isAccepted) {
135  pL.remove(*p);
136  accepted++;
137  }
138  }
139  ++loopCounter;
140  } while(loopCounter<maxLoopCounter && accepted > 0); /* Loop checking, 10.07.2015, D.Mancusi */
141  return pL;
142  }
const char * p
Definition: xmltok.h:285
bool G4bool
Definition: G4Types.hh:79
ParticleList::const_iterator ParticleIter

Here is the call graph for this function:

ParticleList G4INCL::ProjectileRemnant::addMostDynamicalSpectators ( ParticleList  pL)

Add back dynamical spectators to the projectile remnant.

Try as hard as possible to add back all the dynamical spectators. Don't add spectators that lead to negative excitation energies. Start by adding all of them, and repeatedly remove the most troublesome one until the excitation energy becomes non-negative.

Return a list of rejected dynamical spectators.

Definition at line 191 of file G4INCLProjectileRemnant.cc.

191  {
192  // Try as hard as possible to add back all the dynamical spectators.
193  // Don't add spectators that lead to negative excitation energies. Start by
194  // adding all of them, and repeatedly remove the most troublesome one until
195  // the excitation energy becomes non-negative.
196 
197  // Put all the spectators in the projectile
198  ThreeVector theNewMomentum = theMomentum;
199  G4double theNewEnergy = theEnergy;
200  G4int theNewA = theA;
201  G4int theNewZ = theZ;
202  for(ParticleIter p=pL.begin(), e=pL.end(); p!=e; ++p) {
203 // assert((*p)->isNucleon());
204  // Add the initial (off-shell) momentum and energy to the projectile remnant
205  theNewMomentum += getStoredMomentum(*p);
206  theNewEnergy += (*p)->getEnergy();
207  theNewA += (*p)->getA();
208  theNewZ += (*p)->getZ();
209  }
210 
211  // Check that the excitation energy of the new projectile remnant is non-negative
212  const G4double theNewMass = ParticleTable::getTableMass(theNewA,theNewZ);
213  const G4double theNewInvariantMassSquared = theNewEnergy*theNewEnergy-theNewMomentum.mag2();
214 
215  G4bool positiveExcitationEnergy = false;
216  if(theNewInvariantMassSquared>=0.) {
217  const G4double theNewInvariantMass = std::sqrt(theNewInvariantMassSquared);
218  positiveExcitationEnergy = (theNewInvariantMass-theNewMass>-1.e-5);
219  }
220 
221  // Keep removing nucleons from the projectile remnant until we achieve a
222  // non-negative excitation energy.
223  ParticleList rejected;
224  while(!positiveExcitationEnergy && !pL.empty()) { /* Loop checking, 10.07.2015, D.Mancusi */
225  G4double maxExcitationEnergy = -1.E30;
226  ParticleMutableIter best = pL.end();
227  ThreeVector bestMomentum;
228  G4double bestEnergy = -1.;
229  G4int bestA = -1, bestZ = -1;
230  for(ParticleList::iterator p=pL.begin(), e=pL.end(); p!=e; ++p) {
231  // Subtract the initial (off-shell) momentum and energy from the new
232  // projectile remnant
233  const ThreeVector theNewerMomentum = theNewMomentum - getStoredMomentum(*p);
234  const G4double theNewerEnergy = theNewEnergy - (*p)->getEnergy();
235  const G4int theNewerA = theNewA - (*p)->getA();
236  const G4int theNewerZ = theNewZ - (*p)->getZ();
237 
238  const G4double theNewerMass = ParticleTable::getTableMass(theNewerA,theNewerZ);
239  const G4double theNewerInvariantMassSquared = theNewerEnergy*theNewerEnergy-theNewerMomentum.mag2();
240 
241  if(theNewerInvariantMassSquared>=-1.e-5) {
242  const G4double theNewerInvariantMass = std::sqrt(std::max(0.,theNewerInvariantMassSquared));
243  const G4double theNewerExcitationEnergy = ((theNewerA>1) ? theNewerInvariantMass-theNewerMass : 0.);
244  // Pick the nucleon that maximises the excitation energy of the
245  // ProjectileRemnant
246  if(theNewerExcitationEnergy>maxExcitationEnergy) {
247  best = p;
248  maxExcitationEnergy = theNewerExcitationEnergy;
249  bestMomentum = theNewerMomentum;
250  bestEnergy = theNewerEnergy;
251  bestA = theNewerA;
252  bestZ = theNewerZ;
253  }
254  }
255  }
256 
257  // If we couldn't even calculate the excitation energy, fail miserably
258  if(best==pL.end())
259  return pL;
260 
261  rejected.push_back(*best);
262  pL.erase(best);
263  theNewMomentum = bestMomentum;
264  theNewEnergy = bestEnergy;
265  theNewA = bestA;
266  theNewZ = bestZ;
267 
268  if(maxExcitationEnergy>0.) {
269  // Stop here
270  positiveExcitationEnergy = true;
271  }
272  }
273 
274  // Add the accepted participants to the projectile remnant
275  for(ParticleIter p=pL.begin(), e=pL.end(); p!=e; ++p) {
276  particles.push_back(*p);
277  }
278  theA = theNewA;
279  theZ = theNewZ;
280  theMomentum = theNewMomentum;
281  theEnergy = theNewEnergy;
282 
283  return rejected;
284  }
const char * p
Definition: xmltok.h:285
int G4int
Definition: G4Types.hh:78
bool G4bool
Definition: G4Types.hh:79
ParticleList particles
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4ThreadLocal NuclearMassFn getTableMass
Static pointer to the mass function for nuclei.
G4INCL::ThreeVector theMomentum
double G4double
Definition: G4Types.hh:76
ParticleList::iterator ParticleMutableIter
ParticleList::const_iterator ParticleIter

Here is the call graph for this function:

Here is the caller graph for this function:

void G4INCL::ProjectileRemnant::clearEnergyLevels ( )
inline

Clear the stored energy levels.

Definition at line 153 of file G4INCLProjectileRemnant.hh.

153  {
154  theInitialEnergyLevels.clear();
155  theGroundStateEnergies.clear();
156  }

Here is the caller graph for this function:

void G4INCL::ProjectileRemnant::clearStoredComponents ( )
inline

Clear the stored projectile components.

Definition at line 148 of file G4INCLProjectileRemnant.hh.

148  {
149  storedComponents.clear();
150  }

Here is the caller graph for this function:

G4double G4INCL::ProjectileRemnant::computeExcitationEnergyExcept ( const long  exceptID) const

Compute the excitation energy when a nucleon is removed.

Compute the excitation energy of the projectile-like remnant as the difference between the initial and the present configuration. This follows the algorithm proposed by A. Boudard in INCL4.2-HI, as implemented in Geant4.

Returns
the excitation energy

Definition at line 316 of file G4INCLProjectileRemnant.cc.

316  {
317  const EnergyLevels theEnergyLevels = getPresentEnergyLevelsExcept(exceptID);
318  return computeExcitationEnergy(theEnergyLevels);
319  }
std::vector< G4double > EnergyLevels

Here is the caller graph for this function:

G4double G4INCL::ProjectileRemnant::computeExcitationEnergyWith ( const ParticleList pL) const

Compute the excitation energy if some nucleons are put back.

Returns
the excitation energy

Definition at line 321 of file G4INCLProjectileRemnant.cc.

321  {
322  const EnergyLevels theEnergyLevels = getPresentEnergyLevelsWith(pL);
323  return computeExcitationEnergy(theEnergyLevels);
324  }
std::vector< G4double > EnergyLevels

Here is the caller graph for this function:

void G4INCL::ProjectileRemnant::deleteStoredComponents ( )
inline

Clear the stored projectile components and delete the particles.

Definition at line 141 of file G4INCLProjectileRemnant.hh.

141  {
142  for(std::map<long,Particle*>::const_iterator p=storedComponents.begin(), e=storedComponents.end(); p!=e; ++p)
143  delete p->second;
145  }
const char * p
Definition: xmltok.h:285
void clearStoredComponents()
Clear the stored projectile components.

Here is the call graph for this function:

Here is the caller graph for this function:

EnergyLevels const& G4INCL::ProjectileRemnant::getGroundStateEnergies ( ) const
inline

Definition at line 207 of file G4INCLProjectileRemnant.hh.

207  {
208  return theGroundStateEnergies;
209  }
G4int G4INCL::ProjectileRemnant::getNumberStoredComponents ( ) const
inline

Get the number of the stored components.

Definition at line 184 of file G4INCLProjectileRemnant.hh.

184  {
185  return storedComponents.size();
186  }
void G4INCL::ProjectileRemnant::removeParticle ( Particle *const  p,
const G4double  theProjectileCorrection 
)

Remove a nucleon from the projectile remnant.

Parameters
pparticle to be removed
theProjectileCorrectioncorrection to be given to the projectile total energy

Definition at line 76 of file G4INCLProjectileRemnant.cc.

76  {
77 // assert(p->isNucleon());
78 
79  INCL_DEBUG("The following Particle is about to be removed from the ProjectileRemnant:"
80  << '\n' << p->print()
81  << "theProjectileCorrection=" << theProjectileCorrection << '\n');
82  // Update A, Z, momentum and energy of the projectile remnant
83  theA -= p->getA();
84  theZ -= p->getZ();
85 
86  ThreeVector const &oldMomentum = p->getMomentum();
87  const G4double oldEnergy = p->getEnergy();
89 
90 #if !defined(NDEBUG) && !defined(INCLXX_IN_GEANT4_MODE)
91  ThreeVector theTotalMomentum;
92  G4double theTotalEnergy = 0.;
93  const G4double theThreshold = 0.1;
94 #endif
95 
96  if(getA()>0) { // if there are any particles left
97 // assert((unsigned int)getA()==particles.size());
98 
99  const G4double theProjectileCorrectionPerNucleon = theProjectileCorrection / particles.size();
100 
101  // Update the kinematics of the components
102  for(ParticleIter i=particles.begin(), e=particles.end(); i!=e; ++i) {
103  (*i)->setEnergy((*i)->getEnergy() + theProjectileCorrectionPerNucleon);
104  (*i)->setMass((*i)->getInvariantMass());
105 #if !defined(NDEBUG) && !defined(INCLXX_IN_GEANT4_MODE)
106  theTotalMomentum += (*i)->getMomentum();
107  theTotalEnergy += (*i)->getEnergy();
108 #endif
109  }
110  }
111 
112  theMomentum -= oldMomentum;
113  theEnergy -= oldEnergy - theProjectileCorrection;
114 
115 // assert(std::abs((theTotalMomentum-theMomentum).mag())<theThreshold);
116 // assert(std::abs(theTotalEnergy-theEnergy)<theThreshold);
117  INCL_DEBUG("After Particle removal, the ProjectileRemnant looks like this:"
118  << '\n' << print());
119  }
G4int getA() const
Returns the baryon number.
const char * p
Definition: xmltok.h:285
std::string print() const
ParticleList particles
void removeParticle(Particle *const p)
Remove a particle from the cluster components.
G4INCL::ThreeVector theMomentum
double G4double
Definition: G4Types.hh:76
#define INCL_DEBUG(x)
ParticleList::const_iterator ParticleIter

Here is the call graph for this function:

Here is the caller graph for this function:

void G4INCL::ProjectileRemnant::reset ( )

Reset the projectile remnant to the state at the beginning of the cascade.

Definition at line 51 of file G4INCLProjectileRemnant.cc.

51  {
53  thePosition = ThreeVector();
54  theMomentum = ThreeVector();
55  theEnergy = 0.0;
56  thePotentialEnergy = 0.0;
57  theA = 0;
58  theZ = 0;
59  nCollisions = 0;
60 
61  for(std::map<long, Particle*>::const_iterator i=storedComponents.begin(); i!=storedComponents.end(); ++i) {
62  Particle *p = new Particle(*(i->second));
63  EnergyLevelMap::iterator energyIter = theInitialEnergyLevels.find(i->first);
64 // assert(energyIter!=theInitialEnergyLevels.end());
65  const G4double energyLevel = energyIter->second;
66  theInitialEnergyLevels.erase(energyIter);
67  theInitialEnergyLevels[p->getID()] = energyLevel;
68  addParticle(p);
69  }
70  if(theA>0)
71  thePosition /= theA;
72  setTableMass();
73  INCL_DEBUG("ProjectileRemnant object was reset:" << '\n' << print());
74  }
const char * p
Definition: xmltok.h:285
std::string print() const
G4INCL::ThreeVector thePosition
void setTableMass()
Set the mass of the Particle to its table mass.
void addParticle(Particle *const p)
G4INCL::ThreeVector theMomentum
double G4double
Definition: G4Types.hh:76
#define INCL_DEBUG(x)
G4double thePotentialEnergy

Here is the call graph for this function:

void G4INCL::ProjectileRemnant::storeComponents ( )
inline

Store the projectile components.

Definition at line 176 of file G4INCLProjectileRemnant.hh.

176  {
177  for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
178  // Store the particles (needed for forced CN)
179  storedComponents[(*p)->getID()]=new Particle(**p);
180  }
181  }
const char * p
Definition: xmltok.h:285
ParticleList particles
ParticleList::const_iterator ParticleIter

Here is the call graph for this function:

Here is the caller graph for this function:

void G4INCL::ProjectileRemnant::storeEnergyLevels ( )
inline

Store the energy levels.

Definition at line 189 of file G4INCLProjectileRemnant.hh.

189  {
190  EnergyLevels energies;
191 
192  for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
193  const G4double theCMEnergy = (*p)->getEnergy();
194  // Store the CM energy in the EnergyLevels map
195  theInitialEnergyLevels[(*p)->getID()] = theCMEnergy;
196  energies.push_back(theCMEnergy);
197  }
198 
199  std::sort(energies.begin(), energies.end());
200 // assert(energies.size()==(unsigned int)theA);
201  theGroundStateEnergies.resize(energies.size());
202  // Compute the partial sums of the CM energies -- they are our reference
203  // ground-state energies for any number of nucleons
204  std::partial_sum(energies.begin(), energies.end(), theGroundStateEnergies.begin());
205  }
const char * p
Definition: xmltok.h:285
ParticleList particles
std::vector< G4double > EnergyLevels
double G4double
Definition: G4Types.hh:76
ParticleList::const_iterator ParticleIter

Here is the caller graph for this function:


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