Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
|
#include <G4INCLCluster.hh>
Public Member Functions | |
Cluster (const G4int Z, const G4int A, const G4int S, const G4bool createParticleSampler=true) | |
Standard Cluster constructor. | |
template<class Iterator > | |
Cluster (Iterator begin, Iterator end) | |
virtual | ~Cluster () |
Cluster (const Cluster &rhs) | |
Copy constructor. | |
Cluster & | operator= (const Cluster &rhs) |
Assignment operator. | |
void | swap (Cluster &rhs) |
Helper method for the assignment operator. | |
ParticleSpecies | getSpecies () const |
Get the particle species. | |
void | deleteParticles () |
void | clearParticles () |
void | setZ (const G4int Z) |
Set the charge number of the cluster. | |
void | setA (const G4int A) |
Set the mass number of the cluster. | |
void | setS (const G4int S) |
Set the strangess number of the cluster. | |
G4double | getExcitationEnergy () const |
Get the excitation energy of the cluster. | |
void | setExcitationEnergy (const G4double e) |
Set the excitation energy of the cluster. | |
virtual G4double | getTableMass () const |
Get the real particle mass. | |
ParticleList const & | getParticles () const |
void | removeParticle (Particle *const p) |
Remove a particle from the cluster components. | |
void | addParticle (Particle *const p) |
void | updateClusterParameters () |
Set total cluster mass, energy, size, etc. from the particles. | |
void | addParticles (ParticleList const &pL) |
Add a list of particles to the cluster. | |
ParticleList | getParticleList () const |
Returns the list of particles that make up the cluster. | |
std::string | print () const |
virtual void | initializeParticles () |
Initialise the NuclearDensity pointer and sample the particles. | |
void | internalBoostToCM () |
Boost to the CM of the component particles. | |
void | putParticlesOffShell () |
Put the cluster components off shell. | |
void | setPosition (const ThreeVector &position) |
Set the position of the cluster. | |
void | boost (const ThreeVector &aBoostVector) |
Boost the cluster with the indicated velocity. | |
void | freezeInternalMotion () |
Freeze the internal motion of the particles. | |
virtual void | rotatePosition (const G4double angle, const ThreeVector &axis) |
Rotate position of all the particles. | |
virtual void | rotateMomentum (const G4double angle, const ThreeVector &axis) |
Rotate momentum of all the particles. | |
virtual void | makeProjectileSpectator () |
Make all the components projectile spectators, too. | |
virtual void | makeTargetSpectator () |
Make all the components target spectators, too. | |
virtual void | makeParticipant () |
Make all the components participants, too. | |
ThreeVector const & | getSpin () const |
Get the spin of the nucleus. | |
void | setSpin (const ThreeVector &j) |
Set the spin of the nucleus. | |
G4INCL::ThreeVector | getAngularMomentum () const |
Get the total angular momentum (orbital + spin) | |
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. | |
Particle & | operator= (const Particle &rhs) |
Assignment operator. | |
G4INCL::ParticleType | getType () const |
virtual G4INCL::ParticleSpecies | getSpecies () const |
Get the particle species. | |
void | setType (ParticleType t) |
G4bool | isNucleon () const |
ParticipantType | getParticipantType () const |
void | setParticipantType (ParticipantType const p) |
G4bool | isParticipant () const |
G4bool | isTargetSpectator () const |
G4bool | isProjectileSpectator () const |
virtual void | makeParticipant () |
virtual void | makeTargetSpectator () |
virtual void | makeProjectileSpectator () |
G4bool | isPion () const |
Is this a pion? | |
G4bool | isEta () const |
Is this an eta? | |
G4bool | isOmega () const |
Is this an omega? | |
G4bool | isEtaPrime () const |
Is this an etaprime? | |
G4bool | isPhoton () const |
Is this a photon? | |
G4bool | isResonance () const |
Is it a resonance? | |
G4bool | isDelta () const |
Is it a Delta? | |
G4bool | isSigma () const |
Is this a Sigma? | |
G4bool | isKaon () const |
Is this a Kaon? | |
G4bool | isAntiKaon () const |
Is this an antiKaon? | |
G4bool | isLambda () const |
Is this a Lambda? | |
G4bool | isNucleonorLambda () const |
Is this a Nucleon or a Lambda? | |
G4bool | isHyperon () const |
Is this an Hyperon? | |
G4bool | isMeson () const |
Is this a Meson? | |
G4bool | isBaryon () const |
Is this a Baryon? | |
G4bool | isStrange () const |
Is this an Strange? | |
G4int | getA () const |
Returns the baryon number. | |
G4int | getZ () const |
Returns the charge number. | |
G4int | getS () const |
Returns the strangeness number. | |
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. | |
G4double | getMass () const |
Get the cached particle mass. | |
G4double | getINCLMass () const |
Get the INCL particle mass. | |
virtual G4double | getTableMass () const |
Get the tabulated particle mass. | |
G4double | getRealMass () const |
Get the real particle mass. | |
void | setRealMass () |
Set the mass of the Particle to its real mass. | |
void | setTableMass () |
Set the mass of the Particle to its table mass. | |
void | setINCLMass () |
Set the mass of the Particle to its table mass. | |
G4double | getEmissionQValueCorrection (const G4int AParent, const G4int ZParent) const |
Computes correction on the emission Q-value. | |
G4double | getTransferQValueCorrection (const G4int AFrom, const G4int ZFrom, const G4int ATo, const G4int ZTo) const |
Computes correction on the transfer Q-value. | |
G4double | getEmissionQValueCorrection (const G4int AParent, const G4int ZParent, const G4int SParent) const |
Computes correction on the emission Q-value for hypernuclei. | |
G4double | getTransferQValueCorrection (const G4int AFrom, const G4int ZFrom, const G4int SFrom, const G4int ATo, const G4int ZTo, const G4int STo) const |
Computes correction on the transfer Q-value for hypernuclei. | |
G4double | getInvariantMass () const |
Get the the particle invariant mass. | |
G4double | getKineticEnergy () const |
Get the particle kinetic energy. | |
G4double | getPotentialEnergy () const |
Get the particle potential energy. | |
void | setPotentialEnergy (G4double v) |
Set the particle potential energy. | |
G4double | getEnergy () const |
void | setMass (G4double mass) |
void | setEnergy (G4double energy) |
const G4INCL::ThreeVector & | getMomentum () const |
virtual G4INCL::ThreeVector | getAngularMomentum () const |
virtual void | setMomentum (const G4INCL::ThreeVector &momentum) |
const G4INCL::ThreeVector & | getPosition () const |
virtual void | setPosition (const G4INCL::ThreeVector &position) |
G4double | getHelicity () |
void | setHelicity (G4double h) |
void | propagate (G4double step) |
G4int | getNumberOfCollisions () const |
Return the number of collisions undergone by the particle. | |
void | setNumberOfCollisions (G4int n) |
Set the number of collisions undergone by the particle. | |
void | incrementNumberOfCollisions () |
Increment the number of collisions undergone by the particle. | |
G4int | getNumberOfDecays () const |
Return the number of decays undergone by the particle. | |
void | setNumberOfDecays (G4int n) |
Set the number of decays undergone by the particle. | |
void | incrementNumberOfDecays () |
Increment the number of decays undergone by the particle. | |
void | setOutOfWell () |
Mark the particle as out of its potential well. | |
G4bool | isOutOfWell () const |
Check if the particle is out of its potential well. | |
void | setEmissionTime (G4double t) |
G4double | getEmissionTime () |
ThreeVector | getTransversePosition () const |
Transverse component of the position w.r.t. the momentum. | |
ThreeVector | getLongitudinalPosition () const |
Longitudinal component of the position w.r.t. the momentum. | |
const ThreeVector & | adjustMomentumFromEnergy () |
Rescale the momentum to match the total energy. | |
G4double | adjustEnergyFromMomentum () |
Recompute the energy to match the momentum. | |
G4bool | isCluster () const |
void | setFrozenMomentum (const ThreeVector &momentum) |
Set the frozen particle momentum. | |
void | setFrozenEnergy (const G4double energy) |
Set the frozen particle momentum. | |
ThreeVector | getFrozenMomentum () const |
Get the frozen particle momentum. | |
G4double | getFrozenEnergy () const |
Get the frozen particle momentum. | |
ThreeVector | getPropagationVelocity () const |
Get the propagation velocity of the particle. | |
void | freezePropagation () |
Freeze particle propagation. | |
void | thawPropagation () |
Unfreeze particle propagation. | |
virtual void | rotatePositionAndMomentum (const G4double angle, const ThreeVector &axis) |
Rotate the particle position and momentum. | |
virtual void | rotatePosition (const G4double angle, const ThreeVector &axis) |
Rotate the particle position. | |
virtual void | rotateMomentum (const G4double angle, const ThreeVector &axis) |
Rotate the particle momentum. | |
std::string | print () const |
std::string | dump () const |
long | getID () const |
ParticleList const * | getParticles () const |
G4double | getReflectionMomentum () const |
Return the reflection momentum. | |
void | setUncorrelatedMomentum (const G4double p) |
Set the uncorrelated momentum. | |
void | rpCorrelate () |
Make the particle follow a strict r-p correlation. | |
void | rpDecorrelate () |
Make the particle not follow a strict r-p correlation. | |
G4double | getCosRPAngle () const |
Get the cosine of the angle between position and momentum. | |
G4double | getParticleBias () const |
Get the particle bias. | |
void | setParticleBias (G4double ParticleBias) |
Set the particle bias. | |
std::vector< G4int > | getBiasCollisionVector () const |
Get the vector list of biased vertices on the particle path. | |
void | setBiasCollisionVector (std::vector< G4int > BiasCollisionVector) |
Set the vector list of biased vertices on the particle path. | |
G4int | getNumberOfKaon () const |
Number of Kaon inside de nucleus. | |
void | setNumberOfKaon (const G4int NK) |
Protected Attributes | |
ParticleList | particles |
G4double | theExcitationEnergy |
ThreeVector | theSpin |
ParticleSampler * | theParticleSampler |
Protected Attributes inherited from G4INCL::Particle | |
G4int | theZ |
G4int | theA |
G4int | theS |
ParticipantType | theParticipantType |
G4INCL::ParticleType | theType |
G4double | theEnergy |
G4double * | thePropagationEnergy |
G4double | theFrozenEnergy |
G4INCL::ThreeVector | theMomentum |
G4INCL::ThreeVector * | thePropagationMomentum |
G4INCL::ThreeVector | theFrozenMomentum |
G4INCL::ThreeVector | thePosition |
G4int | nCollisions |
G4int | nDecays |
G4double | thePotentialEnergy |
long | ID |
G4bool | rpCorrelated |
G4double | uncorrelatedMomentum |
G4double | theParticleBias |
G4int | theNKaon |
The number of Kaons inside the nucleus (update during the cascade) | |
Additional Inherited Members | |
Static Public Member Functions inherited from G4INCL::Particle | |
static G4double | getTotalBias () |
General bias vector function. | |
static void | setINCLBiasVector (std::vector< G4double > NewVector) |
static void | FillINCLBiasVector (G4double newBias) |
static G4double | getBiasFromVector (std::vector< G4int > VectorBias) |
static std::vector< G4int > | MergeVectorBias (Particle const *const p1, Particle const *const p2) |
static std::vector< G4int > | MergeVectorBias (std::vector< G4int > p1, Particle const *const p2) |
Static Public Attributes inherited from G4INCL::Particle | |
static std::vector< G4double > | INCLBiasVector |
Time ordered vector of all bias applied. | |
static G4ThreadLocal G4int | nextBiasedCollisionID = 0 |
Protected Member Functions inherited from G4INCL::Particle | |
void | swap (Particle &rhs) |
Helper method for the assignment operator. | |
Cluster is a particle (inherits from the Particle class) that is actually a collection of elementary particles.
Definition at line 52 of file G4INCLCluster.hh.
|
inline |
Standard Cluster constructor.
This constructor should mainly be used when constructing Nucleus or when constructing Clusters to be used as composite projectiles.
Definition at line 60 of file G4INCLCluster.hh.
|
inline |
A cluster can be directly built from a list of particles.
Definition at line 79 of file G4INCLCluster.hh.
|
inlinevirtual |
Definition at line 94 of file G4INCLCluster.hh.
|
inline |
Copy constructor.
Definition at line 99 of file G4INCLCluster.hh.
|
inline |
Add one particle to the cluster. This updates the cluster mass, energy, size, etc.
Definition at line 178 of file G4INCLCluster.hh.
Referenced by Cluster(), and G4INCL::ProjectileRemnant::reset().
|
inline |
Add a list of particles to the cluster.
Definition at line 213 of file G4INCLCluster.hh.
|
inline |
Boost the cluster with the indicated velocity.
The Cluster is boosted as a whole, just like any Particle object; moreover, the internal components (particles list) are also boosted, according to Alain Boudard's off-shell recipe.
aBoostVector | the velocity to boost to [c] |
Definition at line 344 of file G4INCLCluster.hh.
Referenced by G4INCL::ProjectileRemnant::ProjectileRemnant().
|
inline |
Definition at line 143 of file G4INCLCluster.hh.
Referenced by deleteParticles().
|
inline |
Definition at line 136 of file G4INCLCluster.hh.
Referenced by G4INCL::Store::clearOutgoing(), G4INCL::Nucleus::decayOutgoingClusters(), G4INCL::ProjectileRemnant::reset(), and G4INCL::ProjectileRemnant::~ProjectileRemnant().
|
inline |
Freeze the internal motion of the particles.
Each particle is assigned a frozen momentum four-vector determined by the collective cluster velocity. This is used for propagation, but not for dynamics. Normal propagation is restored by calling the Particle::thawPropagation() method, which should be done in InteractionAvatar::postInteraction.
Definition at line 367 of file G4INCLCluster.hh.
Referenced by G4INCL::ProjectileRemnant::ProjectileRemnant().
|
inlinevirtual |
Get the total angular momentum (orbital + spin)
Reimplemented from G4INCL::Particle.
Definition at line 428 of file G4INCLCluster.hh.
Referenced by G4INCL::Nucleus::computeRecoilKinematics(), and G4INCL::StandardPropagationModel::shootComposite().
|
inline |
Get the excitation energy of the cluster.
Definition at line 155 of file G4INCLCluster.hh.
Referenced by G4INCL::Nucleus::fillEventInfo(), and G4INCL::Nucleus::getConservationBalance().
|
inline |
Returns the list of particles that make up the cluster.
Definition at line 219 of file G4INCLCluster.hh.
|
inline |
Get the list of particles in the cluster.
Definition at line 169 of file G4INCLCluster.hh.
Referenced by G4INCL::Nucleus::applyFinalState(), G4INCL::CoulombNone::bringToSurface(), and G4INCL::SurfaceAvatar::postInteraction().
|
inlinevirtual |
Get the particle species.
Reimplemented from G4INCL::Particle.
Definition at line 132 of file G4INCLCluster.hh.
|
inline |
Get the spin of the nucleus.
Definition at line 422 of file G4INCLCluster.hh.
Referenced by G4INCL::Nucleus::fillEventInfo(), and getAngularMomentum().
|
inlinevirtual |
Get the real particle mass.
Overloads the Particle method.
Reimplemented from G4INCL::Particle.
Definition at line 164 of file G4INCLCluster.hh.
Referenced by G4INCL::Nucleus::useFusionKinematics().
|
virtual |
Initialise the NuclearDensity pointer and sample the particles.
Reimplemented in G4INCL::Nucleus.
Definition at line 43 of file G4INCLCluster.cc.
Referenced by G4INCL::Nucleus::initializeParticles(), and G4INCL::ProjectileRemnant::ProjectileRemnant().
|
inline |
Boost to the CM of the component particles.
The position of all particles in the particles list is shifted so that their centre of mass is in the origin and their total momentum is zero.
Definition at line 254 of file G4INCLCluster.hh.
Referenced by G4INCL::ProjectileRemnant::ProjectileRemnant().
|
inlinevirtual |
Make all the components participants, too.
Reimplemented from G4INCL::Particle.
Definition at line 414 of file G4INCLCluster.hh.
|
inlinevirtual |
Make all the components projectile spectators, too.
Reimplemented from G4INCL::Particle.
Definition at line 398 of file G4INCLCluster.hh.
Referenced by G4INCL::ProjectileRemnant::ProjectileRemnant().
|
inlinevirtual |
Make all the components target spectators, too.
Reimplemented from G4INCL::Particle.
Definition at line 406 of file G4INCLCluster.hh.
Assignment operator.
Definition at line 114 of file G4INCLCluster.hh.
|
inline |
Definition at line 221 of file G4INCLCluster.hh.
Referenced by boost(), G4INCL::SurfaceAvatar::getChannel(), initializeParticles(), internalBoostToCM(), putParticlesOffShell(), G4INCL::ProjectileRemnant::removeParticle(), and G4INCL::ProjectileRemnant::reset().
|
inline |
Put the cluster components off shell.
The Cluster components are put off shell in such a way that their total energy equals the cluster mass.
Definition at line 306 of file G4INCLCluster.hh.
Referenced by G4INCL::ProjectileRemnant::ProjectileRemnant().
|
inline |
Remove a particle from the cluster components.
Definition at line 172 of file G4INCLCluster.hh.
Referenced by G4INCL::ProjectileRemnant::removeParticle().
|
virtual |
Rotate momentum of all the particles.
This includes the cluster components. Overloads Particle::rotateMomentum().
angle | the rotation angle |
axis | a unit vector representing the rotation axis |
Reimplemented from G4INCL::Particle.
Definition at line 64 of file G4INCLCluster.cc.
|
virtual |
Rotate position of all the particles.
This includes the cluster components. Overloads Particle::rotateMomentum().
angle | the rotation angle |
axis | a unit vector representing the rotation axis |
Reimplemented from G4INCL::Particle.
Definition at line 57 of file G4INCLCluster.cc.
|
inline |
Set the mass number of the cluster.
Definition at line 149 of file G4INCLCluster.hh.
|
inline |
Set the excitation energy of the cluster.
Definition at line 158 of file G4INCLCluster.hh.
Referenced by G4INCL::Nucleus::finalizeProjectileRemnant().
|
inlinevirtual |
Set the position of the cluster.
This overloads the Particle method to take into account that the positions of the cluster members must be updated as well.
Reimplemented from G4INCL::Particle.
Definition at line 328 of file G4INCLCluster.hh.
Referenced by G4INCL::StandardPropagationModel::shootComposite().
|
inline |
Set the strangess number of the cluster.
Definition at line 152 of file G4INCLCluster.hh.
|
inline |
Set the spin of the nucleus.
Definition at line 425 of file G4INCLCluster.hh.
Referenced by G4INCL::Nucleus::finalizeProjectileRemnant(), and G4INCL::Nucleus::useFusionKinematics().
|
inline |
Set the charge number of the cluster.
Definition at line 146 of file G4INCLCluster.hh.
|
inline |
Helper method for the assignment operator.
Definition at line 122 of file G4INCLCluster.hh.
Referenced by operator=().
|
inline |
Set total cluster mass, energy, size, etc. from the particles.
Definition at line 191 of file G4INCLCluster.hh.
Referenced by addParticles(), and initializeParticles().
|
protected |
Definition at line 451 of file G4INCLCluster.hh.
Referenced by G4INCL::ProjectileRemnant::addAllDynamicalSpectators(), G4INCL::ProjectileRemnant::addMostDynamicalSpectators(), addParticle(), addParticles(), boost(), clearParticles(), Cluster(), deleteParticles(), freezeInternalMotion(), getParticleList(), getParticles(), initializeParticles(), G4INCL::Nucleus::initializeParticles(), internalBoostToCM(), makeParticipant(), makeProjectileSpectator(), makeTargetSpectator(), print(), putParticlesOffShell(), removeParticle(), G4INCL::ProjectileRemnant::removeParticle(), rotateMomentum(), rotatePosition(), setPosition(), G4INCL::ProjectileRemnant::storeComponents(), G4INCL::ProjectileRemnant::storeEnergyLevels(), swap(), and updateClusterParameters().
|
protected |
Definition at line 452 of file G4INCLCluster.hh.
Referenced by G4INCL::Nucleus::computeRecoilKinematics(), getExcitationEnergy(), G4INCL::Nucleus::getExcitationEnergy(), setExcitationEnergy(), swap(), and G4INCL::Nucleus::useFusionKinematics().
|
protected |
Definition at line 454 of file G4INCLCluster.hh.
Referenced by Cluster(), initializeParticles(), G4INCL::Nucleus::Nucleus(), G4INCL::Nucleus::setDensity(), swap(), and ~Cluster().
|
protected |
Definition at line 453 of file G4INCLCluster.hh.
Referenced by G4INCL::Nucleus::computeRecoilKinematics(), getSpin(), setSpin(), and swap().