|
Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
|
#include <G4INCLParticle.hh>
Public Member Functions | |
| 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 a Strange? | |
| G4bool | isXi () const |
| Is this a Xi? | |
| G4bool | isAntiNucleon () const |
| Is this an antinucleon? | |
| G4bool | isAntiSigma () const |
| Is this an antiSigma? | |
| G4bool | isAntiXi () const |
| Is this an antiXi? | |
| G4bool | isAntiLambda () const |
| Is this an antiLambda? | |
| G4bool | isAntiHyperon () const |
| Is this an antiHyperon? | |
| G4bool | isAntiBaryon () const |
| Is this an antiBaryon? | |
| G4bool | isAntiNucleonorAntiLambda () const |
| Is this an antiNucleon or an antiLambda? | |
| G4int | getA () const |
| Returns the baryon number. | |
| G4int | getZ () const |
| Returns the charge number. | |
| G4int | getS () const |
| Returns the strangeness number. | |
| G4int | getSrcPair () 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 | setNumberOfSrcPair (int n) |
| Set the number of srcpairs. | |
| 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 | setSrcPartner () |
| Set and reset src partner. | |
| void | resetSrcPartner () |
| G4bool | isSrcPartner () const |
| Check if the particle is a src partner. | |
| 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) |
| G4int | getParentResonancePDGCode () const |
| void | setParentResonancePDGCode (const G4int parentPDGCode) |
| G4int | getParentResonanceID () const |
| void | setParentResonanceID (const G4int parentID) |
Static Public Member Functions | |
| 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 | |
| static std::vector< G4double > | INCLBiasVector |
| Time ordered vector of all bias applied. | |
| static G4ThreadLocal G4int | nextBiasedCollisionID = 0 |
Protected Member Functions | |
| void | swap (Particle &rhs) |
| Helper method for the assignment operator. | |
Protected Attributes | |
| 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 |
| G4int | nSrcPair |
| G4double | thePotentialEnergy |
| long | ID |
| G4bool | rpCorrelated |
| G4double | uncorrelatedMomentum |
| G4double | theParticleBias |
| G4int | theNKaon |
| The number of Kaons inside the nucleus (update during the cascade). | |
| G4int | theParentResonancePDGCode |
| G4int | theParentResonanceID |
Definition at line 75 of file G4INCLParticle.hh.
| G4INCL::Particle::Particle | ( | ) |
Definition at line 59 of file G4INCLParticle.cc.
Referenced by G4INCL::Cluster::addParticle(), G4INCL::Cluster::Cluster(), G4INCL::Cluster::Cluster(), G4INCL::Cluster::Cluster(), MergeVectorBias(), MergeVectorBias(), operator=(), Particle(), G4INCL::Cluster::removeParticle(), G4INCL::ProjectileRemnant::reset(), G4INCL::ProjectileRemnant::storeComponents(), and swap().
| G4INCL::Particle::Particle | ( | ParticleType | t, |
| G4double | energy, | ||
| ThreeVector const & | momentum, | ||
| ThreeVector const & | position ) |
Definition at line 92 of file G4INCLParticle.cc.
| G4INCL::Particle::Particle | ( | ParticleType | t, |
| ThreeVector const & | momentum, | ||
| ThreeVector const & | position ) |
Definition at line 124 of file G4INCLParticle.cc.
|
inlinevirtual |
Definition at line 80 of file G4INCLParticle.hh.
|
inline |
Copy constructor.
Does not copy the particle ID.
Definition at line 86 of file G4INCLParticle.hh.
| G4double G4INCL::Particle::adjustEnergyFromMomentum | ( | ) |
Recompute the energy to match the momentum.
Definition at line 170 of file G4INCLParticle.cc.
Referenced by G4INCL::Nucleus::computeRecoilKinematics(), G4INCL::Nucleus::decayOutgoingPionResonances(), G4INCL::Nucleus::decayOutgoingSigmaZero(), G4INCL::NKbElasticChannel::fillFinalState(), G4INCL::NKbToLpiChannel::fillFinalState(), G4INCL::NKbToNKbChannel::fillFinalState(), G4INCL::NKbToSpiChannel::fillFinalState(), G4INCL::NKElasticChannel::fillFinalState(), G4INCL::NKToNKChannel::fillFinalState(), G4INCL::NSToNSChannel::fillFinalState(), G4INCL::StrangeAbsorbtionChannel::fillFinalState(), G4INCL::PhaseSpaceKopylov::generate(), and G4INCL::Nucleus::restoreSrcPartner().
| const ThreeVector & G4INCL::Particle::adjustMomentumFromEnergy | ( | ) |
Rescale the momentum to match the total energy.
Definition at line 157 of file G4INCLParticle.cc.
Referenced by G4INCL::Cluster::Cluster(), G4INCL::Nucleus::emitInsideAntilambda(), G4INCL::Nucleus::emitInsideKaon(), G4INCL::Nucleus::emitInsideLambda(), G4INCL::Nucleus::emitInsidePions(), G4INCL::Nucleus::emitInsideStrangeParticles(), G4INCL::CrossSections::interactionDistanceKbarN(), G4INCL::CrossSections::interactionDistanceKN(), G4INCL::CrossSections::interactionDistanceNbarN(), G4INCL::CrossSections::interactionDistancenbarN(), G4INCL::CrossSections::interactionDistanceNN(), G4INCL::CrossSections::interactionDistancePiN(), G4INCL::CrossSections::interactionDistanceYN(), G4INCL::StandardPropagationModel::shootParticle(), and G4INCL::KinematicsUtils::transformToLocalEnergyFrame().
|
inline |
Boost the particle using a boost vector.
Example (go to the particle rest frame): particle->boost(particle->boostVector());
Definition at line 517 of file G4INCLParticle.hh.
Referenced by G4INCL::Cluster::boost(), G4INCL::Nucleus::decayOutgoingPionResonances(), G4INCL::Nucleus::decayOutgoingSigmaZero(), and G4INCL::PhaseSpaceKopylov::generate().
|
inline |
Returns a three vector we can give to the boost() -method.
In order to go to the particle rest frame you need to multiply the boost vector by -1.0.
Definition at line 507 of file G4INCLParticle.hh.
Referenced by G4INCL::PhaseSpaceKopylov::generate(), G4INCL::StandardPropagationModel::shootComposite(), and G4INCL::StandardPropagationModel::shootParticle().
|
inline |
Definition at line 1114 of file G4INCLParticle.hh.
|
static |
Definition at line 217 of file G4INCLParticle.cc.
Referenced by G4INCL::InteractionAvatar::postInteraction().
|
inline |
Freeze particle propagation.
Make the particle use theFrozenMomentum and theFrozenEnergy for propagation. The normal state can be restored by calling the thawPropagation() method.
Definition at line 1052 of file G4INCLParticle.hh.
|
inline |
Returns the baryon number.
Definition at line 485 of file G4INCLParticle.hh.
Referenced by G4INCL::Cluster::addParticle(), G4INCL::ClusteringModelIntercomparison::clusterCanEscape(), G4INCL::ClusterDecay::decay(), G4INCL::Nucleus::emitInsideAntilambda(), G4INCL::Nucleus::emitInsideLambda(), G4INCL::Nucleus::emitInsideStrangeParticles(), G4INCL::Nucleus::fillEventInfo(), G4INCL::ParticleEntryChannel::fillFinalState(), G4INCL::Nucleus::getConservationBalance(), G4INCL::NuclearDensity::getTransmissionRadius(), G4INCL::Nucleus::insertParticle(), G4INCL::ClusterDecay::isStable(), G4INCL::ProjectileRemnant::removeParticle(), G4INCL::StandardPropagationModel::shootComposite(), and G4INCL::StandardPropagationModel::shootParticle().
|
inlinevirtual |
Get the angular momentum w.r.t. the origin
Reimplemented in G4INCL::Cluster.
Definition at line 934 of file G4INCLParticle.hh.
Referenced by G4INCL::Cluster::getAngularMomentum(), G4INCL::StandardPropagationModel::shootCompositeAtrest(), and G4INCL::StandardPropagationModel::shootParticle().
|
inline |
Definition at line 496 of file G4INCLParticle.hh.
|
inline |
Get the vector list of biased vertices on the particle path.
Definition at line 1185 of file G4INCLParticle.hh.
Referenced by G4INCL::ClusterDecay::decay(), G4INCL::Nucleus::decayMe(), G4INCL::Nucleus::decayOutgoingClusters(), G4INCL::Nucleus::decayOutgoingPionResonances(), G4INCL::Nucleus::decayOutgoingSigmaZero(), MergeVectorBias(), and MergeVectorBias().
Definition at line 226 of file G4INCLParticle.cc.
Referenced by G4INCL::BinaryCollisionAvatar::getChannel(), G4INCL::ParticleList::getParticleListBias(), and setBiasCollisionVector().
|
inline |
Get the cosine of the angle between position and momentum.
Definition at line 1161 of file G4INCLParticle.hh.
|
inline |
Computes correction on the emission Q-value.
Computes the correction that must be applied to INCL particles in order to obtain the correct Q-value for particle emission from a given nucleus. For absorption, the correction is obviously equal to minus the value returned by this function.
| AParent | the mass number of the emitting nucleus |
| ZParent | the charge number of the emitting nucleus |
Definition at line 736 of file G4INCLParticle.hh.
Referenced by G4INCL::Nucleus::emitInsideAntilambda(), G4INCL::Nucleus::emitInsideKaon(), G4INCL::Nucleus::emitInsideLambda(), G4INCL::Nucleus::emitInsidePions(), G4INCL::Nucleus::emitInsideStrangeParticles(), and G4INCL::SurfaceAvatar::getTransmissionProbability().
|
inline |
Computes correction on the emission Q-value for hypernuclei.
Computes the correction that must be applied to INCL particles in order to obtain the correct Q-value for particle emission from a given nucleus. For absorption, the correction is obviously equal to minus the value returned by this function.
| AParent | the mass number of the emitting nucleus |
| ZParent | the charge number of the emitting nucleus |
| SParent | the strangess number of the emitting nucleus |
Definition at line 812 of file G4INCLParticle.hh.
|
inline |
Definition at line 1009 of file G4INCLParticle.hh.
|
inline |
Get the energy of the particle in MeV.
Definition at line 904 of file G4INCLParticle.hh.
Referenced by G4INCL::Cluster::addParticle(), G4INCL::DeltaDecayChannel::computeDecayTime(), G4INCL::PionResonanceDecayChannel::computeDecayTime(), G4INCL::EtaNElasticChannel::fillFinalState(), G4INCL::EtaNToPiNChannel::fillFinalState(), G4INCL::EtaOrOmegaNToLKChannel::fillFinalState(), G4INCL::EtaOrOmegaNToSKChannel::fillFinalState(), G4INCL::OmegaNElasticChannel::fillFinalState(), G4INCL::OmegaNToPiNChannel::fillFinalState(), G4INCL::ParticleEntryChannel::fillFinalState(), G4INCL::ClusteringModelIntercomparison::getCluster(), G4INCL::StandardPropagationModel::getTime(), G4INCL::KinematicsUtils::makeBoostVector(), G4INCL::KinematicsUtils::momentumInCM(), G4INCL::ProjectileRemnant::removeParticle(), G4INCL::StandardPropagationModel::shootAtrest(), G4INCL::StandardPropagationModel::shootComposite(), G4INCL::StandardPropagationModel::shootCompositeAtrest(), G4INCL::StandardPropagationModel::shootParticle(), G4INCL::KinematicsUtils::squareTotalEnergyInCM(), and G4INCL::KinematicsUtils::transformToLocalEnergyFrame().
|
inline |
|
inline |
|
inline |
Definition at line 960 of file G4INCLParticle.hh.
|
inline |
Definition at line 1128 of file G4INCLParticle.hh.
Referenced by G4INCL::PauliStandard::getBlockingProbability(), G4INCL::SurfaceAvatar::getChannel(), G4INCL::ClusteringModelIntercomparison::getCluster(), and G4INCL::ProjectileRemnant::reset().
|
inline |
Get the INCL particle mass.
Definition at line 549 of file G4INCLParticle.hh.
Referenced by getEmissionQValueCorrection(), getEmissionQValueCorrection(), G4INCL::StandardPropagationModel::getTime(), and setINCLMass().
|
inline |
Get the the particle invariant mass.
Uses the relativistic invariant
![\[ m = \sqrt{E^2 - {\vec p}^2}\]](../../form_40.png)
Definition at line 882 of file G4INCLParticle.hh.
Referenced by Particle().
|
inline |
Get the particle kinetic energy.
Definition at line 893 of file G4INCLParticle.hh.
Referenced by G4INCL::NuclearPotential::NuclearPotentialEnergyIsospin::computePotentialEnergy(), G4INCL::NuclearPotential::NuclearPotentialEnergyIsospinSmooth::computePotentialEnergy(), G4INCL::Nucleus::emitInsideAntilambda(), G4INCL::Nucleus::emitInsideKaon(), G4INCL::Nucleus::emitInsideLambda(), G4INCL::Nucleus::emitInsidePions(), G4INCL::Nucleus::emitInsideStrangeParticles(), G4INCL::Nucleus::fillEventInfo(), G4INCL::Nucleus::getConservationBalance(), G4INCL::KinematicsUtils::getLocalEnergy(), G4INCL::StandardPropagationModel::getTime(), G4INCL::SurfaceAvatar::getTransmissionProbability(), G4INCL::CoulombDistortion::maxImpactParameter(), G4INCL::CDPP::processOneParticle(), G4INCL::StandardPropagationModel::shootCompositeAtrest(), and G4INCL::StandardPropagationModel::shootParticle().
|
inline |
Longitudinal component of the position w.r.t. the momentum.
Definition at line 1017 of file G4INCLParticle.hh.
Referenced by getTransversePosition().
|
inline |
Get the cached particle mass.
Definition at line 546 of file G4INCLParticle.hh.
Referenced by G4INCL::DeltaDecayChannel::computeDecayTime(), G4INCL::PionResonanceDecayChannel::computeDecayTime(), G4INCL::SigmaZeroDecayChannel::computeDecayTime(), G4INCL::Nucleus::emitInsideAntilambda(), G4INCL::Nucleus::emitInsideKaon(), G4INCL::Nucleus::emitInsideLambda(), G4INCL::Nucleus::emitInsidePions(), G4INCL::Nucleus::emitInsideStrangeParticles(), G4INCL::InteractionAvatar::enforceEnergyConservation(), G4INCL::CrossSectionsMultiPionsAndResonances::etaNToPiN(), G4INCL::EtaNElasticChannel::fillFinalState(), G4INCL::EtaNToPiNChannel::fillFinalState(), G4INCL::EtaOrOmegaNToLKChannel::fillFinalState(), G4INCL::EtaOrOmegaNToSKChannel::fillFinalState(), G4INCL::NSToNLChannel::fillFinalState(), G4INCL::NSToNSChannel::fillFinalState(), G4INCL::OmegaNElasticChannel::fillFinalState(), G4INCL::OmegaNToPiNChannel::fillFinalState(), G4INCL::Cluster::freezeInternalMotion(), G4INCL::PhaseSpaceKopylov::generate(), G4INCL::NuclearPotential::INuclearPotential::getFermiMomentum(), G4INCL::KinematicsUtils::getLocalEnergy(), G4INCL::SurfaceAvatar::getTransmissionProbability(), G4INCL::CrossSections::interactionDistanceKbarN(), G4INCL::CrossSections::interactionDistanceKN(), G4INCL::CrossSections::interactionDistanceNbarN(), G4INCL::CrossSections::interactionDistancenbarN(), G4INCL::CrossSections::interactionDistanceNN(), G4INCL::CrossSections::interactionDistancePiN(), G4INCL::CrossSections::interactionDistanceYN(), G4INCL::Cluster::internalBoostToCM(), G4INCL::KinematicsUtils::momentumInCM(), G4INCL::KinematicsUtils::momentumInLab(), G4INCL::CrossSectionsINCL46::NDeltaToNN(), G4INCL::CrossSectionsMultiPions::NDeltaToNN(), G4INCL::CrossSectionsStrangeness::NNToNLK2pi(), G4INCL::CrossSectionsStrangeness::NNToNLKpi(), G4INCL::CrossSectionsStrangeness::NNToNSK2pi(), G4INCL::CrossSectionsStrangeness::NNToNSKpi(), G4INCL::CrossSectionsStrangeness::omegaNToLK(), G4INCL::CrossSectionsMultiPionsAndResonances::omegaNToPiN(), G4INCL::CrossSectionsStrangeness::omegaNToSK(), G4INCL::CrossSectionsMultiPionsAndResonances::piMinuspToEtaN(), G4INCL::CrossSectionsMultiPionsAndResonances::piMinuspToOmegaN(), G4INCL::Cluster::print(), G4INCL::ProjectileRemnant::ProjectileRemnant(), G4INCL::StandardPropagationModel::shootCompositeAtrest(), and G4INCL::StandardPropagationModel::shootParticle().
|
inline |
Get the momentum vector.
Definition at line 928 of file G4INCLParticle.hh.
Referenced by G4INCL::Cluster::addParticle(), G4INCL::ClusteringModelIntercomparison::clusterCanEscape(), G4INCL::SigmaZeroDecayChannel::computeDecayTime(), G4INCL::Nucleus::decayOutgoingPionResonances(), G4INCL::Nucleus::decayOutgoingSigmaZero(), G4INCL::Nucleus::fillEventInfo(), G4INCL::NNToNSKpiChannel::fillFinalState(), G4INCL::NpiToLK2piChannel::fillFinalState(), G4INCL::NpiToLKpiChannel::fillFinalState(), G4INCL::NpiToNKKbChannel::fillFinalState(), G4INCL::NpiToSK2piChannel::fillFinalState(), G4INCL::NpiToSKpiChannel::fillFinalState(), G4INCL::ParticleEntryChannel::fillFinalState(), G4INCL::PauliStandard::getBlockingProbability(), G4INCL::ClusteringModelIntercomparison::getCluster(), G4INCL::Nucleus::getConservationBalance(), G4INCL::KinematicsUtils::getLocalEnergy(), G4INCL::SurfaceAvatar::getTransmissionProbability(), G4INCL::NKbElasticChannel::KaonMomentum(), G4INCL::NKbToLpiChannel::KaonMomentum(), G4INCL::NKbToNKbChannel::KaonMomentum(), G4INCL::NKbToSpiChannel::KaonMomentum(), G4INCL::KinematicsUtils::makeBoostVector(), G4INCL::KinematicsUtils::momentumInCM(), G4INCL::ProjectileRemnant::removeParticle(), G4INCL::ParticleSampler::sampleParticlesIntoList(), G4INCL::StandardPropagationModel::shootComposite(), G4INCL::StandardPropagationModel::shootCompositeAtrest(), and G4INCL::StandardPropagationModel::shootParticle().
|
inline |
Return the number of collisions undergone by the particle.
Definition at line 968 of file G4INCLParticle.hh.
Referenced by G4INCL::Cluster::addParticle().
|
inline |
Return the number of decays undergone by the particle.
Definition at line 977 of file G4INCLParticle.hh.
|
inline |
Number of Kaon inside de nucleus.
Put in the Particle class in order to calculate the "correct" mass of composit particle.
Definition at line 1200 of file G4INCLParticle.hh.
|
inline |
Definition at line 1206 of file G4INCLParticle.hh.
|
inline |
Definition at line 1204 of file G4INCLParticle.hh.
|
inline |
Definition at line 378 of file G4INCLParticle.hh.
Referenced by G4INCL::StandardPropagationModel::generateBinaryCollisionAvatar().
|
inline |
|
inline |
Return a NULL pointer
Definition at line 1133 of file G4INCLParticle.hh.
|
inline |
Set the position vector.
Definition at line 950 of file G4INCLParticle.hh.
Referenced by G4INCL::Cluster::addParticle(), G4INCL::InteractionAvatar::bringParticleInside(), G4INCL::CoulombNone::bringToSurface(), G4INCL::CoulombNone::bringToSurfaceAbar(), G4INCL::ClusteringModelIntercomparison::clusterCanEscape(), G4INCL::EtaNToPiPiNChannel::fillFinalState(), G4INCL::NKbToNKb2piChannel::fillFinalState(), G4INCL::NKToNK2piChannel::fillFinalState(), G4INCL::NpiToMissingStrangenessChannel::fillFinalState(), G4INCL::OmegaNToPiPiNChannel::fillFinalState(), G4INCL::StandardPropagationModel::generateBinaryCollisionAvatar(), G4INCL::PauliStandard::getBlockingProbability(), G4INCL::ClusteringModelIntercomparison::getCluster(), G4INCL::KinematicsUtils::getLocalEnergy(), G4INCL::StandardPropagationModel::getReflectionTime(), G4INCL::StandardPropagationModel::getTime(), and G4INCL::ParticleSampler::sampleParticlesIntoList().
|
inline |
Get the particle potential energy.
Definition at line 896 of file G4INCLParticle.hh.
Referenced by G4INCL::Cluster::addParticle(), G4INCL::Nucleus::emitInsideAntilambda(), G4INCL::Nucleus::emitInsideKaon(), G4INCL::Nucleus::emitInsideLambda(), G4INCL::Nucleus::emitInsidePions(), G4INCL::Nucleus::emitInsideStrangeParticles(), G4INCL::ClusteringModelIntercomparison::getCluster(), G4INCL::KinematicsUtils::getLocalEnergy(), G4INCL::StandardPropagationModel::getTime(), G4INCL::SurfaceAvatar::getTransmissionProbability(), and G4INCL::CDPP::processOneParticle().
|
inline |
Get the propagation velocity of the particle.
Definition at line 1044 of file G4INCLParticle.hh.
Referenced by G4INCL::CoulombNone::bringToSurface(), G4INCL::CoulombNone::bringToSurfaceAbar(), G4INCL::StandardPropagationModel::getReflectionTime(), and G4INCL::StandardPropagationModel::getTime().
|
inline |
Get the real particle mass.
Definition at line 661 of file G4INCLParticle.hh.
Referenced by G4INCL::Cluster::getTableMass(), and setRealMass().
|
inline |
Return the reflection momentum.
The reflection momentum is used by calls to getSurfaceRadius to compute the radius of the sphere where the nucleon moves. It is necessary to introduce fuzzy r-p correlations.
Definition at line 1144 of file G4INCLParticle.hh.
Referenced by G4INCL::KinematicsUtils::getLocalEnergy(), and G4INCL::Nucleus::getSurfaceRadius().
|
inline |
Returns the strangeness number.
Definition at line 491 of file G4INCLParticle.hh.
Referenced by G4INCL::Cluster::addParticle(), G4INCL::ClusteringModelIntercomparison::clusterCanEscape(), G4INCL::ClusterDecay::decay(), G4INCL::Nucleus::emitInsideAntilambda(), G4INCL::Nucleus::emitInsideKaon(), G4INCL::Nucleus::emitInsideLambda(), G4INCL::Nucleus::emitInsideStrangeParticles(), G4INCL::Nucleus::fillEventInfo(), G4INCL::ParticleEntryChannel::fillFinalState(), G4INCL::Nucleus::getConservationBalance(), G4INCL::Nucleus::insertParticle(), G4INCL::ClusterDecay::isStable(), and G4INCL::ProjectileRemnant::removeParticle().
|
inlinevirtual |
Get the particle species.
Reimplemented in G4INCL::Cluster.
Definition at line 196 of file G4INCLParticle.hh.
Referenced by G4INCL::CoulombDistortion::maxImpactParameter(), and G4INCL::StandardPropagationModel::shootParticle().
|
inline |
|
inlinevirtual |
Get the tabulated particle mass.
Reimplemented in G4INCL::Cluster.
Definition at line 605 of file G4INCLParticle.hh.
Referenced by G4INCL::Nucleus::decayOutgoingPionResonances(), G4INCL::Nucleus::decayOutgoingSigmaZero(), getEmissionQValueCorrection(), getEmissionQValueCorrection(), getTransferQValueCorrection(), getTransferQValueCorrection(), and setTableMass().
|
static |
General bias vector function.
Definition at line 315 of file G4INCLParticle.cc.
Referenced by G4INCL::Nucleus::emitInsideAntilambda(), G4INCL::Nucleus::emitInsideKaon(), G4INCL::Nucleus::emitInsideLambda(), G4INCL::Nucleus::emitInsidePions(), G4INCL::Nucleus::emitInsideStrangeParticles(), and G4INCL::EventInfo::remnantToParticle().
|
inline |
Computes correction on the transfer Q-value.
Computes the correction that must be applied to INCL particles in order to obtain the correct Q-value for particle transfer from a given nucleus to another.
Assumes that the receving nucleus is INCL's target nucleus, with the INCL separation energy.
| AFrom | the mass number of the donating nucleus |
| ZFrom | the charge number of the donating nucleus |
| ATo | the mass number of the receiving nucleus |
| ZTo | the charge number of the receiving nucleus |
Definition at line 776 of file G4INCLParticle.hh.
|
inline |
Computes correction on the transfer Q-value for hypernuclei.
Computes the correction that must be applied to INCL particles in order to obtain the correct Q-value for particle transfer from a given nucleus to another.
Assumes that the receving nucleus is INCL's target nucleus, with the INCL separation energy.
| AFrom | the mass number of the donating nucleus |
| ZFrom | the charge number of the donating nucleus |
| SFrom | the strangess number of the donating nucleus |
| ATo | the mass number of the receiving nucleus |
| ZTo | the charge number of the receiving nucleus |
| STo | the strangess number of the receiving nucleus |
Definition at line 853 of file G4INCLParticle.hh.
|
inline |
Transverse component of the position w.r.t. the momentum.
Definition at line 1012 of file G4INCLParticle.hh.
Referenced by G4INCL::StandardPropagationModel::shootComposite(), G4INCL::StandardPropagationModel::shootCompositeAtrest(), and G4INCL::StandardPropagationModel::shootParticle().
|
inline |
Get the particle type.
Definition at line 191 of file G4INCLParticle.hh.
Referenced by G4INCL::PionResonanceDecayChannel::computeDecayTime(), G4INCL::NuclearPotential::INuclearPotential::computeKaonPotentialEnergy(), G4INCL::NuclearPotential::INuclearPotential::computePionPotentialEnergy(), G4INCL::NuclearPotential::INuclearPotential::computePionResonancePotentialEnergy(), G4INCL::NuclearPotential::NuclearPotentialConstant::computePotentialEnergy(), G4INCL::NuclearPotential::NuclearPotentialIsospin::computePotentialEnergy(), G4INCL::CrossSectionsINCL46::elasticNNLegacy(), G4INCL::NDeltaToDeltaLKChannel::fillFinalState(), G4INCL::NDeltaToDeltaSKChannel::fillFinalState(), G4INCL::NDeltaToNSKChannel::fillFinalState(), G4INCL::NKbToL2piChannel::fillFinalState(), G4INCL::NKbToLpiChannel::fillFinalState(), G4INCL::NKbToNKb2piChannel::fillFinalState(), G4INCL::NKbToNKbChannel::fillFinalState(), G4INCL::NKbToNKbpiChannel::fillFinalState(), G4INCL::NKbToS2piChannel::fillFinalState(), G4INCL::NKbToSpiChannel::fillFinalState(), G4INCL::NKToNK2piChannel::fillFinalState(), G4INCL::NKToNKChannel::fillFinalState(), G4INCL::NKToNKpiChannel::fillFinalState(), G4INCL::NpiToMissingStrangenessChannel::fillFinalState(), G4INCL::NSToNLChannel::fillFinalState(), G4INCL::NSToNSChannel::fillFinalState(), G4INCL::PauliStandard::getBlockingProbability(), G4INCL::NuclearPotential::INuclearPotential::getFermiEnergy(), G4INCL::NuclearPotential::INuclearPotential::getFermiMomentum(), G4INCL::KinematicsUtils::getLocalEnergy(), G4INCL::NuclearPotential::INuclearPotential::getSeparationEnergy(), G4INCL::Nucleus::getSurfaceRadius(), G4INCL::StandardPropagationModel::getTime(), G4INCL::NuclearDensity::getTransmissionRadius(), G4INCL::Nucleus::insertParticle(), G4INCL::ParticleConfig::isPair(), G4INCL::CrossSectionsStrangeness::NDeltaToDeltaLK(), G4INCL::CrossSectionsStrangeness::NDeltaToDeltaSK(), G4INCL::CrossSectionsStrangeness::NDeltaToNLK(), G4INCL::CrossSectionsINCL46::NDeltaToNN(), G4INCL::CrossSectionsMultiPions::NDeltaToNN(), G4INCL::CrossSectionsStrangeness::NDeltaToNNKKb(), G4INCL::CrossSectionsStrangeness::NDeltaToNSK(), G4INCL::CrossSectionsStrangeness::NKbToL2pi(), G4INCL::CrossSectionsStrangeness::NKbToLpi(), G4INCL::CrossSectionsStrangeness::NKbToNKb(), G4INCL::CrossSectionsStrangeness::NKbToNKb2pi(), G4INCL::CrossSectionsStrangeness::NKbToNKbpi(), G4INCL::CrossSectionsStrangeness::NKbToS2pi(), G4INCL::CrossSectionsStrangeness::NKbToSpi(), G4INCL::CrossSectionsStrangeness::NKToNK(), G4INCL::CrossSectionsStrangeness::NKToNK2pi(), G4INCL::CrossSectionsStrangeness::NKToNKpi(), G4INCL::CrossSectionsAntiparticles::NNbarCEX(), G4INCL::CrossSectionsAntiparticles::NNbarElastic(), G4INCL::CrossSectionsAntiparticles::NNbarToAnnihilation(), G4INCL::CrossSectionsAntiparticles::NNbarToLLbar(), G4INCL::CrossSectionsAntiparticles::NNbarToNNbar2pi(), G4INCL::CrossSectionsAntiparticles::NNbarToNNbar3pi(), G4INCL::CrossSectionsAntiparticles::NNbarToNNbarpi(), G4INCL::CrossSectionsMultiPions::NNElastic(), G4INCL::CrossSectionsMultiPions::NNOnePi(), G4INCL::CrossSectionsMultiPions::NNOnePiOrDelta(), G4INCL::CrossSectionsMultiPions::NNThreePi(), G4INCL::CrossSectionsStrangeness::NNToMissingStrangeness(), G4INCL::CrossSectionsINCL46::NNToNDelta(), G4INCL::CrossSectionsMultiPions::NNToNDelta(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNDeltaEta(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNDeltaOmega(), G4INCL::CrossSectionsStrangeness::NNToNLK(), G4INCL::CrossSectionsStrangeness::NNToNLK2pi(), G4INCL::CrossSectionsStrangeness::NNToNLKpi(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNEta(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNEtaExclu(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNEtaFourPi(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNEtaOnePi(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNEtaOnePiOrDelta(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNEtaThreePi(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNEtaTwoPi(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNEtaxPi(), G4INCL::CrossSectionsStrangeness::NNToNNKKb(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNOmega(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNOmegaExclu(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNOmegaFourPi(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNOmegaOnePi(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNOmegaOnePiOrDelta(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNOmegaThreePi(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNOmegaTwoPi(), G4INCL::CrossSectionsMultiPionsAndResonances::NNToNNOmegaxPi(), G4INCL::CrossSectionsStrangeness::NNToNSK(), G4INCL::CrossSectionsStrangeness::NNToNSK2pi(), G4INCL::CrossSectionsStrangeness::NNToNSKpi(), G4INCL::CrossSectionsMultiPions::NNTot(), G4INCL::CrossSectionsMultiPions::NNTwoPi(), G4INCL::CrossSectionsStrangeness::NpiToLK(), G4INCL::CrossSectionsStrangeness::NpiToLK2pi(), G4INCL::CrossSectionsStrangeness::NpiToLKpi(), G4INCL::CrossSectionsStrangeness::NpiToNKKb(), G4INCL::CrossSectionsStrangeness::NpiToSK(), G4INCL::CrossSectionsStrangeness::NpiToSK2pi(), G4INCL::CrossSectionsStrangeness::NpiToSKpi(), G4INCL::CrossSectionsStrangeness::NSToNL(), G4INCL::CrossSectionsStrangeness::NSToNS(), G4INCL::CrossSectionsStrangeness::p_pimToSzKz(), G4INCL::CrossSectionsINCL46::piNToDelta(), G4INCL::CrossSectionsMultiPions::piNToDelta(), G4INCL::CrossSectionsMultiPionsAndResonances::piNToEtaN(), G4INCL::CrossSectionsMultiPionsAndResonances::piNToOmegaN(), G4INCL::CrossSectionsMultiPions::piNTot(), and G4INCL::CrossSectionsAntiparticles::total().
|
inline |
Returns the charge number.
Definition at line 488 of file G4INCLParticle.hh.
Referenced by G4INCL::Cluster::addParticle(), G4INCL::CoulombNonRelativistic::bringToSurface(), G4INCL::CoulombNonRelativistic::bringToSurfaceAbar(), G4INCL::ClusterDecay::decay(), G4INCL::Nucleus::emitInsideAnnihilationProducts(), G4INCL::Nucleus::emitInsideKaon(), G4INCL::Nucleus::emitInsidePions(), G4INCL::Nucleus::emitInsideStrangeParticles(), G4INCL::Nucleus::fillEventInfo(), G4INCL::ParticleEntryChannel::fillFinalState(), G4INCL::ClusteringModelIntercomparison::getCluster(), G4INCL::Nucleus::getConservationBalance(), G4INCL::Nucleus::getTransmissionBarrier(), G4INCL::SurfaceAvatar::getTransmissionProbability(), G4INCL::NuclearDensity::getTransmissionRadius(), G4INCL::Nucleus::insertParticle(), G4INCL::ClusterDecay::isStable(), G4INCL::ProjectileRemnant::removeParticle(), and G4INCL::StandardPropagationModel::shootComposite().
|
inline |
Increment the number of collisions undergone by the particle.
Definition at line 974 of file G4INCLParticle.hh.
|
inline |
Increment the number of decays undergone by the particle.
Definition at line 983 of file G4INCLParticle.hh.
|
inline |
Is this an antiBaryon?
Definition at line 479 of file G4INCLParticle.hh.
|
inline |
Is this an antiHyperon?
Definition at line 476 of file G4INCLParticle.hh.
Referenced by isAntiBaryon().
|
inline |
Is this an antiKaon?
Definition at line 440 of file G4INCLParticle.hh.
Referenced by G4INCL::CrossSectionsAntiparticles::elastic(), G4INCL::CrossSectionsStrangeness::elastic(), G4INCL::Nucleus::insertParticle(), isMeson(), isStrange(), G4INCL::CrossSectionsStrangeness::NKbelastic(), G4INCL::CrossSectionsStrangeness::NKbToL2pi(), G4INCL::CrossSectionsStrangeness::NKbToLpi(), G4INCL::CrossSectionsStrangeness::NKbToNKb(), G4INCL::CrossSectionsStrangeness::NKbToNKb2pi(), G4INCL::CrossSectionsStrangeness::NKbToNKbpi(), G4INCL::CrossSectionsStrangeness::NKbToS2pi(), G4INCL::CrossSectionsStrangeness::NKbToSpi(), G4INCL::CrossSectionsAntiparticles::total(), and G4INCL::CrossSectionsStrangeness::total().
|
inline |
Is this an antiLambda?
Definition at line 473 of file G4INCLParticle.hh.
Referenced by isAntiHyperon(), and isAntiNucleonorAntiLambda().
|
inline |
Is this an antinucleon?
Definition at line 464 of file G4INCLParticle.hh.
Referenced by G4INCL::CrossSectionsAntiparticles::elastic(), G4INCL::StandardPropagationModel::generateBinaryCollisionAvatar(), G4INCL::Nucleus::insertParticle(), isAntiBaryon(), isAntiNucleonorAntiLambda(), G4INCL::CrossSectionsAntiparticles::NNbarCEX(), G4INCL::CrossSectionsAntiparticles::NNbarElastic(), G4INCL::CrossSectionsAntiparticles::NNbarToAnnihilation(), G4INCL::CrossSectionsAntiparticles::NNbarToLLbar(), G4INCL::CrossSectionsAntiparticles::NNbarToNNbar2pi(), G4INCL::CrossSectionsAntiparticles::NNbarToNNbar3pi(), G4INCL::CrossSectionsAntiparticles::NNbarToNNbarpi(), G4INCL::InteractionAvatar::preInteractionLocalEnergy(), G4INCL::StandardPropagationModel::shootParticle(), and G4INCL::CrossSectionsAntiparticles::total().
|
inline |
|
inline |
Is this an antiSigma?
Definition at line 467 of file G4INCLParticle.hh.
Referenced by isAntiHyperon().
|
inline |
|
inline |
Is this a Baryon?
Definition at line 455 of file G4INCLParticle.hh.
Referenced by G4INCL::CDPP::processOneParticle().
|
inline |
Definition at line 1027 of file G4INCLParticle.hh.
Referenced by getEmissionQValueCorrection(), getEmissionQValueCorrection(), and G4INCL::SurfaceAvatar::postInteraction().
|
inline |
Is it a Delta?
Definition at line 429 of file G4INCLParticle.hh.
Referenced by G4INCL::CrossSectionsAntiparticles::elastic(), G4INCL::CrossSectionsINCL46::elastic(), G4INCL::CrossSectionsMultiPions::elastic(), G4INCL::CrossSectionsMultiPionsAndResonances::elastic(), G4INCL::CrossSectionsStrangeness::elastic(), G4INCL::CrossSectionsTruncatedMultiPions::elastic(), G4INCL::NuclearPotential::INuclearPotential::getFermiMomentum(), isResonance(), G4INCL::CrossSectionsINCL46::NDeltaToNN(), G4INCL::CrossSectionsMultiPions::NDeltaToNN(), G4INCL::CrossSectionsMultiPions::NNTot(), G4INCL::RecombinationChannel::RecombinationChannel(), G4INCL::CrossSectionsAntiparticles::total(), G4INCL::CrossSectionsINCL46::total(), G4INCL::CrossSectionsMultiPions::total(), G4INCL::CrossSectionsMultiPionsAndResonances::total(), and G4INCL::CrossSectionsStrangeness::total().
|
inline |
Is this an eta?
Definition at line 414 of file G4INCLParticle.hh.
Referenced by G4INCL::CrossSectionsAntiparticles::elastic(), G4INCL::CrossSectionsMultiPionsAndResonances::elastic(), G4INCL::CrossSectionsStrangeness::elastic(), G4INCL::CrossSectionsMultiPionsAndResonances::etaNElastic(), G4INCL::CrossSectionsStrangeness::etaNToLK(), G4INCL::CrossSectionsMultiPionsAndResonances::etaNToPiN(), G4INCL::CrossSectionsMultiPionsAndResonances::etaNToPiPiN(), G4INCL::CrossSectionsStrangeness::etaNToSK(), isMeson(), G4INCL::CrossSectionsAntiparticles::total(), G4INCL::CrossSectionsMultiPionsAndResonances::total(), and G4INCL::CrossSectionsStrangeness::total().
|
inline |
Is this an etaprime?
Definition at line 420 of file G4INCLParticle.hh.
Referenced by isMeson(), G4INCL::CrossSectionsAntiparticles::total(), G4INCL::CrossSectionsMultiPionsAndResonances::total(), and G4INCL::CrossSectionsStrangeness::total().
|
inline |
Is this an Hyperon?
Definition at line 449 of file G4INCLParticle.hh.
Referenced by G4INCL::CrossSectionsAntiparticles::elastic(), G4INCL::CrossSectionsStrangeness::elastic(), isBaryon(), isStrange(), and G4INCL::CrossSectionsStrangeness::NYelastic().
|
inline |
Is this a Kaon?
Definition at line 437 of file G4INCLParticle.hh.
Referenced by G4INCL::CrossSectionsAntiparticles::elastic(), G4INCL::CrossSectionsStrangeness::elastic(), G4INCL::Nucleus::insertParticle(), isMeson(), isStrange(), G4INCL::CrossSectionsStrangeness::NKelastic(), G4INCL::CrossSectionsStrangeness::NKToNK(), G4INCL::CrossSectionsStrangeness::NKToNK2pi(), G4INCL::CrossSectionsStrangeness::NKToNKpi(), G4INCL::CrossSectionsAntiparticles::total(), and G4INCL::CrossSectionsStrangeness::total().
|
inline |
Is this a Lambda?
Definition at line 443 of file G4INCLParticle.hh.
Referenced by G4INCL::NpiToMissingStrangenessChannel::fillFinalState(), G4INCL::Nucleus::getSurfaceRadius(), G4INCL::Nucleus::insertParticle(), isHyperon(), isNucleonorLambda(), G4INCL::CrossSectionsStrangeness::NLToNS(), G4INCL::CrossSectionsAntiparticles::total(), and G4INCL::CrossSectionsStrangeness::total().
|
inline |
Is this a Meson?
Definition at line 452 of file G4INCLParticle.hh.
Referenced by G4INCL::StandardPropagationModel::generateBinaryCollisionAvatar(), G4INCL::InteractionAvatar::preInteractionLocalEnergy(), G4INCL::CDPP::processOneParticle(), and G4INCL::StandardPropagationModel::shootParticle().
|
inline |
Is this a nucleon?
Definition at line 371 of file G4INCLParticle.hh.
Referenced by G4INCL::NuclearPotential::NuclearPotentialEnergyIsospin::computePotentialEnergy(), G4INCL::NuclearPotential::NuclearPotentialEnergyIsospinSmooth::computePotentialEnergy(), G4INCL::CrossSectionsAntiparticles::elastic(), G4INCL::CrossSectionsINCL46::elastic(), G4INCL::CrossSectionsMultiPions::elastic(), G4INCL::CrossSectionsMultiPionsAndResonances::elastic(), G4INCL::CrossSectionsStrangeness::elastic(), G4INCL::CrossSectionsTruncatedMultiPions::elastic(), G4INCL::CrossSectionsINCL46::elasticNNLegacy(), G4INCL::StandardPropagationModel::generateBinaryCollisionAvatar(), G4INCL::Nucleus::getSurfaceRadius(), G4INCL::Nucleus::insertParticle(), isBaryon(), isNucleonorLambda(), G4INCL::CrossSectionsMultiPions::NNElastic(), G4INCL::CrossSectionsMultiPions::NNTot(), G4INCL::CrossSectionsStrangeness::omegaNToLK(), G4INCL::CrossSectionsStrangeness::omegaNToSK(), G4INCL::CrossSectionsMultiPions::piMinuspIne(), G4INCL::CrossSectionsMultiPions::piMinuspOnePi(), G4INCL::CrossSectionsMultiPions::piMinuspTwoPi(), G4INCL::CrossSectionsMultiPions::piNIne(), G4INCL::CrossSectionsMultiPions::piNOnePi(), G4INCL::CrossSectionsINCL46::piNToDelta(), G4INCL::CrossSectionsMultiPions::piNToxPiN(), G4INCL::CrossSectionsMultiPions::piNTwoPi(), G4INCL::CrossSectionsMultiPions::piPluspIne(), G4INCL::CrossSectionsMultiPions::piPluspOnePi(), G4INCL::CrossSectionsMultiPions::piPluspTwoPi(), G4INCL::CrossSectionsAntiparticles::total(), G4INCL::CrossSectionsINCL46::total(), G4INCL::CrossSectionsMultiPions::total(), G4INCL::CrossSectionsMultiPionsAndResonances::total(), and G4INCL::CrossSectionsStrangeness::total().
|
inline |
Is this a Nucleon or a Lambda?
Definition at line 446 of file G4INCLParticle.hh.
|
inline |
Is this an omega?
Definition at line 417 of file G4INCLParticle.hh.
Referenced by G4INCL::CrossSectionsMultiPionsAndResonances::elastic(), isMeson(), G4INCL::CrossSectionsMultiPionsAndResonances::omegaNElastic(), G4INCL::CrossSectionsMultiPionsAndResonances::omegaNInelastic(), G4INCL::CrossSectionsStrangeness::omegaNToLK(), G4INCL::CrossSectionsMultiPionsAndResonances::omegaNToPiN(), G4INCL::CrossSectionsStrangeness::omegaNToSK(), G4INCL::CrossSectionsAntiparticles::total(), G4INCL::CrossSectionsMultiPionsAndResonances::total(), and G4INCL::CrossSectionsStrangeness::total().
|
inline |
Check if the particle is out of its potential well.
Definition at line 999 of file G4INCLParticle.hh.
Referenced by G4INCL::NuclearPotential::INuclearPotential::computeKaonPotentialEnergy(), G4INCL::NuclearPotential::INuclearPotential::computePionPotentialEnergy(), and G4INCL::NuclearPotential::INuclearPotential::computePionResonancePotentialEnergy().
|
inline |
Definition at line 386 of file G4INCLParticle.hh.
Referenced by G4INCL::StandardPropagationModel::generateBinaryCollisionAvatar().
|
inline |
Is this a photon?
Definition at line 423 of file G4INCLParticle.hh.
Referenced by G4INCL::StandardPropagationModel::generateBinaryCollisionAvatar(), and G4INCL::InteractionAvatar::preInteractionLocalEnergy().
|
inline |
Is this a pion?
Definition at line 411 of file G4INCLParticle.hh.
Referenced by G4INCL::CrossSectionsAntiparticles::elastic(), G4INCL::CrossSectionsMultiPions::elastic(), G4INCL::CrossSectionsMultiPionsAndResonances::elastic(), G4INCL::CrossSectionsStrangeness::elastic(), G4INCL::CrossSectionsTruncatedMultiPions::elastic(), G4INCL::StandardPropagationModel::generateBinaryCollisionAvatar(), G4INCL::Nucleus::insertParticle(), isMeson(), G4INCL::CrossSectionsStrangeness::NpiToLK(), G4INCL::CrossSectionsStrangeness::NpiToLK2pi(), G4INCL::CrossSectionsStrangeness::NpiToLKpi(), G4INCL::CrossSectionsStrangeness::NpiToMissingStrangeness(), G4INCL::CrossSectionsStrangeness::NpiToNKKb(), G4INCL::CrossSectionsStrangeness::NpiToSK(), G4INCL::CrossSectionsStrangeness::NpiToSK2pi(), G4INCL::CrossSectionsStrangeness::NpiToSKpi(), G4INCL::CrossSectionsMultiPionsAndResonances::piMinuspToEtaN(), G4INCL::CrossSectionsMultiPionsAndResonances::piMinuspToOmegaN(), G4INCL::CrossSectionsINCL46::piNToDelta(), G4INCL::CrossSectionsMultiPions::piNToDelta(), G4INCL::CrossSectionsMultiPions::piNTot(), G4INCL::CrossSectionsAntiparticles::total(), G4INCL::CrossSectionsINCL46::total(), G4INCL::CrossSectionsMultiPions::total(), G4INCL::CrossSectionsMultiPionsAndResonances::total(), and G4INCL::CrossSectionsStrangeness::total().
|
inline |
Definition at line 394 of file G4INCLParticle.hh.
|
inline |
Is it a resonance?
Definition at line 426 of file G4INCLParticle.hh.
Referenced by G4INCL::StandardPropagationModel::generateBinaryCollisionAvatar(), G4INCL::Nucleus::getSurfaceRadius(), isBaryon(), Particle(), G4INCL::CDPP::processOneParticle(), and setType().
|
inline |
Is this a Sigma?
Definition at line 434 of file G4INCLParticle.hh.
Referenced by isHyperon(), G4INCL::CrossSectionsStrangeness::NSToNL(), G4INCL::CrossSectionsStrangeness::NSToNS(), G4INCL::CrossSectionsAntiparticles::total(), and G4INCL::CrossSectionsStrangeness::total().
|
inline |
Check if the particle is a src partner.
Definition at line 1006 of file G4INCLParticle.hh.
|
inline |
|
inline |
Definition at line 390 of file G4INCLParticle.hh.
Referenced by G4INCL::Nucleus::insertParticle().
|
inline |
|
inline |
Lorentz-contract the particle position around some center.
Apply Lorentz contraction to the position component along the direction of the boost vector.
| aBoostVector | the boost vector (velocity) [c] |
| refPos | the reference position |
Definition at line 535 of file G4INCLParticle.hh.
|
inlinevirtual |
Reimplemented in G4INCL::Cluster.
Definition at line 398 of file G4INCLParticle.hh.
Referenced by G4INCL::Store::loadParticles(), and G4INCL::Cluster::makeParticipant().
|
inlinevirtual |
Reimplemented in G4INCL::Cluster.
Definition at line 406 of file G4INCLParticle.hh.
Referenced by G4INCL::Cluster::makeProjectileSpectator(), G4INCL::StandardPropagationModel::shootCompositeAtrest(), and G4INCL::StandardPropagationModel::shootParticle().
|
inlinevirtual |
Reimplemented in G4INCL::Cluster.
Definition at line 402 of file G4INCLParticle.hh.
Referenced by G4INCL::Cluster::makeTargetSpectator().
|
static |
Definition at line 238 of file G4INCLParticle.cc.
Referenced by G4INCL::BinaryCollisionAvatar::getChannel(), G4INCL::ParticleList::getParticleListBias(), and G4INCL::ParticleList::getParticleListBiasVector().
|
static |
Definition at line 277 of file G4INCLParticle.cc.
Assignment operator.
Does not copy the particle ID.
Definition at line 181 of file G4INCLParticle.hh.
Referenced by G4INCL::Cluster::operator=().
|
inline |
Definition at line 1097 of file G4INCLParticle.hh.
Referenced by adjustMomentumFromEnergy(), G4INCL::Nucleus::emitInsideAnnihilationProducts(), G4INCL::Nucleus::emitInsideAntilambda(), G4INCL::Nucleus::emitInsideKaon(), G4INCL::Nucleus::emitInsideLambda(), G4INCL::Nucleus::emitInsidePions(), G4INCL::Nucleus::emitInsideStrangeParticles(), G4INCL::KinematicsUtils::getLocalEnergy(), G4INCL::StandardPropagationModel::getReflectionTime(), G4INCL::CrossSectionsMultiPions::piNTot(), G4INCL::ProjectileRemnant::removeParticle(), and G4INCL::Nucleus::restoreSrcPartner().
|
inline |
Definition at line 963 of file G4INCLParticle.hh.
Referenced by G4INCL::StandardPropagationModel::generateBinaryCollisionAvatar().
|
inline |
Definition at line 1003 of file G4INCLParticle.hh.
|
inlinevirtual |
Rotate the particle momentum.
| angle | the rotation angle |
| axis | a unit vector representing the rotation axis |
Reimplemented in G4INCL::Cluster.
Definition at line 1092 of file G4INCLParticle.hh.
Referenced by G4INCL::Cluster::rotateMomentum(), and rotatePositionAndMomentum().
|
inlinevirtual |
Rotate the particle position.
| angle | the rotation angle |
| axis | a unit vector representing the rotation axis |
Reimplemented in G4INCL::Cluster.
Definition at line 1083 of file G4INCLParticle.hh.
Referenced by G4INCL::Cluster::rotatePosition(), and rotatePositionAndMomentum().
|
inlinevirtual |
Rotate the particle position and momentum.
| angle | the rotation angle |
| axis | a unit vector representing the rotation axis |
Definition at line 1073 of file G4INCLParticle.hh.
|
inline |
Make the particle follow a strict r-p correlation.
Definition at line 1155 of file G4INCLParticle.hh.
Referenced by G4INCL::InteractionAvatar::bringParticleInside(), and G4INCL::SurfaceAvatar::postInteraction().
|
inline |
Make the particle not follow a strict r-p correlation.
Definition at line 1158 of file G4INCLParticle.hh.
|
inline |
Set the vector list of biased vertices on the particle path.
Definition at line 1188 of file G4INCLParticle.hh.
Referenced by G4INCL::Nucleus::decayOutgoingPionResonances(), G4INCL::Nucleus::decayOutgoingSigmaZero(), and G4INCL::SurfaceAvatar::postInteraction().
|
inline |
Definition at line 1008 of file G4INCLParticle.hh.
Referenced by G4INCL::Nucleus::decayOutgoingPionResonances(), G4INCL::Nucleus::emitInsideAntilambda(), G4INCL::Nucleus::emitInsideKaon(), G4INCL::Nucleus::emitInsideLambda(), G4INCL::Nucleus::emitInsidePions(), and G4INCL::Nucleus::emitInsideStrangeParticles().
|
inline |
Set the energy of the particle in MeV.
Definition at line 920 of file G4INCLParticle.hh.
Referenced by G4INCL::Nucleus::emitInsideAntilambda(), G4INCL::Nucleus::emitInsideKaon(), G4INCL::Nucleus::emitInsideLambda(), G4INCL::Nucleus::emitInsidePions(), G4INCL::Nucleus::emitInsideStrangeParticles(), G4INCL::EtaNElasticChannel::fillFinalState(), G4INCL::EtaNToPiNChannel::fillFinalState(), G4INCL::EtaOrOmegaNToLKChannel::fillFinalState(), G4INCL::EtaOrOmegaNToSKChannel::fillFinalState(), G4INCL::NSToNLChannel::fillFinalState(), G4INCL::NSToNSChannel::fillFinalState(), G4INCL::OmegaNElasticChannel::fillFinalState(), G4INCL::OmegaNToPiNChannel::fillFinalState(), G4INCL::CrossSections::interactionDistanceKbarN(), G4INCL::CrossSections::interactionDistanceKN(), G4INCL::CrossSections::interactionDistanceNbarN(), G4INCL::CrossSections::interactionDistancenbarN(), G4INCL::CrossSections::interactionDistanceNN(), G4INCL::CrossSections::interactionDistancePiN(), G4INCL::CrossSections::interactionDistanceYN(), G4INCL::StandardPropagationModel::shootParticle(), G4INCL::KinematicsUtils::transformToLocalEnergyFrame(), and G4INCL::Nucleus::useFusionKinematics().
|
inline |
|
inline |
Set the frozen particle momentum.
Definition at line 1032 of file G4INCLParticle.hh.
|
inline |
Definition at line 961 of file G4INCLParticle.hh.
|
static |
Definition at line 321 of file G4INCLParticle.cc.
|
inline |
Set the mass of the Particle to its table mass.
Definition at line 723 of file G4INCLParticle.hh.
Referenced by G4INCL::Cluster::Cluster(), G4INCL::Cluster::Cluster(), setType(), and G4INCL::StandardPropagationModel::shootParticle().
|
inline |
Set the mass of the particle in MeV/c^2.
Definition at line 912 of file G4INCLParticle.hh.
Referenced by G4INCL::Nucleus::computeRecoilKinematics(), G4INCL::NDeltaToDeltaLKChannel::fillFinalState(), G4INCL::PhaseSpaceKopylov::generate(), Particle(), setINCLMass(), setRealMass(), setTableMass(), and G4INCL::Nucleus::useFusionKinematics().
|
inlinevirtual |
Set the momentum vector.
Definition at line 942 of file G4INCLParticle.hh.
Referenced by G4INCL::Nucleus::decayOutgoingPionResonances(), G4INCL::Nucleus::decayOutgoingSigmaZero(), G4INCL::EtaNElasticChannel::fillFinalState(), G4INCL::EtaNToPiNChannel::fillFinalState(), G4INCL::EtaOrOmegaNToLKChannel::fillFinalState(), G4INCL::EtaOrOmegaNToSKChannel::fillFinalState(), G4INCL::NKbElasticChannel::fillFinalState(), G4INCL::NKbToLpiChannel::fillFinalState(), G4INCL::NKbToNKbChannel::fillFinalState(), G4INCL::NKbToSpiChannel::fillFinalState(), G4INCL::NKElasticChannel::fillFinalState(), G4INCL::NKToNKChannel::fillFinalState(), G4INCL::NSToNLChannel::fillFinalState(), G4INCL::NSToNSChannel::fillFinalState(), G4INCL::OmegaNElasticChannel::fillFinalState(), G4INCL::OmegaNToPiNChannel::fillFinalState(), G4INCL::StrangeAbsorbtionChannel::fillFinalState(), G4INCL::PhaseSpaceKopylov::generate(), G4INCL::Nucleus::restoreSrcPartner(), and G4INCL::Nucleus::useFusionKinematics().
|
inline |
Set the number of collisions undergone by the particle.
Definition at line 971 of file G4INCLParticle.hh.
|
inline |
Set the number of decays undergone by the particle.
Definition at line 980 of file G4INCLParticle.hh.
|
inline |
Definition at line 1201 of file G4INCLParticle.hh.
|
inline |
|
inline |
Mark the particle as out of its potential well.
This flag is used to control pions created outside their potential well in delta decay. The pion potential checks it and returns zero if it is true (necessary in order to correctly enforce energy conservation). The Nucleus::applyFinalState() method uses it to determine whether new avatars should be generated for the particle.
Definition at line 996 of file G4INCLParticle.hh.
|
inline |
Definition at line 1207 of file G4INCLParticle.hh.
|
inline |
Definition at line 1205 of file G4INCLParticle.hh.
|
inline |
Definition at line 382 of file G4INCLParticle.hh.
|
inline |
Set the particle bias.
Definition at line 1182 of file G4INCLParticle.hh.
Referenced by setBiasCollisionVector().
|
inlinevirtual |
Reimplemented in G4INCL::Cluster.
Definition at line 955 of file G4INCLParticle.hh.
Referenced by G4INCL::InteractionAvatar::bringParticleInside(), G4INCL::CoulombNone::bringToSurface(), G4INCL::CoulombNone::bringToSurfaceAbar(), G4INCL::ClusteringModelIntercomparison::getCluster(), G4INCL::ParticleSampler::sampleParticlesIntoList(), G4INCL::Cluster::setPosition(), and G4INCL::StandardPropagationModel::shootParticle().
|
inline |
Set the particle potential energy.
Definition at line 899 of file G4INCLParticle.hh.
Referenced by G4INCL::Nucleus::emitInsideAntilambda(), G4INCL::Nucleus::emitInsideKaon(), G4INCL::Nucleus::emitInsideLambda(), G4INCL::Nucleus::emitInsidePions(), G4INCL::Nucleus::emitInsideStrangeParticles(), G4INCL::Store::loadParticles(), and G4INCL::Nucleus::updatePotentialEnergy().
|
inline |
Set the mass of the Particle to its real mass.
Definition at line 717 of file G4INCLParticle.hh.
Referenced by G4INCL::ClusterDecay::decay().
|
inline |
Set and reset src partner.
Definition at line 1002 of file G4INCLParticle.hh.
|
inline |
Set the mass of the Particle to its table mass.
Definition at line 720 of file G4INCLParticle.hh.
Referenced by G4INCL::Nucleus::decayOutgoingPionResonances(), G4INCL::Nucleus::decayOutgoingSigmaZero(), G4INCL::Nucleus::emitInsideAntilambda(), G4INCL::Nucleus::emitInsideKaon(), G4INCL::Nucleus::emitInsideLambda(), G4INCL::Nucleus::emitInsidePions(), G4INCL::Nucleus::emitInsideStrangeParticles(), G4INCL::ProjectileRemnant::ProjectileRemnant(), G4INCL::ProjectileRemnant::reset(), and G4INCL::StandardPropagationModel::shootCompositeAtrest().
|
inline |
Definition at line 200 of file G4INCLParticle.hh.
Referenced by G4INCL::Cluster::Cluster(), G4INCL::Cluster::Cluster(), G4INCL::ClusterDecay::decay(), G4INCL::EtaNToPiNChannel::fillFinalState(), G4INCL::EtaNToPiPiNChannel::fillFinalState(), G4INCL::EtaOrOmegaNToLKChannel::fillFinalState(), G4INCL::EtaOrOmegaNToSKChannel::fillFinalState(), G4INCL::NDeltaToDeltaLKChannel::fillFinalState(), G4INCL::NKbToL2piChannel::fillFinalState(), G4INCL::NKbToLpiChannel::fillFinalState(), G4INCL::NKbToNKb2piChannel::fillFinalState(), G4INCL::NKbToNKbChannel::fillFinalState(), G4INCL::NKbToNKbpiChannel::fillFinalState(), G4INCL::NKbToS2piChannel::fillFinalState(), G4INCL::NKbToSpiChannel::fillFinalState(), G4INCL::NKToNK2piChannel::fillFinalState(), G4INCL::NKToNKChannel::fillFinalState(), G4INCL::NKToNKpiChannel::fillFinalState(), G4INCL::NNbarToNNbar2piChannel::fillFinalState(), G4INCL::NNToMissingStrangenessChannel::fillFinalState(), G4INCL::NpiToMissingStrangenessChannel::fillFinalState(), G4INCL::NSToNLChannel::fillFinalState(), G4INCL::NSToNSChannel::fillFinalState(), G4INCL::OmegaNToPiNChannel::fillFinalState(), G4INCL::OmegaNToPiPiNChannel::fillFinalState(), G4INCL::PiNToMultiPionsChannel::fillFinalState(), G4INCL::StrangeAbsorbtionChannel::fillFinalState(), Particle(), and Particle().
|
inline |
Set the uncorrelated momentum.
Definition at line 1152 of file G4INCLParticle.hh.
|
inlineprotected |
Helper method for the assignment operator.
Definition at line 131 of file G4INCLParticle.hh.
Referenced by operator=(), and G4INCL::Cluster::swap().
|
inline |
Unfreeze particle propagation.
Make the particle use theMomentum and theEnergy for propagation. Call this method to restore the normal propagation if the freezePropagation() method has been called.
Definition at line 1063 of file G4INCLParticle.hh.
|
protected |
Definition at line 1248 of file G4INCLParticle.hh.
Referenced by dump(), getID(), Particle(), Particle(), Particle(), Particle(), G4INCL::Cluster::print(), and print().
|
static |
Time ordered vector of all bias applied.
/!\ Caution /!\ methods Assotiated to G4VectorCache<T> are: Push_back(…), operator[], Begin(), End(), Clear(), Size() and Pop_back()
Definition at line 1225 of file G4INCLParticle.hh.
Referenced by FillINCLBiasVector(), getBiasFromVector(), getTotalBias(), G4INCL::INCL::processEvent(), and setINCLBiasVector().
|
protected |
Definition at line 1244 of file G4INCLParticle.hh.
Referenced by G4INCL::Cluster::addParticle(), getNumberOfCollisions(), incrementNumberOfCollisions(), Particle(), Particle(), Particle(), Particle(), G4INCL::ProjectileRemnant::reset(), setNumberOfCollisions(), swap(), and G4INCL::Cluster::updateClusterParameters().
|
protected |
Definition at line 1245 of file G4INCLParticle.hh.
Referenced by getNumberOfDecays(), incrementNumberOfDecays(), Particle(), Particle(), Particle(), Particle(), setNumberOfDecays(), and swap().
|
static |
Definition at line 1231 of file G4INCLParticle.hh.
Referenced by FillINCLBiasVector(), G4INCL::InteractionAvatar::postInteraction(), and G4INCL::INCL::processEvent().
|
protected |
Definition at line 1246 of file G4INCLParticle.hh.
Referenced by dump(), getSrcPair(), Particle(), Particle(), Particle(), Particle(), print(), resetSrcPartner(), setNumberOfSrcPair(), and swap().
|
protected |
Definition at line 1250 of file G4INCLParticle.hh.
Referenced by getReflectionMomentum(), Particle(), Particle(), Particle(), Particle(), rpCorrelate(), rpDecorrelate(), and swap().
|
protected |
Definition at line 1234 of file G4INCLParticle.hh.
Referenced by G4INCL::ProjectileRemnant::addAllDynamicalSpectators(), G4INCL::ProjectileRemnant::addMostDynamicalSpectators(), G4INCL::Cluster::addParticle(), G4INCL::Nucleus::applyFinalState(), G4INCL::Cluster::Cluster(), G4INCL::Cluster::Cluster(), G4INCL::Cluster::Cluster(), G4INCL::Nucleus::computeRecoilKinematics(), G4INCL::Nucleus::decayInsideDeltas(), G4INCL::Nucleus::decayInsideStrangeParticles(), G4INCL::Nucleus::decayMe(), G4INCL::Nucleus::emitInsideAnnihilationProducts(), G4INCL::Nucleus::emitInsideAntilambda(), G4INCL::Nucleus::emitInsideKaon(), G4INCL::Nucleus::emitInsideLambda(), G4INCL::Nucleus::emitInsidePions(), G4INCL::Nucleus::emitInsideStrangeParticles(), getA(), getEmissionQValueCorrection(), getEmissionQValueCorrection(), getINCLMass(), getRealMass(), G4INCL::Cluster::getSpecies(), getTableMass(), getTransferQValueCorrection(), getTransferQValueCorrection(), G4INCL::Cluster::initializeParticles(), G4INCL::Nucleus::insertParticle(), G4INCL::Cluster::internalBoostToCM(), G4INCL::Nucleus::Nucleus(), Particle(), Particle(), G4INCL::Cluster::print(), G4INCL::ProjectileRemnant::ProjectileRemnant(), G4INCL::ProjectileRemnant::removeParticle(), G4INCL::ProjectileRemnant::reset(), G4INCL::Cluster::setA(), setType(), swap(), and G4INCL::Cluster::updateClusterParameters().
|
protected |
Definition at line 1237 of file G4INCLParticle.hh.
Referenced by G4INCL::ProjectileRemnant::addAllDynamicalSpectators(), G4INCL::ProjectileRemnant::addMostDynamicalSpectators(), G4INCL::Cluster::addParticle(), adjustEnergyFromMomentum(), adjustMomentumFromEnergy(), boost(), boostVector(), dump(), getBeta(), getEnergy(), getInvariantMass(), getKineticEnergy(), G4INCL::Cluster::internalBoostToCM(), Particle(), Particle(), Particle(), Particle(), G4INCL::Cluster::print(), print(), G4INCL::ProjectileRemnant::removeParticle(), G4INCL::ProjectileRemnant::reset(), setEnergy(), swap(), thawPropagation(), G4INCL::Cluster::updateClusterParameters(), and G4INCL::Nucleus::useFusionKinematics().
|
protected |
Definition at line 1239 of file G4INCLParticle.hh.
Referenced by freezePropagation(), getFrozenEnergy(), Particle(), Particle(), Particle(), Particle(), setFrozenEnergy(), and swap().
|
protected |
Definition at line 1242 of file G4INCLParticle.hh.
Referenced by freezePropagation(), getFrozenMomentum(), Particle(), Particle(), Particle(), Particle(), rotateMomentum(), setFrozenMomentum(), and swap().
|
protected |
Definition at line 1240 of file G4INCLParticle.hh.
Referenced by G4INCL::ProjectileRemnant::addAllDynamicalSpectators(), G4INCL::ProjectileRemnant::addMostDynamicalSpectators(), G4INCL::Cluster::addParticle(), adjustEnergyFromMomentum(), adjustMomentumFromEnergy(), boost(), boostVector(), G4INCL::Nucleus::computeRecoilKinematics(), dump(), G4INCL::Cluster::freezeInternalMotion(), getAngularMomentum(), getBeta(), getInvariantMass(), getMomentum(), getReflectionMomentum(), G4INCL::Cluster::internalBoostToCM(), Particle(), Particle(), Particle(), Particle(), G4INCL::Cluster::print(), print(), G4INCL::ProjectileRemnant::removeParticle(), G4INCL::ProjectileRemnant::reset(), rotateMomentum(), setMomentum(), swap(), thawPropagation(), G4INCL::Cluster::updateClusterParameters(), and G4INCL::Nucleus::useFusionKinematics().
|
protected |
The number of Kaons inside the nucleus (update during the cascade).
Definition at line 1255 of file G4INCLParticle.hh.
Referenced by G4INCL::Nucleus::emitInsideKaon(), getNumberOfKaon(), Particle(), Particle(), Particle(), Particle(), and setNumberOfKaon().
|
protected |
Definition at line 1259 of file G4INCLParticle.hh.
Referenced by getParentResonanceID(), Particle(), Particle(), Particle(), Particle(), setParentResonanceID(), and swap().
|
protected |
Definition at line 1258 of file G4INCLParticle.hh.
Referenced by getParentResonancePDGCode(), Particle(), Particle(), Particle(), Particle(), setParentResonancePDGCode(), and swap().
|
protected |
Definition at line 1235 of file G4INCLParticle.hh.
Referenced by getParticipantType(), isParticipant(), isProjectileSpectator(), isTargetSpectator(), makeParticipant(), makeProjectileSpectator(), makeTargetSpectator(), Particle(), Particle(), Particle(), Particle(), setParticipantType(), and swap().
|
protected |
Definition at line 1253 of file G4INCLParticle.hh.
Referenced by getParticleBias(), Particle(), Particle(), Particle(), Particle(), setParticleBias(), and swap().
|
protected |
Definition at line 1243 of file G4INCLParticle.hh.
Referenced by G4INCL::Cluster::addParticle(), G4INCL::Cluster::boost(), G4INCL::Cluster::Cluster(), G4INCL::Nucleus::computeRecoilKinematics(), dump(), getAngularMomentum(), getCosRPAngle(), getLongitudinalPosition(), getPosition(), getTransversePosition(), G4INCL::Cluster::initializeParticles(), G4INCL::Nucleus::initializeParticles(), G4INCL::Cluster::internalBoostToCM(), lorentzContract(), Particle(), Particle(), Particle(), Particle(), G4INCL::Cluster::print(), print(), propagate(), G4INCL::ProjectileRemnant::reset(), rotatePosition(), G4INCL::Cluster::setPosition(), setPosition(), swap(), and G4INCL::Cluster::updateClusterParameters().
|
protected |
Definition at line 1247 of file G4INCLParticle.hh.
Referenced by G4INCL::Cluster::addParticle(), getPotentialEnergy(), Particle(), Particle(), Particle(), Particle(), print(), G4INCL::ProjectileRemnant::reset(), setPotentialEnergy(), swap(), and G4INCL::Cluster::updateClusterParameters().
|
protected |
Definition at line 1238 of file G4INCLParticle.hh.
Referenced by freezePropagation(), Particle(), Particle(), Particle(), Particle(), swap(), and thawPropagation().
|
protected |
Definition at line 1241 of file G4INCLParticle.hh.
Referenced by freezePropagation(), getCosRPAngle(), getLongitudinalPosition(), getPropagationVelocity(), Particle(), Particle(), Particle(), Particle(), swap(), and thawPropagation().
|
protected |
Definition at line 1234 of file G4INCLParticle.hh.
Referenced by G4INCL::ProjectileRemnant::addAllDynamicalSpectators(), G4INCL::ProjectileRemnant::addMostDynamicalSpectators(), G4INCL::Cluster::addParticle(), G4INCL::Nucleus::applyFinalState(), G4INCL::Cluster::Cluster(), G4INCL::Cluster::Cluster(), G4INCL::Nucleus::computeRecoilKinematics(), G4INCL::Nucleus::decayMe(), G4INCL::Nucleus::emitInsideAnnihilationProducts(), G4INCL::Nucleus::emitInsideAntilambda(), G4INCL::Nucleus::emitInsideKaon(), G4INCL::Nucleus::emitInsideLambda(), G4INCL::Nucleus::emitInsidePions(), G4INCL::Nucleus::emitInsideStrangeParticles(), getEmissionQValueCorrection(), getEmissionQValueCorrection(), getINCLMass(), getRealMass(), getS(), G4INCL::Cluster::getSpecies(), getTableMass(), getTransferQValueCorrection(), G4INCL::Nucleus::insertParticle(), G4INCL::Nucleus::Nucleus(), Particle(), Particle(), G4INCL::Cluster::print(), G4INCL::ProjectileRemnant::ProjectileRemnant(), G4INCL::ProjectileRemnant::removeParticle(), G4INCL::Cluster::setS(), setType(), swap(), and G4INCL::Cluster::updateClusterParameters().
|
protected |
Definition at line 1236 of file G4INCLParticle.hh.
Referenced by dump(), getINCLMass(), getRealMass(), getSpecies(), getTableMass(), getType(), isAntiKaon(), isAntiLambda(), isAntiNucleon(), isAntiSigma(), isAntiXi(), isCluster(), isDelta(), isEta(), isEtaPrime(), isKaon(), isLambda(), isNucleon(), isOmega(), isPhoton(), isPion(), isSigma(), isXi(), Particle(), Particle(), G4INCL::Cluster::print(), print(), setType(), and swap().
|
protected |
Definition at line 1234 of file G4INCLParticle.hh.
Referenced by G4INCL::ProjectileRemnant::addAllDynamicalSpectators(), G4INCL::ProjectileRemnant::addMostDynamicalSpectators(), G4INCL::Cluster::addParticle(), G4INCL::Nucleus::applyFinalState(), G4INCL::Cluster::Cluster(), G4INCL::Cluster::Cluster(), G4INCL::Nucleus::computeRecoilKinematics(), G4INCL::Nucleus::decayInsideDeltas(), G4INCL::Nucleus::decayInsideStrangeParticles(), G4INCL::Nucleus::decayMe(), G4INCL::Nucleus::emitInsideAnnihilationProducts(), G4INCL::Nucleus::emitInsideAntilambda(), G4INCL::Nucleus::emitInsideKaon(), G4INCL::Nucleus::emitInsideLambda(), G4INCL::Nucleus::emitInsidePions(), G4INCL::Nucleus::emitInsideStrangeParticles(), getEmissionQValueCorrection(), getEmissionQValueCorrection(), getINCLMass(), getRealMass(), G4INCL::Cluster::getSpecies(), getTableMass(), getTransferQValueCorrection(), getTransferQValueCorrection(), G4INCL::Nucleus::getTransmissionBarrier(), getZ(), G4INCL::Cluster::initializeParticles(), G4INCL::Nucleus::insertParticle(), G4INCL::Nucleus::Nucleus(), Particle(), Particle(), G4INCL::Cluster::print(), G4INCL::ProjectileRemnant::ProjectileRemnant(), G4INCL::ProjectileRemnant::removeParticle(), G4INCL::ProjectileRemnant::reset(), setType(), G4INCL::Cluster::setZ(), swap(), and G4INCL::Cluster::updateClusterParameters().
|
protected |
Definition at line 1251 of file G4INCLParticle.hh.
Referenced by getReflectionMomentum(), Particle(), Particle(), Particle(), Particle(), setUncorrelatedMomentum(), and swap().