Geant4-11
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Private Attributes | Static Private Attributes
G4INCL::Particle Class Reference

#include <G4INCLParticle.hh>

Inheritance diagram for G4INCL::Particle:
G4INCL::Cluster G4INCL::Nucleus G4INCL::ProjectileRemnant

Public Member Functions

G4double adjustEnergyFromMomentum ()
 Recompute the energy to match the momentum. More...
 
const ThreeVectoradjustMomentumFromEnergy ()
 Rescale the momentum to match the total energy. More...
 
void boost (const ThreeVector &aBoostVector)
 
ThreeVector boostVector () const
 
std::string dump () const
 
void freezePropagation ()
 Freeze particle propagation. More...
 
G4int getA () const
 Returns the baryon number. More...
 
virtual G4INCL::ThreeVector getAngularMomentum () const
 
G4double getBeta () const
 
std::vector< G4intgetBiasCollisionVector () const
 Get the vector list of biased vertices on the particle path. More...
 
G4double getCosRPAngle () const
 Get the cosine of the angle between position and momentum. More...
 
G4double getEmissionQValueCorrection (const G4int AParent, const G4int ZParent) const
 Computes correction on the emission Q-value. More...
 
G4double getEmissionQValueCorrection (const G4int AParent, const G4int ZParent, const G4int SParent) const
 Computes correction on the emission Q-value for hypernuclei. More...
 
G4double getEmissionTime ()
 
G4double getEnergy () const
 
G4double getFrozenEnergy () const
 Get the frozen particle momentum. More...
 
ThreeVector getFrozenMomentum () const
 Get the frozen particle momentum. More...
 
G4double getHelicity ()
 
long getID () const
 
G4double getINCLMass () const
 Get the INCL particle mass. More...
 
G4double getInvariantMass () const
 Get the the particle invariant mass. More...
 
G4double getKineticEnergy () const
 Get the particle kinetic energy. More...
 
ThreeVector getLongitudinalPosition () const
 Longitudinal component of the position w.r.t. the momentum. More...
 
G4double getMass () const
 Get the cached particle mass. More...
 
const G4INCL::ThreeVectorgetMomentum () const
 
G4int getNumberOfCollisions () const
 Return the number of collisions undergone by the particle. More...
 
G4int getNumberOfDecays () const
 Return the number of decays undergone by the particle. More...
 
G4int getNumberOfKaon () const
 Number of Kaon inside de nucleus. More...
 
ParticipantType getParticipantType () const
 
G4double getParticleBias () const
 Get the particle bias. More...
 
ParticleList const * getParticles () const
 
const G4INCL::ThreeVectorgetPosition () const
 
G4double getPotentialEnergy () const
 Get the particle potential energy. More...
 
ThreeVector getPropagationVelocity () const
 Get the propagation velocity of the particle. More...
 
G4double getRealMass () const
 Get the real particle mass. More...
 
G4double getReflectionMomentum () const
 Return the reflection momentum. More...
 
G4int getS () const
 Returns the strangeness number. More...
 
virtual G4INCL::ParticleSpecies getSpecies () const
 Get the particle species. More...
 
virtual G4double getTableMass () const
 Get the tabulated particle mass. 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 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. More...
 
ThreeVector getTransversePosition () const
 Transverse component of the position w.r.t. the momentum. More...
 
G4INCL::ParticleType getType () const
 
G4int getZ () const
 Returns the charge number. More...
 
void incrementNumberOfCollisions ()
 Increment the number of collisions undergone by the particle. More...
 
void incrementNumberOfDecays ()
 Increment the number of decays undergone by the particle. More...
 
G4bool isAntiKaon () const
 Is this an antiKaon? More...
 
G4bool isBaryon () const
 Is this a Baryon? More...
 
G4bool isCluster () const
 
G4bool isDelta () const
 Is it a Delta? More...
 
G4bool isEta () const
 Is this an eta? More...
 
G4bool isEtaPrime () const
 Is this an etaprime? More...
 
G4bool isHyperon () const
 Is this an Hyperon? More...
 
G4bool isKaon () const
 Is this a Kaon? More...
 
G4bool isLambda () const
 Is this a Lambda? More...
 
G4bool isMeson () const
 Is this a Meson? More...
 
G4bool isNucleon () const
 
G4bool isNucleonorLambda () const
 Is this a Nucleon or a Lambda? More...
 
G4bool isOmega () const
 Is this an omega? More...
 
G4bool isOutOfWell () const
 Check if the particle is out of its potential well. More...
 
G4bool isParticipant () const
 
G4bool isPhoton () const
 Is this a photon? More...
 
G4bool isPion () const
 Is this a pion? More...
 
G4bool isProjectileSpectator () const
 
G4bool isResonance () const
 Is it a resonance? More...
 
G4bool isSigma () const
 Is this a Sigma? More...
 
G4bool isStrange () const
 Is this an Strange? More...
 
G4bool isTargetSpectator () const
 
void lorentzContract (const ThreeVector &aBoostVector, const ThreeVector &refPos)
 Lorentz-contract the particle position around some center. More...
 
virtual void makeParticipant ()
 
virtual void makeProjectileSpectator ()
 
virtual void makeTargetSpectator ()
 
Particleoperator= (const Particle &rhs)
 Assignment operator. More...
 
 Particle ()
 
 Particle (const Particle &rhs)
 Copy constructor. More...
 
 Particle (ParticleType t, G4double energy, ThreeVector const &momentum, ThreeVector const &position)
 
 Particle (ParticleType t, ThreeVector const &momentum, ThreeVector const &position)
 
std::string print () const
 
void propagate (G4double step)
 
virtual void rotateMomentum (const G4double angle, const ThreeVector &axis)
 Rotate the particle momentum. More...
 
virtual void rotatePosition (const G4double angle, const ThreeVector &axis)
 Rotate the particle position. More...
 
virtual void rotatePositionAndMomentum (const G4double angle, const ThreeVector &axis)
 Rotate the particle position and 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...
 
void setBiasCollisionVector (std::vector< G4int > BiasCollisionVector)
 Set the vector list of biased vertices on the particle path. More...
 
void setEmissionTime (G4double t)
 
void setEnergy (G4double energy)
 
void setFrozenEnergy (const G4double energy)
 Set the frozen particle momentum. More...
 
void setFrozenMomentum (const ThreeVector &momentum)
 Set the frozen particle momentum. More...
 
void setHelicity (G4double h)
 
void setINCLMass ()
 Set the mass of the Particle to its table mass. More...
 
void setMass (G4double mass)
 
virtual void setMomentum (const G4INCL::ThreeVector &momentum)
 
void setNumberOfCollisions (G4int n)
 Set the number of collisions undergone by the particle. More...
 
void setNumberOfDecays (G4int n)
 Set the number of decays undergone by the particle. More...
 
void setNumberOfKaon (const G4int NK)
 
void setOutOfWell ()
 Mark the particle as out of its potential well. More...
 
void setParticipantType (ParticipantType const p)
 
void setParticleBias (G4double ParticleBias)
 Set the particle bias. More...
 
virtual void setPosition (const G4INCL::ThreeVector &position)
 
void setPotentialEnergy (G4double v)
 Set the particle potential energy. 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 setType (ParticleType t)
 
void setUncorrelatedMomentum (const G4double p)
 Set the uncorrelated momentum. More...
 
void thawPropagation ()
 Unfreeze particle propagation. More...
 
virtual ~Particle ()
 

Static Public Member Functions

static void FillINCLBiasVector (G4double newBias)
 
static G4double getBiasFromVector (std::vector< G4int > VectorBias)
 
static G4double getTotalBias ()
 General bias vector function. More...
 
static std::vector< G4intMergeVectorBias (Particle const *const p1, Particle const *const p2)
 
static std::vector< G4intMergeVectorBias (std::vector< G4int > p1, Particle const *const p2)
 
static void setINCLBiasVector (std::vector< G4double > NewVector)
 

Static Public Attributes

static std::vector< G4doubleINCLBiasVector
 Time ordered vector of all bias applied. More...
 
static G4ThreadLocal G4int nextBiasedCollisionID = 0
 

Protected Member Functions

void swap (Particle &rhs)
 Helper method for the assignment operator. More...
 

Protected Attributes

long ID
 
G4int nCollisions
 
G4int nDecays
 
G4bool rpCorrelated
 
G4int theA
 
G4double theEnergy
 
G4double theFrozenEnergy
 
G4INCL::ThreeVector theFrozenMomentum
 
G4INCL::ThreeVector theMomentum
 
G4int theNKaon
 The number of Kaons inside the nucleus (update during the cascade) More...
 
ParticipantType theParticipantType
 
G4double theParticleBias
 
G4INCL::ThreeVector thePosition
 
G4double thePotentialEnergy
 
G4doublethePropagationEnergy
 
G4INCL::ThreeVectorthePropagationMomentum
 
G4int theS
 
G4INCL::ParticleType theType
 
G4int theZ
 
G4double uncorrelatedMomentum
 

Private Attributes

G4double emissionTime
 
G4bool outOfWell
 
std::vector< G4inttheBiasCollisionVector
 Time ordered vector of all biased vertices on the particle path. More...
 
G4double theHelicity
 
G4double theMass
 

Static Private Attributes

static G4ThreadLocal long nextID = 1
 

Detailed Description

Definition at line 75 of file G4INCLParticle.hh.

Constructor & Destructor Documentation

◆ Particle() [1/4]

G4INCL::Particle::Particle ( )

Definition at line 59 of file G4INCLParticle.cc.

60 : theZ(0), theA(0), theS(0),
63 theEnergy(0.0),
66 theMomentum(ThreeVector(0.,0.,0.)),
69 thePosition(ThreeVector(0.,0.,0.)),
70 nCollisions(0),
71 nDecays(0),
73 rpCorrelated(false),
76 theNKaon(0),
77 theHelicity(0.0),
78 emissionTime(0.0),
79 outOfWell(false),
80 theMass(0.)
81 {
82 ID = nextID;
83 nextID++;
84 }
G4INCL::ThreeVector * thePropagationMomentum
G4INCL::ThreeVector theMomentum
G4double thePotentialEnergy
ParticipantType theParticipantType
G4INCL::ParticleType theType
G4int theNKaon
The number of Kaons inside the nucleus (update during the cascade)
static G4ThreadLocal long nextID
G4double uncorrelatedMomentum
G4double * thePropagationEnergy
G4INCL::ThreeVector thePosition
G4INCL::ThreeVector theFrozenMomentum

References ID, and nextID.

Referenced by G4INCL::Cluster::Cluster(), G4INCL::ProjectileRemnant::reset(), and G4INCL::ProjectileRemnant::storeComponents().

◆ Particle() [2/4]

G4INCL::Particle::Particle ( ParticleType  t,
G4double  energy,
ThreeVector const &  momentum,
ThreeVector const &  position 
)

Definition at line 86 of file G4INCLParticle.cc.

91 theMomentum(momentum),
95 nCollisions(0), nDecays(0),
97 rpCorrelated(false),
100 theNKaon(0),
101 theHelicity(0.0),
102 emissionTime(0.0), outOfWell(false)
103 {
105 ID = nextID;
106 nextID++;
107 if(theEnergy <= 0.0) {
108 INCL_WARN("Particle with energy " << theEnergy << " created." << '\n');
109 }
110 setType(t);
112 }
#define INCL_WARN(x)
void setMass(G4double mass)
G4double getInvariantMass() const
Get the the particle invariant mass.
void setType(ParticleType t)
G4double mag() const
G4double energy(const ThreeVector &p, const G4double m)

References getInvariantMass(), ID, INCL_WARN, nextID, setMass(), setType(), G4INCL::TargetSpectator, theEnergy, and theParticipantType.

◆ Particle() [3/4]

G4INCL::Particle::Particle ( ParticleType  t,
ThreeVector const &  momentum,
ThreeVector const &  position 
)

Definition at line 114 of file G4INCLParticle.cc.

117 theMomentum(momentum),
121 nCollisions(0), nDecays(0),
123 rpCorrelated(false),
125 theParticleBias(1.),
126 theNKaon(0),
127 theHelicity(0.0),
128 emissionTime(0.0), outOfWell(false)
129 {
131 ID = nextID;
132 nextID++;
133 setType(t);
134 if( isResonance() ) {
135 INCL_ERROR("Cannot create resonance without specifying its momentum four-vector." << '\n');
136 }
140 }
#define INCL_ERROR(x)
double G4double
Definition: G4Types.hh:83
G4bool isResonance() const
Is it a resonance?
G4double mag2() const

References G4INCL::KinematicsUtils::energy(), ID, INCL_ERROR, isResonance(), G4INCL::ThreeVector::mag2(), nextID, setType(), G4INCL::TargetSpectator, theEnergy, theFrozenEnergy, theMass, theMomentum, and theParticipantType.

◆ ~Particle()

virtual G4INCL::Particle::~Particle ( )
inlinevirtual

Definition at line 80 of file G4INCLParticle.hh.

80{}

◆ Particle() [4/4]

G4INCL::Particle::Particle ( const Particle rhs)
inline

Copy constructor.

Does not copy the particle ID.

Definition at line 86 of file G4INCLParticle.hh.

86 :
87 theZ(rhs.theZ),
88 theA(rhs.theA),
89 theS(rhs.theS),
90 theParticipantType(rhs.theParticipantType),
91 theType(rhs.theType),
92 theEnergy(rhs.theEnergy),
93 theFrozenEnergy(rhs.theFrozenEnergy),
94 theMomentum(rhs.theMomentum),
95 theFrozenMomentum(rhs.theFrozenMomentum),
96 thePosition(rhs.thePosition),
97 nCollisions(rhs.nCollisions),
98 nDecays(rhs.nDecays),
99 thePotentialEnergy(rhs.thePotentialEnergy),
100 rpCorrelated(rhs.rpCorrelated),
101 uncorrelatedMomentum(rhs.uncorrelatedMomentum),
102 theParticleBias(rhs.theParticleBias),
103 theNKaon(rhs.theNKaon),
104 theHelicity(rhs.theHelicity),
105 emissionTime(rhs.emissionTime),
106 outOfWell(rhs.outOfWell),
107 theMass(rhs.theMass)
108 {
109 if(rhs.thePropagationEnergy == &(rhs.theFrozenEnergy))
111 else
113 if(rhs.thePropagationMomentum == &(rhs.theFrozenMomentum))
115 else
117 // ID intentionally not copied
118 ID = nextID++;
119
120 theBiasCollisionVector = rhs.theBiasCollisionVector;
121 }
std::vector< G4int > theBiasCollisionVector
Time ordered vector of all biased vertices on the particle path.

References ID, nextID, theBiasCollisionVector, theEnergy, theFrozenEnergy, theFrozenMomentum, theMomentum, thePropagationEnergy, and thePropagationMomentum.

Member Function Documentation

◆ adjustEnergyFromMomentum()

G4double G4INCL::Particle::adjustEnergyFromMomentum ( )

◆ adjustMomentumFromEnergy()

const ThreeVector & G4INCL::Particle::adjustMomentumFromEnergy ( )

◆ boost()

void G4INCL::Particle::boost ( const ThreeVector aBoostVector)
inline

Boost the particle using a boost vector.

Example (go to the particle rest frame): particle->boost(particle->boostVector());

Definition at line 422 of file G4INCLParticle.hh.

422 {
423 const G4double beta2 = aBoostVector.mag2();
424 const G4double gamma = 1.0 / std::sqrt(1.0 - beta2);
425 const G4double bp = theMomentum.dot(aBoostVector);
426 const G4double alpha = (gamma*gamma)/(1.0 + gamma);
427
428 theMomentum = theMomentum + aBoostVector * (alpha * bp - gamma * theEnergy);
429 theEnergy = gamma * (theEnergy - bp);
430 }
static const G4double bp
static const G4double alpha
G4double dot(const ThreeVector &v) const

References alpha, bp, G4INCL::ThreeVector::dot(), G4INCL::ThreeVector::mag2(), theEnergy, and theMomentum.

Referenced by G4INCL::Cluster::boost(), G4INCL::Nucleus::computeOneNucleonRecoilKinematics(), G4INCL::Nucleus::decayOutgoingPionResonances(), G4INCL::Nucleus::decayOutgoingSigmaZero(), G4INCL::PhaseSpaceKopylov::generate(), G4INCL::InteractionAvatar::preInteraction(), G4INCL::ClusterDecay::anonymous_namespace{G4INCLClusterDecay.cc}::threeBodyDecay(), and G4INCL::ClusterDecay::anonymous_namespace{G4INCLClusterDecay.cc}::twoBodyDecay().

◆ boostVector()

ThreeVector G4INCL::Particle::boostVector ( ) const
inline

◆ dump()

std::string G4INCL::Particle::dump ( ) const
inline

Definition at line 968 of file G4INCLParticle.hh.

968 {
969 std::stringstream ss;
970 ss << "(particle " << ID << " ";
972 ss << '\n'
973 << thePosition.dump()
974 << '\n'
975 << theMomentum.dump()
976 << '\n'
977 << theEnergy << ")" << '\n';
978 return ss.str();
979 };
std::string dump() const
std::string getName(const ParticleType t)
Get the native INCL name of the particle.

References G4INCL::ThreeVector::dump(), G4INCL::ParticleTable::getName(), ID, theEnergy, theMomentum, thePosition, and theType.

Referenced by G4INCL::BinaryCollisionAvatar::dump(), G4INCL::DecayAvatar::dump(), G4INCL::ParticleEntryAvatar::dump(), and G4INCL::SurfaceAvatar::dump().

◆ FillINCLBiasVector()

void G4INCL::Particle::FillINCLBiasVector ( G4double  newBias)
static

Definition at line 202 of file G4INCLParticle.cc.

202 {
203// assert(G4int(Particle::INCLBiasVector.size())==nextBiasedCollisionID);
204 //assert(G4int(Particle::INCLBiasVector.Size())==nextBiasedCollisionID);
205// assert(std::fabs(newBias - 1.) > 1E-6);
206 Particle::INCLBiasVector.push_back(newBias);
207 //Particle::INCLBiasVector.Push_back(newBias);
209 }
static std::vector< G4double > INCLBiasVector
Time ordered vector of all bias applied.
static G4ThreadLocal G4int nextBiasedCollisionID

References INCLBiasVector, and nextBiasedCollisionID.

Referenced by G4INCL::InteractionAvatar::postInteraction().

◆ freezePropagation()

void G4INCL::Particle::freezePropagation ( )
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 908 of file G4INCLParticle.hh.

References theFrozenEnergy, theFrozenMomentum, thePropagationEnergy, and thePropagationMomentum.

◆ getA()

G4int G4INCL::Particle::getA ( ) const
inline

Returns the baryon number.

Definition at line 393 of file G4INCLParticle.hh.

393{ return theA; }

References theA.

Referenced by G4INCL::ProjectileRemnant::addDynamicalSpectator(), G4INCL::Cluster::addParticle(), G4INCL::ClusteringModelIntercomparison::clusterCanEscape(), G4INCL::Nucleus::computeOneNucleonRecoilKinematics(), G4INCL::INCL::continueCascade(), G4INCL::ClusterDecay::decay(), G4INCL::Nucleus::emitInsideLambda(), G4INCL::Nucleus::emitInsideStrangeParticles(), G4INCL::Nucleus::fillEventInfo(), G4INCL::ParticleEntryChannel::fillFinalState(), G4INCL::TransmissionChannel::fillFinalState(), G4INCL::Nucleus::finalizeProjectileRemnant(), G4INCL::ClusteringModelIntercomparison::findClusterStartingFrom(), G4INCL::PauliStandard::getBlockingProbability(), G4INCL::ClusteringModelIntercomparison::getCluster(), G4INCL::Nucleus::getConservationBalance(), G4INCL::SurfaceAvatar::getTransmissionProbability(), G4INCL::NuclearDensity::getTransmissionRadius(), G4INCL::TransmissionChannel::initializeKineticEnergyOutside(), G4INCL::Nucleus::insertParticle(), G4INCL::ClusterDecay::isStable(), G4INCL::ClusterDecay::anonymous_namespace{G4INCLClusterDecay.cc}::phaseSpaceDecay(), G4INCL::INCL::postCascade(), G4INCL::INCL::preCascade(), G4INCL::INCL::RecoilCMFunctor::RecoilCMFunctor(), G4INCL::INCL::RecoilFunctor::RecoilFunctor(), G4INCL::ClusterDecay::anonymous_namespace{G4INCLClusterDecay.cc}::recursiveDecay(), G4INCL::ProjectileRemnant::removeParticle(), G4INCL::INCL::RecoilCMFunctor::scaleParticleCMMomenta(), G4INCL::INCL::RecoilFunctor::scaleParticleEnergies(), G4INCL::InteractionAvatar::ViolationEMomentumFunctor::scaleParticleMomenta(), G4INCL::StandardPropagationModel::shootComposite(), G4INCL::StandardPropagationModel::shootParticle(), G4INCL::ClusterDecay::anonymous_namespace{G4INCLClusterDecay.cc}::threeBodyDecay(), and G4INCL::ClusterDecay::anonymous_namespace{G4INCLClusterDecay.cc}::twoBodyDecay().

◆ getAngularMomentum()

virtual G4INCL::ThreeVector G4INCL::Particle::getAngularMomentum ( ) const
inlinevirtual

Get the angular momentum w.r.t. the origin

Reimplemented in G4INCL::Cluster.

Definition at line 800 of file G4INCLParticle.hh.

801 {
803 };
ThreeVector vector(const ThreeVector &v) const

References theMomentum, thePosition, and G4INCL::ThreeVector::vector().

Referenced by G4INCL::Cluster::getAngularMomentum(), and G4INCL::StandardPropagationModel::shootParticle().

◆ getBeta()

G4double G4INCL::Particle::getBeta ( ) const
inline

Definition at line 401 of file G4INCLParticle.hh.

401 {
402 const G4double P = theMomentum.mag();
403 return P/theEnergy;
404 }
static double P[]

References G4INCL::ThreeVector::mag(), P, theEnergy, and theMomentum.

◆ getBiasCollisionVector()

std::vector< G4int > G4INCL::Particle::getBiasCollisionVector ( ) const
inline

◆ getBiasFromVector()

G4double G4INCL::Particle::getBiasFromVector ( std::vector< G4int VectorBias)
static

Definition at line 211 of file G4INCLParticle.cc.

211 {
212 if(VectorBias.empty()) return 1.;
213
214 G4double ParticleBias = 1.;
215
216 for(G4int i=0; i<G4int(VectorBias.size()); i++){
217 ParticleBias *= Particle::INCLBiasVector[G4int(VectorBias[i])];
218 }
219
220 return ParticleBias;
221 }
int G4int
Definition: G4Types.hh:85

References INCLBiasVector.

Referenced by G4INCL::BinaryCollisionAvatar::getChannel(), G4INCL::ParticleList::getParticleListBias(), and setBiasCollisionVector().

◆ getCosRPAngle()

G4double G4INCL::Particle::getCosRPAngle ( ) const
inline

Get the cosine of the angle between position and momentum.

Definition at line 1014 of file G4INCLParticle.hh.

1014 {
1016 if(norm>0.)
1017 return thePosition.dot(*thePropagationMomentum) / std::sqrt(norm);
1018 else
1019 return 1.;
1020 }

References G4INCL::ThreeVector::dot(), G4INCL::ThreeVector::mag2(), thePosition, and thePropagationMomentum.

Referenced by G4INCL::SurfaceAvatar::initializeRefractionVariables(), and G4INCL::ParticleEntryChannel::particleEnters().

◆ getEmissionQValueCorrection() [1/2]

G4double G4INCL::Particle::getEmissionQValueCorrection ( const G4int  AParent,
const G4int  ZParent 
) const
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.

Parameters
AParentthe mass number of the emitting nucleus
ZParentthe charge number of the emitting nucleus
Returns
the correction

Definition at line 602 of file G4INCLParticle.hh.

602 {
603 const G4int SParent = 0;
604 const G4int ADaughter = AParent - theA;
605 const G4int ZDaughter = ZParent - theZ;
606 const G4int SDaughter = 0;
607
608 // Note the minus sign here
609 G4double theQValue;
610 if(isCluster())
611 theQValue = -ParticleTable::getTableQValue(theA, theZ, theS, ADaughter, ZDaughter, SDaughter);
612 else {
613 const G4double massTableParent = ParticleTable::getTableMass(AParent,ZParent,SParent);
614 const G4double massTableDaughter = ParticleTable::getTableMass(ADaughter,ZDaughter,SDaughter);
615 const G4double massTableParticle = getTableMass();
616 theQValue = massTableParent - massTableDaughter - massTableParticle;
617 }
618
619 const G4double massINCLParent = ParticleTable::getINCLMass(AParent,ZParent,SParent);
620 const G4double massINCLDaughter = ParticleTable::getINCLMass(ADaughter,ZDaughter,SDaughter);
621 const G4double massINCLParticle = getINCLMass();
622
623 // The rhs corresponds to the INCL Q-value
624 return theQValue - (massINCLParent-massINCLDaughter-massINCLParticle);
625 }
G4bool isCluster() const
G4double getINCLMass() const
Get the INCL particle mass.
virtual G4double getTableMass() const
Get the tabulated particle mass.
G4ThreadLocal NuclearMassFn getTableMass
Static pointer to the mass function for nuclei.
G4double getTableQValue(const G4int A1, const G4int Z1, const G4int S1, const G4int A2, const G4int Z2, const G4int S2)
Get Q-value (in MeV/c^2)
G4double getINCLMass(const G4int A, const G4int Z, const G4int S)
Get INCL nuclear mass (in MeV/c^2)

References getINCLMass(), G4INCL::ParticleTable::getINCLMass(), getTableMass(), G4INCL::ParticleTable::getTableMass, G4INCL::ParticleTable::getTableQValue(), isCluster(), theA, theS, and theZ.

Referenced by G4INCL::Nucleus::emitInsideKaon(), G4INCL::Nucleus::emitInsideLambda(), G4INCL::Nucleus::emitInsidePions(), G4INCL::Nucleus::emitInsideStrangeParticles(), G4INCL::ParticleEntryChannel::fillFinalState(), G4INCL::TransmissionChannel::fillFinalState(), G4INCL::SurfaceAvatar::getTransmissionProbability(), and G4INCL::TransmissionChannel::initializeKineticEnergyOutside().

◆ getEmissionQValueCorrection() [2/2]

G4double G4INCL::Particle::getEmissionQValueCorrection ( const G4int  AParent,
const G4int  ZParent,
const G4int  SParent 
) const
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.

Parameters
AParentthe mass number of the emitting nucleus
ZParentthe charge number of the emitting nucleus
SParentthe strangess number of the emitting nucleus
Returns
the correction

Definition at line 678 of file G4INCLParticle.hh.

678 {
679 const G4int ADaughter = AParent - theA;
680 const G4int ZDaughter = ZParent - theZ;
681 const G4int SDaughter = SParent - theS;
682
683 // Note the minus sign here
684 G4double theQValue;
685 if(isCluster())
686 theQValue = -ParticleTable::getTableQValue(theA, theZ, theS, ADaughter, ZDaughter, SDaughter);
687 else {
688 const G4double massTableParent = ParticleTable::getTableMass(AParent,ZParent,SParent);
689 const G4double massTableDaughter = ParticleTable::getTableMass(ADaughter,ZDaughter,SDaughter);
690 const G4double massTableParticle = getTableMass();
691 theQValue = massTableParent - massTableDaughter - massTableParticle;
692 }
693
694 const G4double massINCLParent = ParticleTable::getINCLMass(AParent,ZParent,SParent);
695 const G4double massINCLDaughter = ParticleTable::getINCLMass(ADaughter,ZDaughter,SDaughter);
696 const G4double massINCLParticle = getINCLMass();
697
698 // The rhs corresponds to the INCL Q-value
699 return theQValue - (massINCLParent-massINCLDaughter-massINCLParticle);
700 }

References getINCLMass(), G4INCL::ParticleTable::getINCLMass(), getTableMass(), G4INCL::ParticleTable::getTableMass, G4INCL::ParticleTable::getTableQValue(), isCluster(), theA, theS, and theZ.

◆ getEmissionTime()

G4double G4INCL::Particle::getEmissionTime ( )
inline

◆ getEnergy()

G4double G4INCL::Particle::getEnergy ( ) const
inline

◆ getFrozenEnergy()

G4double G4INCL::Particle::getFrozenEnergy ( ) const
inline

Get the frozen particle momentum.

Definition at line 897 of file G4INCLParticle.hh.

897{ return theFrozenEnergy; }

References theFrozenEnergy.

◆ getFrozenMomentum()

ThreeVector G4INCL::Particle::getFrozenMomentum ( ) const
inline

Get the frozen particle momentum.

Definition at line 894 of file G4INCLParticle.hh.

894{ return theFrozenMomentum; }

References theFrozenMomentum.

◆ getHelicity()

G4double G4INCL::Particle::getHelicity ( )
inline

Definition at line 826 of file G4INCLParticle.hh.

826{ return theHelicity; };

References theHelicity.

Referenced by G4INCL::DeltaDecayChannel::sampleAngles().

◆ getID()

long G4INCL::Particle::getID ( ) const
inline

◆ getINCLMass()

G4double G4INCL::Particle::getINCLMass ( ) const
inline

Get the INCL particle mass.

Definition at line 454 of file G4INCLParticle.hh.

454 {
455 switch(theType) {
456 case Proton:
457 case Neutron:
458 case PiPlus:
459 case PiMinus:
460 case PiZero:
461 case Lambda:
462 case SigmaPlus:
463 case SigmaZero:
464 case SigmaMinus:
465 case KPlus:
466 case KZero:
467 case KZeroBar:
468 case KShort:
469 case KLong:
470 case KMinus:
471 case Eta:
472 case Omega:
473 case EtaPrime:
474 case Photon:
476 break;
477
478 case DeltaPlusPlus:
479 case DeltaPlus:
480 case DeltaZero:
481 case DeltaMinus:
482 return theMass;
483 break;
484
485 case Composite:
487 break;
488
489 default:
490 INCL_ERROR("Particle::getINCLMass: Unknown particle type." << '\n');
491 return 0.0;
492 break;
493 }
494 }

References G4INCL::Composite, G4INCL::DeltaMinus, G4INCL::DeltaPlus, G4INCL::DeltaPlusPlus, G4INCL::DeltaZero, G4INCL::Eta, G4INCL::EtaPrime, G4INCL::ParticleTable::getINCLMass(), INCL_ERROR, G4INCL::KLong, G4INCL::KMinus, G4INCL::KPlus, G4INCL::KShort, G4INCL::KZero, G4INCL::KZeroBar, G4INCL::Lambda, G4INCL::Neutron, G4INCL::Omega, G4INCL::Photon, G4INCL::PiMinus, G4INCL::PiPlus, G4INCL::PiZero, G4INCL::Proton, G4INCL::SigmaMinus, G4INCL::SigmaPlus, G4INCL::SigmaZero, theA, theMass, theS, theType, and theZ.

Referenced by G4INCL::ParticleEntryChannel::fillFinalState(), getEmissionQValueCorrection(), and setINCLMass().

◆ getInvariantMass()

G4double G4INCL::Particle::getInvariantMass ( ) const
inline

Get the the particle invariant mass.

Uses the relativistic invariant

\[ m = \sqrt{E^2 - {\vec p}^2}\]

Definition at line 748 of file G4INCLParticle.hh.

748 {
749 const G4double mass = std::pow(theEnergy, 2) - theMomentum.dot(theMomentum);
750 if(mass < 0.0) {
751 INCL_ERROR("E*E - p*p is negative." << '\n');
752 return 0.0;
753 } else {
754 return std::sqrt(mass);
755 }
756 };

References G4INCL::ThreeVector::dot(), INCL_ERROR, theEnergy, and theMomentum.

Referenced by G4INCL::Nucleus::finalizeProjectileRemnant(), and Particle().

◆ getKineticEnergy()

G4double G4INCL::Particle::getKineticEnergy ( ) const
inline

◆ getLongitudinalPosition()

ThreeVector G4INCL::Particle::getLongitudinalPosition ( ) const
inline

Longitudinal component of the position w.r.t. the momentum.

Definition at line 873 of file G4INCLParticle.hh.

References G4INCL::ThreeVector::dot(), G4INCL::ThreeVector::mag2(), thePosition, and thePropagationMomentum.

Referenced by G4INCL::CoulombNonRelativistic::coulombDeviation(), and getTransversePosition().

◆ getMass()

G4double G4INCL::Particle::getMass ( ) const
inline

Get the cached particle mass.

Definition at line 451 of file G4INCLParticle.hh.

451{ return theMass; }

References theMass.

Referenced by G4INCL::DeltaDecayChannel::computeDecayTime(), G4INCL::PionResonanceDecayChannel::computeDecayTime(), G4INCL::SigmaZeroDecayChannel::computeDecayTime(), G4INCL::Nucleus::computeOneNucleonRecoilKinematics(), G4INCL::Nucleus::emitInsideKaon(), G4INCL::Nucleus::emitInsideLambda(), G4INCL::Nucleus::emitInsidePions(), G4INCL::Nucleus::emitInsideStrangeParticles(), G4INCL::InteractionAvatar::enforceEnergyConservation(), G4INCL::CrossSectionsMultiPionsAndResonances::etaNToPiN(), G4INCL::DeltaDecayChannel::fillFinalState(), G4INCL::EtaNElasticChannel::fillFinalState(), G4INCL::EtaNToPiNChannel::fillFinalState(), G4INCL::NSToNLChannel::fillFinalState(), G4INCL::NSToNSChannel::fillFinalState(), G4INCL::OmegaNElasticChannel::fillFinalState(), G4INCL::OmegaNToPiNChannel::fillFinalState(), G4INCL::PionResonanceDecayChannel::fillFinalState(), G4INCL::RecombinationChannel::fillFinalState(), G4INCL::SigmaZeroDecayChannel::fillFinalState(), G4INCL::TransmissionChannel::fillFinalState(), G4INCL::Cluster::freezeInternalMotion(), G4INCL::PhaseSpaceKopylov::generate(), G4INCL::NuclearPotential::INuclearPotential::getFermiMomentum(), G4INCL::KinematicsUtils::getLocalEnergy(), G4INCL::SurfaceAvatar::getTransmissionProbability(), G4INCL::PhaseSpaceRauboldLynch::initialize(), G4INCL::TransmissionChannel::initializeKineticEnergyOutside(), G4INCL::CrossSections::interactionDistanceKbarN(), G4INCL::CrossSections::interactionDistanceKN(), 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::CrossSectionsMultiPionsAndResonances::omegaNToPiN(), G4INCL::ParticleEntryChannel::particleEnters(), G4INCL::TransmissionChannel::particleLeaves(), G4INCL::ClusterDecay::anonymous_namespace{G4INCLClusterDecay.cc}::phaseSpaceDecay(), G4INCL::CrossSectionsMultiPionsAndResonances::piMinuspToEtaN(), G4INCL::CrossSectionsMultiPionsAndResonances::piMinuspToOmegaN(), G4INCL::Cluster::print(), G4INCL::ProjectileRemnant::ProjectileRemnant(), G4INCL::StandardPropagationModel::shootParticle(), G4INCL::ClusterDecay::anonymous_namespace{G4INCLClusterDecay.cc}::threeBodyDecay(), and G4INCL::ClusterDecay::anonymous_namespace{G4INCLClusterDecay.cc}::twoBodyDecay().

◆ getMomentum()

const G4INCL::ThreeVector & G4INCL::Particle::getMomentum ( ) const
inline

Get the momentum vector.

Definition at line 794 of file G4INCLParticle.hh.

795 {
796 return theMomentum;
797 };

References theMomentum.

Referenced by G4INCL::Cluster::addParticle(), G4INCL::anonymous_namespace{G4INCLPhaseSpaceGenerator.cc}::bias(), G4INCL::ClusteringModelIntercomparison::clusterCanEscape(), G4INCL::SigmaZeroDecayChannel::computeDecayTime(), G4INCL::Nucleus::computeOneNucleonRecoilKinematics(), G4INCL::Nucleus::computeRecoilKinematics(), G4INCL::CoulombNonRelativistic::coulombDeviation(), G4INCL::Nucleus::decayOutgoingPionResonances(), G4INCL::Nucleus::decayOutgoingSigmaZero(), G4INCL::Nucleus::fillEventInfo(), G4INCL::DeltaProductionChannel::fillFinalState(), G4INCL::ElasticChannel::fillFinalState(), G4INCL::NDeltaEtaProductionChannel::fillFinalState(), G4INCL::NDeltaOmegaProductionChannel::fillFinalState(), G4INCL::NNToNSKpiChannel::fillFinalState(), G4INCL::NpiToLK2piChannel::fillFinalState(), G4INCL::NpiToLKpiChannel::fillFinalState(), G4INCL::NpiToNKKbChannel::fillFinalState(), G4INCL::NpiToSK2piChannel::fillFinalState(), G4INCL::NpiToSKpiChannel::fillFinalState(), G4INCL::ParticleEntryChannel::fillFinalState(), G4INCL::ReflectionChannel::fillFinalState(), G4INCL::StrangeAbsorbtionChannel::fillFinalState(), G4INCL::PhaseSpaceGenerator::generateBiased(), G4INCL::PauliStandard::getBlockingProbability(), G4INCL::ClusteringModelIntercomparison::getCluster(), G4INCL::Nucleus::getConservationBalance(), G4INCL::KinematicsUtils::getLocalEnergy(), G4INCL::ProjectileRemnant::getStoredMomentum(), G4INCL::SurfaceAvatar::getTransmissionProbability(), G4INCL::NKbElasticChannel::KaonMomentum(), G4INCL::NKbToLpiChannel::KaonMomentum(), G4INCL::NKbToNKbChannel::KaonMomentum(), G4INCL::NKbToSpiChannel::KaonMomentum(), G4INCL::KinematicsUtils::makeBoostVector(), G4INCL::INCL::makeCompoundNucleus(), G4INCL::KinematicsUtils::momentumInCM(), G4INCL::ParticleEntryChannel::particleEnters(), G4INCL::TransmissionChannel::particleLeaves(), G4INCL::BinaryCollisionAvatar::postInteraction(), G4INCL::InteractionAvatar::preInteraction(), G4INCL::INCL::RecoilCMFunctor::RecoilCMFunctor(), G4INCL::INCL::RecoilFunctor::RecoilFunctor(), G4INCL::ProjectileRemnant::removeParticle(), G4INCL::ParticleSampler::sampleParticlesIntoList(), G4INCL::StandardPropagationModel::shootComposite(), G4INCL::StandardPropagationModel::shootParticle(), and G4INCL::ClusterDecay::anonymous_namespace{G4INCLClusterDecay.cc}::twoBodyDecay().

◆ getNumberOfCollisions()

G4int G4INCL::Particle::getNumberOfCollisions ( ) const
inline

Return the number of collisions undergone by the particle.

Definition at line 834 of file G4INCLParticle.hh.

834{ return nCollisions; }

References nCollisions.

Referenced by G4INCL::Cluster::addParticle().

◆ getNumberOfDecays()

G4int G4INCL::Particle::getNumberOfDecays ( ) const
inline

Return the number of decays undergone by the particle.

Definition at line 843 of file G4INCLParticle.hh.

843{ return nDecays; }

References nDecays.

◆ getNumberOfKaon()

G4int G4INCL::Particle::getNumberOfKaon ( ) const
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 1053 of file G4INCLParticle.hh.

1053{ return theNKaon; };

References theNKaon.

Referenced by G4INCL::ParticleEntryChannel::fillFinalState(), G4INCL::SurfaceAvatar::getChannel(), and G4INCL::BinaryCollisionAvatar::postInteraction().

◆ getParticipantType()

ParticipantType G4INCL::Particle::getParticipantType ( ) const
inline

Definition at line 310 of file G4INCLParticle.hh.

310 {
311 return theParticipantType;
312 }

References theParticipantType.

Referenced by G4INCL::StandardPropagationModel::generateBinaryCollisionAvatar().

◆ getParticleBias()

G4double G4INCL::Particle::getParticleBias ( ) const
inline

Get the particle bias.

Definition at line 1032 of file G4INCLParticle.hh.

1032{ return theParticleBias; };

References theParticleBias.

◆ getParticles()

ParticleList const * G4INCL::Particle::getParticles ( ) const
inline

Return a NULL pointer

Definition at line 986 of file G4INCLParticle.hh.

986 {
987 INCL_WARN("Particle::getParticles() method was called on a Particle object" << '\n');
988 return 0;
989 }

References INCL_WARN.

◆ getPosition()

const G4INCL::ThreeVector & G4INCL::Particle::getPosition ( ) const
inline

Set the position vector.

Definition at line 816 of file G4INCLParticle.hh.

817 {
818 return thePosition;
819 };

References thePosition.

Referenced by G4INCL::Cluster::addParticle(), G4INCL::InteractionAvatar::bringParticleInside(), G4INCL::CoulombNone::bringToSurface(), G4INCL::ClusteringModelIntercomparison::clusterCanEscape(), G4INCL::DeltaDecayChannel::fillFinalState(), G4INCL::EtaNToPiPiNChannel::fillFinalState(), G4INCL::NDeltaEtaProductionChannel::fillFinalState(), G4INCL::NDeltaOmegaProductionChannel::fillFinalState(), G4INCL::NDeltaToDeltaSKChannel::fillFinalState(), G4INCL::NDeltaToNLKChannel::fillFinalState(), G4INCL::NDeltaToNNKKbChannel::fillFinalState(), G4INCL::NDeltaToNSKChannel::fillFinalState(), G4INCL::NKbToNKb2piChannel::fillFinalState(), G4INCL::NKToNK2piChannel::fillFinalState(), G4INCL::NNEtaToMultiPionsChannel::fillFinalState(), G4INCL::NNOmegaToMultiPionsChannel::fillFinalState(), G4INCL::NNToMissingStrangenessChannel::fillFinalState(), G4INCL::NNToMultiPionsChannel::fillFinalState(), G4INCL::NNToNLK2piChannel::fillFinalState(), G4INCL::NNToNLKChannel::fillFinalState(), G4INCL::NNToNLKpiChannel::fillFinalState(), G4INCL::NNToNNEtaChannel::fillFinalState(), G4INCL::NNToNNKKbChannel::fillFinalState(), G4INCL::NNToNNOmegaChannel::fillFinalState(), G4INCL::NNToNSK2piChannel::fillFinalState(), G4INCL::NNToNSKChannel::fillFinalState(), G4INCL::NNToNSKpiChannel::fillFinalState(), G4INCL::NpiToMissingStrangenessChannel::fillFinalState(), G4INCL::OmegaNToPiPiNChannel::fillFinalState(), G4INCL::PionResonanceDecayChannel::fillFinalState(), G4INCL::ReflectionChannel::fillFinalState(), G4INCL::SigmaZeroDecayChannel::fillFinalState(), G4INCL::StandardPropagationModel::generateBinaryCollisionAvatar(), G4INCL::PauliStandard::getBlockingProbability(), G4INCL::BinaryCollisionAvatar::getChannel(), G4INCL::ClusteringModelIntercomparison::getCluster(), G4INCL::KinematicsUtils::getLocalEnergy(), G4INCL::StandardPropagationModel::getReflectionTime(), G4INCL::StandardPropagationModel::getTime(), G4INCL::ParticleEntryChannel::particleEnters(), G4INCL::TransmissionChannel::particleLeaves(), G4INCL::ClusterDecay::anonymous_namespace{G4INCLClusterDecay.cc}::phaseSpaceDecay(), G4INCL::BinaryCollisionAvatar::postInteraction(), G4INCL::ParticleSampler::sampleParticlesIntoList(), G4INCL::ClusterDecay::anonymous_namespace{G4INCLClusterDecay.cc}::threeBodyDecay(), and G4INCL::ClusterDecay::anonymous_namespace{G4INCLClusterDecay.cc}::twoBodyDecay().

◆ getPotentialEnergy()

G4double G4INCL::Particle::getPotentialEnergy ( ) const
inline

◆ getPropagationVelocity()

ThreeVector G4INCL::Particle::getPropagationVelocity ( ) const
inline

Get the propagation velocity of the particle.

Definition at line 900 of file G4INCLParticle.hh.

900{ return (*thePropagationMomentum)/(*thePropagationEnergy); }

References thePropagationMomentum.

Referenced by G4INCL::CoulombNone::bringToSurface(), G4INCL::StandardPropagationModel::getReflectionTime(), and G4INCL::StandardPropagationModel::getTime().

◆ getRealMass()

G4double G4INCL::Particle::getRealMass ( ) const
inline

Get the real particle mass.

Definition at line 540 of file G4INCLParticle.hh.

540 {
541 switch(theType) {
542 case Proton:
543 case Neutron:
544 case PiPlus:
545 case PiMinus:
546 case PiZero:
547 case Lambda:
548 case SigmaPlus:
549 case SigmaZero:
550 case SigmaMinus:
551 case KPlus:
552 case KZero:
553 case KZeroBar:
554 case KShort:
555 case KLong:
556 case KMinus:
557 case Eta:
558 case Omega:
559 case EtaPrime:
560 case Photon:
562 break;
563
564 case DeltaPlusPlus:
565 case DeltaPlus:
566 case DeltaZero:
567 case DeltaMinus:
568 return theMass;
569 break;
570
571 case Composite:
573 break;
574
575 default:
576 INCL_ERROR("Particle::getRealMass: Unknown particle type." << '\n');
577 return 0.0;
578 break;
579 }
580 }
G4double getRealMass(const G4INCL::ParticleType t)
Get particle mass (in MeV/c^2)

References G4INCL::Composite, G4INCL::DeltaMinus, G4INCL::DeltaPlus, G4INCL::DeltaPlusPlus, G4INCL::DeltaZero, G4INCL::Eta, G4INCL::EtaPrime, G4INCL::ParticleTable::getRealMass(), INCL_ERROR, G4INCL::KLong, G4INCL::KMinus, G4INCL::KPlus, G4INCL::KShort, G4INCL::KZero, G4INCL::KZeroBar, G4INCL::Lambda, G4INCL::Neutron, G4INCL::Omega, G4INCL::Photon, G4INCL::PiMinus, G4INCL::PiPlus, G4INCL::PiZero, G4INCL::Proton, G4INCL::SigmaMinus, G4INCL::SigmaPlus, G4INCL::SigmaZero, theA, theMass, theS, theType, and theZ.

Referenced by G4INCL::Cluster::getTableMass(), and setRealMass().

◆ getReflectionMomentum()

G4double G4INCL::Particle::getReflectionMomentum ( ) const
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 997 of file G4INCLParticle.hh.

997 {
998 if(rpCorrelated)
999 return theMomentum.mag();
1000 else
1001 return uncorrelatedMomentum;
1002 }

References G4INCL::ThreeVector::mag(), rpCorrelated, theMomentum, and uncorrelatedMomentum.

Referenced by G4INCL::KinematicsUtils::getLocalEnergy(), and G4INCL::Nucleus::getSurfaceRadius().

◆ getS()

G4int G4INCL::Particle::getS ( ) const
inline

Returns the strangeness number.

Definition at line 399 of file G4INCLParticle.hh.

399{ return theS; }

References theS.

Referenced by G4INCL::ProjectileRemnant::addDynamicalSpectator(), G4INCL::Cluster::addParticle(), G4INCL::ClusteringModelIntercomparison::clusterCanEscape(), G4INCL::Nucleus::computeOneNucleonRecoilKinematics(), G4INCL::ClusterDecay::decay(), G4INCL::Nucleus::emitInsideKaon(), G4INCL::Nucleus::emitInsideLambda(), G4INCL::Nucleus::emitInsideStrangeParticles(), G4INCL::Nucleus::fillEventInfo(), G4INCL::ParticleEntryChannel::fillFinalState(), G4INCL::TransmissionChannel::fillFinalState(), G4INCL::Nucleus::finalizeProjectileRemnant(), G4INCL::Nucleus::getConservationBalance(), G4INCL::SurfaceAvatar::getTransmissionProbability(), G4INCL::TransmissionChannel::initializeKineticEnergyOutside(), G4INCL::Nucleus::insertParticle(), G4INCL::ClusterDecay::isStable(), G4INCL::ClusterDecay::anonymous_namespace{G4INCLClusterDecay.cc}::phaseSpaceDecay(), G4INCL::INCL::preCascade(), G4INCL::ClusterDecay::anonymous_namespace{G4INCLClusterDecay.cc}::recursiveDecay(), G4INCL::ProjectileRemnant::removeParticle(), G4INCL::INCL::RecoilCMFunctor::scaleParticleCMMomenta(), G4INCL::INCL::RecoilFunctor::scaleParticleEnergies(), G4INCL::StandardPropagationModel::shootComposite(), G4INCL::StandardPropagationModel::shootParticle(), G4INCL::ClusterDecay::anonymous_namespace{G4INCLClusterDecay.cc}::threeBodyDecay(), and G4INCL::ClusterDecay::anonymous_namespace{G4INCLClusterDecay.cc}::twoBodyDecay().

◆ getSpecies()

virtual G4INCL::ParticleSpecies G4INCL::Particle::getSpecies ( ) const
inlinevirtual

◆ getTableMass()

virtual G4double G4INCL::Particle::getTableMass ( ) const
inlinevirtual

Get the tabulated particle mass.

Reimplemented in G4INCL::Cluster.

Definition at line 497 of file G4INCLParticle.hh.

497 {
498 switch(theType) {
499 case Proton:
500 case Neutron:
501 case PiPlus:
502 case PiMinus:
503 case PiZero:
504 case Lambda:
505 case SigmaPlus:
506 case SigmaZero:
507 case SigmaMinus:
508 case KPlus:
509 case KZero:
510 case KZeroBar:
511 case KShort:
512 case KLong:
513 case KMinus:
514 case Eta:
515 case Omega:
516 case EtaPrime:
517 case Photon:
519 break;
520
521 case DeltaPlusPlus:
522 case DeltaPlus:
523 case DeltaZero:
524 case DeltaMinus:
525 return theMass;
526 break;
527
528 case Composite:
530 break;
531
532 default:
533 INCL_ERROR("Particle::getTableMass: Unknown particle type." << '\n');
534 return 0.0;
535 break;
536 }
537 }
G4ThreadLocal ParticleMassFn getTableParticleMass
Static pointer to the mass function for particles.

References G4INCL::Composite, G4INCL::DeltaMinus, G4INCL::DeltaPlus, G4INCL::DeltaPlusPlus, G4INCL::DeltaZero, G4INCL::Eta, G4INCL::EtaPrime, G4INCL::ParticleTable::getTableMass, G4INCL::ParticleTable::getTableParticleMass, INCL_ERROR, G4INCL::KLong, G4INCL::KMinus, G4INCL::KPlus, G4INCL::KShort, G4INCL::KZero, G4INCL::KZeroBar, G4INCL::Lambda, G4INCL::Neutron, G4INCL::Omega, G4INCL::Photon, G4INCL::PiMinus, G4INCL::PiPlus, G4INCL::PiZero, G4INCL::Proton, G4INCL::SigmaMinus, G4INCL::SigmaPlus, G4INCL::SigmaZero, theA, theMass, theS, theType, and theZ.

Referenced by G4INCL::Nucleus::decayOutgoingPionResonances(), G4INCL::Nucleus::decayOutgoingSigmaZero(), G4INCL::ParticleEntryChannel::fillFinalState(), G4INCL::TransmissionChannel::fillFinalState(), getEmissionQValueCorrection(), getTransferQValueCorrection(), and setTableMass().

◆ getTotalBias()

G4double G4INCL::Particle::getTotalBias ( )
static

General bias vector function.

Definition at line 300 of file G4INCLParticle.cc.

300 {
301 G4double TotalBias = 1.;
302 for(G4int i=0; i<G4int(INCLBiasVector.size());i++) TotalBias *= Particle::INCLBiasVector[i];
303 return TotalBias;
304 }

References INCLBiasVector.

Referenced by G4INCL::Nucleus::emitInsideKaon(), G4INCL::Nucleus::emitInsideLambda(), G4INCL::Nucleus::emitInsidePions(), G4INCL::Nucleus::emitInsideStrangeParticles(), G4INCL::INCL::postCascade(), and G4INCL::EventInfo::remnantToParticle().

◆ getTransferQValueCorrection() [1/2]

G4double G4INCL::Particle::getTransferQValueCorrection ( const G4int  AFrom,
const G4int  ZFrom,
const G4int  ATo,
const G4int  ZTo 
) const
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.

Parameters
AFromthe mass number of the donating nucleus
ZFromthe charge number of the donating nucleus
ATothe mass number of the receiving nucleus
ZTothe charge number of the receiving nucleus
Returns
the correction

Definition at line 642 of file G4INCLParticle.hh.

642 {
643 const G4int SFrom = 0;
644 const G4int STo = 0;
645 const G4int AFromDaughter = AFrom - theA;
646 const G4int ZFromDaughter = ZFrom - theZ;
647 const G4int SFromDaughter = 0;
648 const G4int AToDaughter = ATo + theA;
649 const G4int ZToDaughter = ZTo + theZ;
650 const G4int SToDaughter = 0;
651 const G4double theQValue = ParticleTable::getTableQValue(AToDaughter,ZToDaughter,SToDaughter,AFromDaughter,ZFromDaughter,SFromDaughter,AFrom,ZFrom,SFrom);
652
653 const G4double massINCLTo = ParticleTable::getINCLMass(ATo,ZTo,STo);
654 const G4double massINCLToDaughter = ParticleTable::getINCLMass(AToDaughter,ZToDaughter,SToDaughter);
655 /* Note that here we have to use the table mass in the INCL Q-value. We
656 * cannot use theMass, because at this stage the particle is probably
657 * still off-shell; and we cannot use getINCLMass(), because it leads to
658 * violations of global energy conservation.
659 */
660 const G4double massINCLParticle = getTableMass();
661
662 // The rhs corresponds to the INCL Q-value for particle absorption
663 return theQValue - (massINCLToDaughter-massINCLTo-massINCLParticle);
664 }

References G4INCL::ParticleTable::getINCLMass(), getTableMass(), G4INCL::ParticleTable::getTableQValue(), theA, and theZ.

◆ getTransferQValueCorrection() [2/2]

G4double G4INCL::Particle::getTransferQValueCorrection ( const G4int  AFrom,
const G4int  ZFrom,
const G4int  SFrom,
const G4int  ATo,
const G4int  ZTo,
const G4int  STo 
) const
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.

Parameters
AFromthe mass number of the donating nucleus
ZFromthe charge number of the donating nucleus
SFromthe strangess number of the donating nucleus
ATothe mass number of the receiving nucleus
ZTothe charge number of the receiving nucleus
STothe strangess number of the receiving nucleus
Returns
the correction

Definition at line 719 of file G4INCLParticle.hh.

719 {
720 const G4int AFromDaughter = AFrom - theA;
721 const G4int ZFromDaughter = ZFrom - theZ;
722 const G4int SFromDaughter = SFrom - theS;
723 const G4int AToDaughter = ATo + theA;
724 const G4int ZToDaughter = ZTo + theZ;
725 const G4int SToDaughter = STo + theS;
726 const G4double theQValue = ParticleTable::getTableQValue(AToDaughter,ZToDaughter,SFromDaughter,AFromDaughter,ZFromDaughter,SToDaughter,AFrom,ZFrom,SFrom);
727
728 const G4double massINCLTo = ParticleTable::getINCLMass(ATo,ZTo,STo);
729 const G4double massINCLToDaughter = ParticleTable::getINCLMass(AToDaughter,ZToDaughter,SToDaughter);
730 /* Note that here we have to use the table mass in the INCL Q-value. We
731 * cannot use theMass, because at this stage the particle is probably
732 * still off-shell; and we cannot use getINCLMass(), because it leads to
733 * violations of global energy conservation.
734 */
735 const G4double massINCLParticle = getTableMass();
736
737 // The rhs corresponds to the INCL Q-value for particle absorption
738 return theQValue - (massINCLToDaughter-massINCLTo-massINCLParticle);
739 }

References G4INCL::ParticleTable::getINCLMass(), getTableMass(), G4INCL::ParticleTable::getTableQValue(), theA, theS, and theZ.

◆ getTransversePosition()

ThreeVector G4INCL::Particle::getTransversePosition ( ) const
inline

Transverse component of the position w.r.t. the momentum.

Definition at line 868 of file G4INCLParticle.hh.

868 {
870 }
ThreeVector getLongitudinalPosition() const
Longitudinal component of the position w.r.t. the momentum.

References getLongitudinalPosition(), and thePosition.

Referenced by G4INCL::CoulombNonRelativistic::coulombDeviation(), G4INCL::StandardPropagationModel::shootComposite(), and G4INCL::StandardPropagationModel::shootParticle().

◆ getType()

G4INCL::ParticleType G4INCL::Particle::getType ( ) const
inline

Get the particle type.

See also
G4INCL::ParticleType

Definition at line 178 of file G4INCLParticle.hh.

178 {
179 return theType;
180 };

References theType.

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::DeltaDecayChannel::fillFinalState(), G4INCL::DeltaProductionChannel::fillFinalState(), G4INCL::ElasticChannel::fillFinalState(), G4INCL::NDeltaEtaProductionChannel::fillFinalState(), G4INCL::NDeltaOmegaProductionChannel::fillFinalState(), G4INCL::NDeltaToDeltaLKChannel::fillFinalState(), G4INCL::NDeltaToDeltaSKChannel::fillFinalState(), G4INCL::NDeltaToNLKChannel::fillFinalState(), G4INCL::NDeltaToNNKKbChannel::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::NNEtaToMultiPionsChannel::fillFinalState(), G4INCL::NNOmegaToMultiPionsChannel::fillFinalState(), G4INCL::NNToMissingStrangenessChannel::fillFinalState(), G4INCL::NNToMultiPionsChannel::fillFinalState(), G4INCL::NNToNLK2piChannel::fillFinalState(), G4INCL::NNToNLKChannel::fillFinalState(), G4INCL::NNToNLKpiChannel::fillFinalState(), G4INCL::NNToNNEtaChannel::fillFinalState(), G4INCL::NNToNNKKbChannel::fillFinalState(), G4INCL::NNToNNOmegaChannel::fillFinalState(), G4INCL::NNToNSK2piChannel::fillFinalState(), G4INCL::NNToNSKChannel::fillFinalState(), G4INCL::NNToNSKpiChannel::fillFinalState(), G4INCL::NpiToMissingStrangenessChannel::fillFinalState(), G4INCL::NSToNLChannel::fillFinalState(), G4INCL::NSToNSChannel::fillFinalState(), G4INCL::PionResonanceDecayChannel::fillFinalState(), G4INCL::RecombinationChannel::fillFinalState(), G4INCL::PauliStandard::getBlockingProbability(), G4INCL::DecayAvatar::getChannel(), G4INCL::NuclearPotential::INuclearPotential::getFermiEnergy(), G4INCL::NuclearPotential::INuclearPotential::getFermiMomentum(), G4INCL::KinematicsUtils::getLocalEnergy(), G4INCL::NuclearPotential::INuclearPotential::getSeparationEnergy(), G4INCL::Nucleus::getSurfaceRadius(), 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::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(), and G4INCL::CrossSectionsMultiPions::piNTot().

◆ getZ()

G4int G4INCL::Particle::getZ ( ) const
inline

Returns the charge number.

Definition at line 396 of file G4INCLParticle.hh.

396{ return theZ; }

References theZ.

Referenced by G4INCL::ProjectileRemnant::addDynamicalSpectator(), G4INCL::Cluster::addParticle(), G4INCL::CoulombNonRelativistic::bringToSurface(), G4INCL::Nucleus::computeOneNucleonRecoilKinematics(), G4INCL::ClusterDecay::decay(), G4INCL::CoulombNonRelativistic::distortOut(), G4INCL::Nucleus::emitInsideKaon(), G4INCL::Nucleus::emitInsidePions(), G4INCL::Nucleus::emitInsideStrangeParticles(), G4INCL::Nucleus::fillEventInfo(), G4INCL::ParticleEntryChannel::fillFinalState(), G4INCL::TransmissionChannel::fillFinalState(), G4INCL::Nucleus::finalizeProjectileRemnant(), G4INCL::ClusteringModelIntercomparison::findClusterStartingFrom(), G4INCL::PauliStandard::getBlockingProbability(), G4INCL::ClusteringModelIntercomparison::getCluster(), G4INCL::Nucleus::getConservationBalance(), G4INCL::Nucleus::getTransmissionBarrier(), G4INCL::SurfaceAvatar::getTransmissionProbability(), G4INCL::NuclearDensity::getTransmissionRadius(), G4INCL::TransmissionChannel::initializeKineticEnergyOutside(), G4INCL::Nucleus::insertParticle(), G4INCL::ClusterDecay::isStable(), G4INCL::ClusterDecay::anonymous_namespace{G4INCLClusterDecay.cc}::phaseSpaceDecay(), G4INCL::INCL::preCascade(), G4INCL::ClusterDecay::anonymous_namespace{G4INCLClusterDecay.cc}::recursiveDecay(), G4INCL::ProjectileRemnant::removeParticle(), G4INCL::INCL::RecoilCMFunctor::scaleParticleCMMomenta(), G4INCL::INCL::RecoilFunctor::scaleParticleEnergies(), G4INCL::StandardPropagationModel::shootComposite(), G4INCL::StandardPropagationModel::shootParticle(), G4INCL::ClusterDecay::anonymous_namespace{G4INCLClusterDecay.cc}::threeBodyDecay(), and G4INCL::ClusterDecay::anonymous_namespace{G4INCLClusterDecay.cc}::twoBodyDecay().

◆ incrementNumberOfCollisions()

void G4INCL::Particle::incrementNumberOfCollisions ( )
inline

Increment the number of collisions undergone by the particle.

Definition at line 840 of file G4INCLParticle.hh.

840{ nCollisions++; }

References nCollisions.

◆ incrementNumberOfDecays()

void G4INCL::Particle::incrementNumberOfDecays ( )
inline

Increment the number of decays undergone by the particle.

Definition at line 849 of file G4INCLParticle.hh.

849{ nDecays++; }

References nDecays.

◆ isAntiKaon()

G4bool G4INCL::Particle::isAntiKaon ( ) const
inline

◆ isBaryon()

G4bool G4INCL::Particle::isBaryon ( ) const
inline

Is this a Baryon?

Definition at line 387 of file G4INCLParticle.hh.

387{ return (isNucleon() || isResonance() || isHyperon()); }
G4bool isHyperon() const
Is this an Hyperon?
G4bool isNucleon() const

References isHyperon(), isNucleon(), and isResonance().

Referenced by G4INCL::CDPP::processOneParticle().

◆ isCluster()

G4bool G4INCL::Particle::isCluster ( ) const
inline

Definition at line 883 of file G4INCLParticle.hh.

883 {
884 return (theType == Composite);
885 }

References G4INCL::Composite, and theType.

Referenced by getEmissionQValueCorrection(), and G4INCL::SurfaceAvatar::postInteraction().

◆ isDelta()

G4bool G4INCL::Particle::isDelta ( ) const
inline

◆ isEta()

G4bool G4INCL::Particle::isEta ( ) const
inline

◆ isEtaPrime()

G4bool G4INCL::Particle::isEtaPrime ( ) const
inline

Is this an etaprime?

Definition at line 352 of file G4INCLParticle.hh.

352{ return (theType == EtaPrime); }

References G4INCL::EtaPrime, and theType.

Referenced by isMeson(), G4INCL::CrossSectionsMultiPionsAndResonances::total(), and G4INCL::CrossSectionsStrangeness::total().

◆ isHyperon()

G4bool G4INCL::Particle::isHyperon ( ) const
inline

Is this an Hyperon?

Definition at line 381 of file G4INCLParticle.hh.

381{ return (isLambda() || isSigma()); }
G4bool isLambda() const
Is this a Lambda?
G4bool isSigma() const
Is this a Sigma?

References isLambda(), and isSigma().

Referenced by G4INCL::CrossSectionsStrangeness::elastic(), isBaryon(), isStrange(), and G4INCL::CrossSectionsStrangeness::NYelastic().

◆ isKaon()

G4bool G4INCL::Particle::isKaon ( ) const
inline

◆ isLambda()

G4bool G4INCL::Particle::isLambda ( ) const
inline

◆ isMeson()

G4bool G4INCL::Particle::isMeson ( ) const
inline

Is this a Meson?

Definition at line 384 of file G4INCLParticle.hh.

384{ return (isPion() || isKaon() || isAntiKaon() || isEta() || isEtaPrime() || isOmega()); }
G4bool isEtaPrime() const
Is this an etaprime?
G4bool isOmega() const
Is this an omega?
G4bool isEta() const
Is this an eta?
G4bool isPion() const
Is this a pion?
G4bool isAntiKaon() const
Is this an antiKaon?
G4bool isKaon() const
Is this a Kaon?

References isAntiKaon(), isEta(), isEtaPrime(), isKaon(), isOmega(), and isPion().

Referenced by G4INCL::StandardPropagationModel::generateBinaryCollisionAvatar(), G4INCL::InteractionAvatar::preInteractionLocalEnergy(), G4INCL::CDPP::processOneParticle(), and G4INCL::StandardPropagationModel::shootParticle().

◆ isNucleon()

G4bool G4INCL::Particle::isNucleon ( ) const
inline

Is this a nucleon?

Definition at line 303 of file G4INCLParticle.hh.

303 {
305 return true;
306 else
307 return false;
308 };

References G4INCL::Neutron, G4INCL::Proton, and theType.

Referenced by G4INCL::NuclearPotential::NuclearPotentialEnergyIsospin::computePotentialEnergy(), G4INCL::NuclearPotential::NuclearPotentialEnergyIsospinSmooth::computePotentialEnergy(), G4INCL::CrossSectionsINCL46::elastic(), G4INCL::CrossSectionsMultiPions::elastic(), G4INCL::CrossSectionsMultiPionsAndResonances::elastic(), G4INCL::CrossSectionsStrangeness::elastic(), G4INCL::CrossSectionsTruncatedMultiPions::elastic(), G4INCL::CrossSectionsINCL46::elasticNNLegacy(), G4INCL::EtaNElasticChannel::fillFinalState(), G4INCL::EtaNToPiNChannel::fillFinalState(), G4INCL::EtaNToPiPiNChannel::fillFinalState(), G4INCL::NDeltaToNNKKbChannel::fillFinalState(), G4INCL::NKbElasticChannel::fillFinalState(), G4INCL::NKbToL2piChannel::fillFinalState(), G4INCL::NKbToLpiChannel::fillFinalState(), G4INCL::NKbToNKb2piChannel::fillFinalState(), G4INCL::NKbToNKbChannel::fillFinalState(), G4INCL::NKbToNKbpiChannel::fillFinalState(), G4INCL::NKbToS2piChannel::fillFinalState(), G4INCL::NKbToSpiChannel::fillFinalState(), G4INCL::NKElasticChannel::fillFinalState(), G4INCL::NKToNK2piChannel::fillFinalState(), G4INCL::NKToNKChannel::fillFinalState(), G4INCL::NKToNKpiChannel::fillFinalState(), G4INCL::NLToNSChannel::fillFinalState(), G4INCL::NpiToLK2piChannel::fillFinalState(), G4INCL::NpiToLKChannel::fillFinalState(), G4INCL::NpiToLKpiChannel::fillFinalState(), G4INCL::NpiToNKKbChannel::fillFinalState(), G4INCL::NpiToSK2piChannel::fillFinalState(), G4INCL::NpiToSKChannel::fillFinalState(), G4INCL::NpiToSKpiChannel::fillFinalState(), G4INCL::NSToNLChannel::fillFinalState(), G4INCL::NSToNSChannel::fillFinalState(), G4INCL::NYElasticChannel::fillFinalState(), G4INCL::OmegaNElasticChannel::fillFinalState(), G4INCL::OmegaNToPiNChannel::fillFinalState(), G4INCL::OmegaNToPiPiNChannel::fillFinalState(), G4INCL::PiNElasticChannel::fillFinalState(), G4INCL::PiNToDeltaChannel::fillFinalState(), G4INCL::PiNToEtaChannel::fillFinalState(), G4INCL::PiNToMultiPionsChannel::fillFinalState(), G4INCL::PiNToOmegaChannel::fillFinalState(), G4INCL::StrangeAbsorbtionChannel::fillFinalState(), G4INCL::StandardPropagationModel::generateBinaryCollisionAvatar(), G4INCL::BinaryCollisionAvatar::getChannel(), G4INCL::DecayAvatar::getChannel(), G4INCL::Nucleus::getSurfaceRadius(), G4INCL::Nucleus::insertParticle(), isBaryon(), isNucleonorLambda(), G4INCL::CrossSectionsMultiPions::NNElastic(), G4INCL::CrossSectionsMultiPions::NNTot(), 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::CrossSectionsINCL46::total(), G4INCL::CrossSectionsMultiPions::total(), G4INCL::CrossSectionsMultiPionsAndResonances::total(), and G4INCL::CrossSectionsStrangeness::total().

◆ isNucleonorLambda()

G4bool G4INCL::Particle::isNucleonorLambda ( ) const
inline

Is this a Nucleon or a Lambda?

Definition at line 378 of file G4INCLParticle.hh.

378{ return (isNucleon() || isLambda()); }

References isLambda(), and isNucleon().

Referenced by G4INCL::ParticleEntryChannel::fillFinalState(), and G4INCL::SurfaceAvatar::getChannel().

◆ isOmega()

G4bool G4INCL::Particle::isOmega ( ) const
inline

◆ isOutOfWell()

G4bool G4INCL::Particle::isOutOfWell ( ) const
inline

◆ isParticipant()

G4bool G4INCL::Particle::isParticipant ( ) const
inline

◆ isPhoton()

G4bool G4INCL::Particle::isPhoton ( ) const
inline

Is this a photon?

Definition at line 355 of file G4INCLParticle.hh.

355{ return (theType == Photon); }

References G4INCL::Photon, and theType.

◆ isPion()

G4bool G4INCL::Particle::isPion ( ) const
inline

◆ isProjectileSpectator()

G4bool G4INCL::Particle::isProjectileSpectator ( ) const
inline

Definition at line 326 of file G4INCLParticle.hh.

326 {
328 }

References G4INCL::ProjectileSpectator, and theParticipantType.

Referenced by G4INCL::SurfaceAvatar::getChannel().

◆ isResonance()

G4bool G4INCL::Particle::isResonance ( ) const
inline

◆ isSigma()

G4bool G4INCL::Particle::isSigma ( ) const
inline

◆ isStrange()

G4bool G4INCL::Particle::isStrange ( ) const
inline

Is this an Strange?

Definition at line 390 of file G4INCLParticle.hh.

390{ return (isKaon() || isAntiKaon() || isHyperon()); }

References isAntiKaon(), isHyperon(), and isKaon().

◆ isTargetSpectator()

G4bool G4INCL::Particle::isTargetSpectator ( ) const
inline

◆ lorentzContract()

void G4INCL::Particle::lorentzContract ( const ThreeVector aBoostVector,
const ThreeVector refPos 
)
inline

Lorentz-contract the particle position around some center.

Apply Lorentz contraction to the position component along the direction of the boost vector.

Parameters
aBoostVectorthe boost vector (velocity) [c]
refPosthe reference position

Definition at line 440 of file G4INCLParticle.hh.

440 {
441 const G4double beta2 = aBoostVector.mag2();
442 const G4double gamma = 1.0 / std::sqrt(1.0 - beta2);
443 const ThreeVector theRelativePosition = thePosition - refPos;
444 const ThreeVector transversePosition = theRelativePosition - aBoostVector * (theRelativePosition.dot(aBoostVector) / aBoostVector.mag2());
445 const ThreeVector longitudinalPosition = theRelativePosition - transversePosition;
446
447 thePosition = refPos + transversePosition + longitudinalPosition / gamma;
448 }

References G4INCL::ThreeVector::dot(), G4INCL::ThreeVector::mag2(), and thePosition.

◆ makeParticipant()

virtual void G4INCL::Particle::makeParticipant ( )
inlinevirtual

◆ makeProjectileSpectator()

virtual void G4INCL::Particle::makeProjectileSpectator ( )
inlinevirtual

◆ makeTargetSpectator()

virtual void G4INCL::Particle::makeTargetSpectator ( )
inlinevirtual

Reimplemented in G4INCL::Cluster.

Definition at line 334 of file G4INCLParticle.hh.

334 {
336 }

References G4INCL::TargetSpectator, and theParticipantType.

Referenced by G4INCL::Cluster::makeTargetSpectator().

◆ MergeVectorBias() [1/2]

std::vector< G4int > G4INCL::Particle::MergeVectorBias ( Particle const *const  p1,
Particle const *const  p2 
)
static

Definition at line 223 of file G4INCLParticle.cc.

223 {
224 std::vector<G4int> MergedVectorBias;
225 std::vector<G4int> VectorBias1 = p1->getBiasCollisionVector();
226 std::vector<G4int> VectorBias2 = p2->getBiasCollisionVector();
227 G4int i = 0;
228 G4int j = 0;
229 if(VectorBias1.size()==0 && VectorBias2.size()==0) return MergedVectorBias;
230 else if(VectorBias1.size()==0) return VectorBias2;
231 else if(VectorBias2.size()==0) return VectorBias1;
232
233 while(i < G4int(VectorBias1.size()) || j < G4int(VectorBias2.size())){
234 if(VectorBias1[i]==VectorBias2[j]){
235 MergedVectorBias.push_back(VectorBias1[i]);
236 i++;
237 j++;
238 if(i == G4int(VectorBias1.size())){
239 for(;j<G4int(VectorBias2.size());j++) MergedVectorBias.push_back(VectorBias2[j]);
240 }
241 else if(j == G4int(VectorBias2.size())){
242 for(;i<G4int(VectorBias1.size());i++) MergedVectorBias.push_back(VectorBias1[i]);
243 }
244 } else if(VectorBias1[i]<VectorBias2[j]){
245 MergedVectorBias.push_back(VectorBias1[i]);
246 i++;
247 if(i == G4int(VectorBias1.size())){
248 for(;j<G4int(VectorBias2.size());j++) MergedVectorBias.push_back(VectorBias2[j]);
249 }
250 }
251 else {
252 MergedVectorBias.push_back(VectorBias2[j]);
253 j++;
254 if(j == G4int(VectorBias2.size())){
255 for(;i<G4int(VectorBias1.size());i++) MergedVectorBias.push_back(VectorBias1[i]);
256 }
257 }
258 }
259 return MergedVectorBias;
260 }

References getBiasCollisionVector().

Referenced by G4INCL::BinaryCollisionAvatar::getChannel(), G4INCL::ParticleList::getParticleListBias(), and G4INCL::ParticleList::getParticleListBiasVector().

◆ MergeVectorBias() [2/2]

std::vector< G4int > G4INCL::Particle::MergeVectorBias ( std::vector< G4int p1,
Particle const *const  p2 
)
static

Definition at line 262 of file G4INCLParticle.cc.

262 {
263 std::vector<G4int> MergedVectorBias;
264 std::vector<G4int> VectorBias = p2->getBiasCollisionVector();
265 G4int i = 0;
266 G4int j = 0;
267 if(p1.size()==0 && VectorBias.size()==0) return MergedVectorBias;
268 else if(p1.size()==0) return VectorBias;
269 else if(VectorBias.size()==0) return p1;
270
271 while(i < G4int(p1.size()) || j < G4int(VectorBias.size())){
272 if(p1[i]==VectorBias[j]){
273 MergedVectorBias.push_back(p1[i]);
274 i++;
275 j++;
276 if(i == G4int(p1.size())){
277 for(;j<G4int(VectorBias.size());j++) MergedVectorBias.push_back(VectorBias[j]);
278 }
279 else if(j == G4int(VectorBias.size())){
280 for(;i<G4int(p1.size());i++) MergedVectorBias.push_back(p1[i]);
281 }
282 } else if(p1[i]<VectorBias[j]){
283 MergedVectorBias.push_back(p1[i]);
284 i++;
285 if(i == G4int(p1.size())){
286 for(;j<G4int(VectorBias.size());j++) MergedVectorBias.push_back(VectorBias[j]);
287 }
288 }
289 else {
290 MergedVectorBias.push_back(VectorBias[j]);
291 j++;
292 if(j == G4int(VectorBias.size())){
293 for(;i<G4int(p1.size());i++) MergedVectorBias.push_back(p1[i]);
294 }
295 }
296 }
297 return MergedVectorBias;
298 }

References getBiasCollisionVector().

◆ operator=()

Particle & G4INCL::Particle::operator= ( const Particle rhs)
inline

Assignment operator.

Does not copy the particle ID.

Definition at line 168 of file G4INCLParticle.hh.

168 {
169 Particle temporaryParticle(rhs);
170 swap(temporaryParticle);
171 return *this;
172 }
void swap(Particle &rhs)
Helper method for the assignment operator.

References swap().

Referenced by G4INCL::Cluster::operator=().

◆ print()

std::string G4INCL::Particle::print ( ) const
inline

◆ propagate()

void G4INCL::Particle::propagate ( G4double  step)
inline

Definition at line 829 of file G4INCLParticle.hh.

829 {
830 thePosition += ((*thePropagationMomentum)*(step/(*thePropagationEnergy)));
831 };

References thePosition.

Referenced by G4INCL::StandardPropagationModel::generateBinaryCollisionAvatar().

◆ rotateMomentum()

virtual void G4INCL::Particle::rotateMomentum ( const G4double  angle,
const ThreeVector axis 
)
inlinevirtual

Rotate the particle momentum.

Parameters
anglethe rotation angle
axisa unit vector representing the rotation axis

Reimplemented in G4INCL::Cluster.

Definition at line 948 of file G4INCLParticle.hh.

948 {
949 theMomentum.rotate(angle, axis);
951 }
static const G4double angle[DIMMOTT]
void rotate(const G4double angle, const ThreeVector &axis)
Rotate the vector by a given angle around a given axis.

References angle, G4INCL::ThreeVector::rotate(), theFrozenMomentum, and theMomentum.

Referenced by G4INCL::Cluster::rotateMomentum(), and rotatePositionAndMomentum().

◆ rotatePosition()

virtual void G4INCL::Particle::rotatePosition ( const G4double  angle,
const ThreeVector axis 
)
inlinevirtual

Rotate the particle position.

Parameters
anglethe rotation angle
axisa unit vector representing the rotation axis

Reimplemented in G4INCL::Cluster.

Definition at line 939 of file G4INCLParticle.hh.

939 {
940 thePosition.rotate(angle, axis);
941 }

References angle, G4INCL::ThreeVector::rotate(), and thePosition.

Referenced by G4INCL::Cluster::rotatePosition(), and rotatePositionAndMomentum().

◆ rotatePositionAndMomentum()

virtual void G4INCL::Particle::rotatePositionAndMomentum ( const G4double  angle,
const ThreeVector axis 
)
inlinevirtual

Rotate the particle position and momentum.

Parameters
anglethe rotation angle
axisa unit vector representing the rotation axis

Definition at line 929 of file G4INCLParticle.hh.

929 {
930 rotatePosition(angle, axis);
931 rotateMomentum(angle, axis);
932 }
virtual void rotateMomentum(const G4double angle, const ThreeVector &axis)
Rotate the particle momentum.
virtual void rotatePosition(const G4double angle, const ThreeVector &axis)
Rotate the particle position.

References angle, rotateMomentum(), and rotatePosition().

Referenced by G4INCL::CoulombNonRelativistic::coulombDeviation().

◆ rpCorrelate()

void G4INCL::Particle::rpCorrelate ( )
inline

Make the particle follow a strict r-p correlation.

Definition at line 1008 of file G4INCLParticle.hh.

1008{ rpCorrelated = true; }

References rpCorrelated.

Referenced by G4INCL::InteractionAvatar::bringParticleInside(), G4INCL::ParticleEntryAvatar::postInteraction(), and G4INCL::SurfaceAvatar::postInteraction().

◆ rpDecorrelate()

void G4INCL::Particle::rpDecorrelate ( )
inline

Make the particle not follow a strict r-p correlation.

Definition at line 1011 of file G4INCLParticle.hh.

1011{ rpCorrelated = false; }

References rpCorrelated.

◆ setBiasCollisionVector()

void G4INCL::Particle::setBiasCollisionVector ( std::vector< G4int BiasCollisionVector)
inline

Set the vector list of biased vertices on the particle path.

Definition at line 1041 of file G4INCLParticle.hh.

1041 {
1042 this->theBiasCollisionVector = BiasCollisionVector;
1043 this->setParticleBias(Particle::getBiasFromVector(BiasCollisionVector));
1044 }
void setParticleBias(G4double ParticleBias)
Set the particle bias.
static G4double getBiasFromVector(std::vector< G4int > VectorBias)

References getBiasFromVector(), setParticleBias(), and theBiasCollisionVector.

Referenced by G4INCL::Nucleus::decayOutgoingPionResonances(), G4INCL::Nucleus::decayOutgoingSigmaZero(), and G4INCL::SurfaceAvatar::postInteraction().

◆ setEmissionTime()

void G4INCL::Particle::setEmissionTime ( G4double  t)
inline

◆ setEnergy()

void G4INCL::Particle::setEnergy ( G4double  energy)
inline

◆ setFrozenEnergy()

void G4INCL::Particle::setFrozenEnergy ( const G4double  energy)
inline

Set the frozen particle momentum.

Definition at line 891 of file G4INCLParticle.hh.

References G4INCL::KinematicsUtils::energy(), and theFrozenEnergy.

◆ setFrozenMomentum()

void G4INCL::Particle::setFrozenMomentum ( const ThreeVector momentum)
inline

Set the frozen particle momentum.

Definition at line 888 of file G4INCLParticle.hh.

888{ theFrozenMomentum = momentum; }

References theFrozenMomentum.

◆ setHelicity()

void G4INCL::Particle::setHelicity ( G4double  h)
inline

◆ setINCLBiasVector()

void G4INCL::Particle::setINCLBiasVector ( std::vector< G4double NewVector)
static

Definition at line 306 of file G4INCLParticle.cc.

306 {
307 Particle::INCLBiasVector = NewVector;
308 }

References INCLBiasVector.

◆ setINCLMass()

void G4INCL::Particle::setINCLMass ( )
inline

◆ setMass()

void G4INCL::Particle::setMass ( G4double  mass)
inline

◆ setMomentum()

virtual void G4INCL::Particle::setMomentum ( const G4INCL::ThreeVector momentum)
inlinevirtual

Set the momentum vector.

Definition at line 808 of file G4INCLParticle.hh.

809 {
810 this->theMomentum = momentum;
811 };

References theMomentum.

Referenced by G4INCL::Nucleus::computeOneNucleonRecoilKinematics(), G4INCL::Nucleus::decayOutgoingPionResonances(), G4INCL::Nucleus::decayOutgoingSigmaZero(), G4INCL::DeltaDecayChannel::fillFinalState(), G4INCL::DeltaProductionChannel::fillFinalState(), G4INCL::ElasticChannel::fillFinalState(), G4INCL::EtaNElasticChannel::fillFinalState(), G4INCL::EtaNToPiNChannel::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::PionResonanceDecayChannel::fillFinalState(), G4INCL::RecombinationChannel::fillFinalState(), G4INCL::ReflectionChannel::fillFinalState(), G4INCL::SigmaZeroDecayChannel::fillFinalState(), G4INCL::StrangeAbsorbtionChannel::fillFinalState(), G4INCL::PhaseSpaceKopylov::generate(), G4INCL::PhaseSpaceRauboldLynch::generateEvent(), G4INCL::INCL::makeCompoundNucleus(), G4INCL::ParticleEntryChannel::particleEnters(), G4INCL::TransmissionChannel::particleLeaves(), G4INCL::ClusterDecay::anonymous_namespace{G4INCLClusterDecay.cc}::phaseSpaceDecay(), G4INCL::INCL::RecoilCMFunctor::scaleParticleCMMomenta(), G4INCL::INCL::RecoilFunctor::scaleParticleEnergies(), G4INCL::ClusterDecay::anonymous_namespace{G4INCLClusterDecay.cc}::threeBodyDecay(), G4INCL::ClusterDecay::anonymous_namespace{G4INCLClusterDecay.cc}::twoBodyDecay(), and G4INCL::Nucleus::useFusionKinematics().

◆ setNumberOfCollisions()

void G4INCL::Particle::setNumberOfCollisions ( G4int  n)
inline

Set the number of collisions undergone by the particle.

Definition at line 837 of file G4INCLParticle.hh.

References CLHEP::detail::n, and nCollisions.

◆ setNumberOfDecays()

void G4INCL::Particle::setNumberOfDecays ( G4int  n)
inline

◆ setNumberOfKaon()

void G4INCL::Particle::setNumberOfKaon ( const G4int  NK)
inline

◆ setOutOfWell()

void G4INCL::Particle::setOutOfWell ( )
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 859 of file G4INCLParticle.hh.

859{ outOfWell = true; }

References outOfWell.

◆ setParticipantType()

void G4INCL::Particle::setParticipantType ( ParticipantType const  p)
inline

Definition at line 314 of file G4INCLParticle.hh.

314 {
316 }

References theParticipantType.

◆ setParticleBias()

void G4INCL::Particle::setParticleBias ( G4double  ParticleBias)
inline

Set the particle bias.

Definition at line 1035 of file G4INCLParticle.hh.

1035{ this->theParticleBias = ParticleBias; }

References theParticleBias.

Referenced by setBiasCollisionVector().

◆ setPosition()

virtual void G4INCL::Particle::setPosition ( const G4INCL::ThreeVector position)
inlinevirtual

◆ setPotentialEnergy()

void G4INCL::Particle::setPotentialEnergy ( G4double  v)
inline

◆ setRealMass()

void G4INCL::Particle::setRealMass ( )
inline

◆ setTableMass()

void G4INCL::Particle::setTableMass ( )
inline

◆ setType()

void G4INCL::Particle::setType ( ParticleType  t)
inline

Definition at line 187 of file G4INCLParticle.hh.

187 {
188 theType = t;
189 switch(theType)
190 {
191 case DeltaPlusPlus:
192 theA = 1;
193 theZ = 2;
194 theS = 0;
195 break;
196 case Proton:
197 case DeltaPlus:
198 theA = 1;
199 theZ = 1;
200 theS = 0;
201 break;
202 case Neutron:
203 case DeltaZero:
204 theA = 1;
205 theZ = 0;
206 theS = 0;
207 break;
208 case DeltaMinus:
209 theA = 1;
210 theZ = -1;
211 theS = 0;
212 break;
213 case PiPlus:
214 theA = 0;
215 theZ = 1;
216 theS = 0;
217 break;
218 case PiZero:
219 case Eta:
220 case Omega:
221 case EtaPrime:
222 case Photon:
223 theA = 0;
224 theZ = 0;
225 theS = 0;
226 break;
227 case PiMinus:
228 theA = 0;
229 theZ = -1;
230 theS = 0;
231 break;
232 case Lambda:
233 theA = 1;
234 theZ = 0;
235 theS = -1;
236 break;
237 case SigmaPlus:
238 theA = 1;
239 theZ = 1;
240 theS = -1;
241 break;
242 case SigmaZero:
243 theA = 1;
244 theZ = 0;
245 theS = -1;
246 break;
247 case SigmaMinus:
248 theA = 1;
249 theZ = -1;
250 theS = -1;
251 break;
252 case KPlus:
253 theA = 0;
254 theZ = 1;
255 theS = 1;
256 break;
257 case KZero:
258 theA = 0;
259 theZ = 0;
260 theS = 1;
261 break;
262 case KZeroBar:
263 theA = 0;
264 theZ = 0;
265 theS = -1;
266 break;
267 case KShort:
268 theA = 0;
269 theZ = 0;
270// theS should not be defined
271 break;
272 case KLong:
273 theA = 0;
274 theZ = 0;
275// theS should not be defined
276 break;
277 case KMinus:
278 theA = 0;
279 theZ = -1;
280 theS = -1;
281 break;
282 case Composite:
283 // INCL_ERROR("Trying to set particle type to Composite! Construct a Cluster object instead" << '\n');
284 theA = 0;
285 theZ = 0;
286 theS = 0;
287 break;
288 case UnknownParticle:
289 theA = 0;
290 theZ = 0;
291 theS = 0;
292 INCL_ERROR("Trying to set particle type to Unknown!" << '\n');
293 break;
294 }
295
296 if( !isResonance() && t!=Composite )
297 setINCLMass();
298 }
void setINCLMass()
Set the mass of the Particle to its table mass.

References G4INCL::Composite, G4INCL::DeltaMinus, G4INCL::DeltaPlus, G4INCL::DeltaPlusPlus, G4INCL::DeltaZero, G4INCL::Eta, G4INCL::EtaPrime, INCL_ERROR, isResonance(), G4INCL::KLong, G4INCL::KMinus, G4INCL::KPlus, G4INCL::KShort, G4INCL::KZero, G4INCL::KZeroBar, G4INCL::Lambda, G4INCL::Neutron, G4INCL::Omega, G4INCL::Photon, G4INCL::PiMinus, G4INCL::PiPlus, G4INCL::PiZero, G4INCL::Proton, setINCLMass(), G4INCL::SigmaMinus, G4INCL::SigmaPlus, G4INCL::SigmaZero, theA, theS, theType, theZ, and G4INCL::UnknownParticle.

Referenced by G4INCL::Cluster::Cluster(), G4INCL::ClusterDecay::decay(), G4INCL::DeltaDecayChannel::fillFinalState(), G4INCL::DeltaProductionChannel::fillFinalState(), G4INCL::ElasticChannel::fillFinalState(), G4INCL::EtaNToPiNChannel::fillFinalState(), G4INCL::EtaNToPiPiNChannel::fillFinalState(), G4INCL::NDeltaEtaProductionChannel::fillFinalState(), G4INCL::NDeltaOmegaProductionChannel::fillFinalState(), G4INCL::NDeltaToDeltaLKChannel::fillFinalState(), G4INCL::NDeltaToDeltaSKChannel::fillFinalState(), G4INCL::NDeltaToNLKChannel::fillFinalState(), G4INCL::NDeltaToNNKKbChannel::fillFinalState(), G4INCL::NDeltaToNSKChannel::fillFinalState(), G4INCL::NeutralKaonDecayChannel::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::NNEtaToMultiPionsChannel::fillFinalState(), G4INCL::NNOmegaToMultiPionsChannel::fillFinalState(), G4INCL::NNToMissingStrangenessChannel::fillFinalState(), G4INCL::NNToMultiPionsChannel::fillFinalState(), G4INCL::NNToNLK2piChannel::fillFinalState(), G4INCL::NNToNLKChannel::fillFinalState(), G4INCL::NNToNLKpiChannel::fillFinalState(), G4INCL::NNToNNEtaChannel::fillFinalState(), G4INCL::NNToNNKKbChannel::fillFinalState(), G4INCL::NNToNNOmegaChannel::fillFinalState(), G4INCL::NNToNSK2piChannel::fillFinalState(), G4INCL::NNToNSKChannel::fillFinalState(), G4INCL::NNToNSKpiChannel::fillFinalState(), G4INCL::NpiToMissingStrangenessChannel::fillFinalState(), G4INCL::NSToNLChannel::fillFinalState(), G4INCL::NSToNSChannel::fillFinalState(), G4INCL::OmegaNToPiNChannel::fillFinalState(), G4INCL::OmegaNToPiPiNChannel::fillFinalState(), G4INCL::PiNToMultiPionsChannel::fillFinalState(), G4INCL::PionResonanceDecayChannel::fillFinalState(), G4INCL::RecombinationChannel::fillFinalState(), G4INCL::SigmaZeroDecayChannel::fillFinalState(), G4INCL::StrangeAbsorbtionChannel::fillFinalState(), and Particle().

◆ setUncorrelatedMomentum()

void G4INCL::Particle::setUncorrelatedMomentum ( const G4double  p)
inline

◆ swap()

void G4INCL::Particle::swap ( Particle rhs)
inlineprotected

Helper method for the assignment operator.

Definition at line 125 of file G4INCLParticle.hh.

125 {
126 std::swap(theZ, rhs.theZ);
127 std::swap(theA, rhs.theA);
128 std::swap(theS, rhs.theS);
129 std::swap(theParticipantType, rhs.theParticipantType);
130 std::swap(theType, rhs.theType);
131 if(rhs.thePropagationEnergy == &(rhs.theFrozenEnergy))
133 else
135 std::swap(theEnergy, rhs.theEnergy);
136 std::swap(theFrozenEnergy, rhs.theFrozenEnergy);
137 if(rhs.thePropagationMomentum == &(rhs.theFrozenMomentum))
139 else
141 std::swap(theMomentum, rhs.theMomentum);
142 std::swap(theFrozenMomentum, rhs.theFrozenMomentum);
143 std::swap(thePosition, rhs.thePosition);
144 std::swap(nCollisions, rhs.nCollisions);
145 std::swap(nDecays, rhs.nDecays);
146 std::swap(thePotentialEnergy, rhs.thePotentialEnergy);
147 // ID intentionally not swapped
148
149 std::swap(theHelicity, rhs.theHelicity);
150 std::swap(emissionTime, rhs.emissionTime);
151 std::swap(outOfWell, rhs.outOfWell);
152
153 std::swap(theMass, rhs.theMass);
154 std::swap(rpCorrelated, rhs.rpCorrelated);
155 std::swap(uncorrelatedMomentum, rhs.uncorrelatedMomentum);
156
157 std::swap(theParticleBias, rhs.theParticleBias);
158 std::swap(theBiasCollisionVector, rhs.theBiasCollisionVector);
159
160 }

References emissionTime, nCollisions, nDecays, outOfWell, rpCorrelated, theA, theBiasCollisionVector, theEnergy, theFrozenEnergy, theFrozenMomentum, theHelicity, theMass, theMomentum, theParticipantType, theParticleBias, thePosition, thePotentialEnergy, thePropagationEnergy, thePropagationMomentum, theS, theType, theZ, and uncorrelatedMomentum.

Referenced by operator=(), and G4INCL::Cluster::swap().

◆ thawPropagation()

void G4INCL::Particle::thawPropagation ( )
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 919 of file G4INCLParticle.hh.

References theEnergy, theMomentum, thePropagationEnergy, and thePropagationMomentum.

Referenced by G4INCL::ReflectionChannel::fillFinalState().

Field Documentation

◆ emissionTime

G4double G4INCL::Particle::emissionTime
private

Definition at line 1104 of file G4INCLParticle.hh.

Referenced by getEmissionTime(), setEmissionTime(), and swap().

◆ ID

long G4INCL::Particle::ID
protected

Definition at line 1093 of file G4INCLParticle.hh.

Referenced by dump(), getID(), Particle(), G4INCL::Cluster::print(), and print().

◆ INCLBiasVector

std::vector< G4double > G4INCL::Particle::INCLBiasVector
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 1071 of file G4INCLParticle.hh.

Referenced by FillINCLBiasVector(), getBiasFromVector(), getTotalBias(), G4INCL::INCL::processEvent(), and setINCLBiasVector().

◆ nCollisions

G4int G4INCL::Particle::nCollisions
protected

◆ nDecays

G4int G4INCL::Particle::nDecays
protected

◆ nextBiasedCollisionID

G4ThreadLocal G4int G4INCL::Particle::nextBiasedCollisionID = 0
static

◆ nextID

G4ThreadLocal long G4INCL::Particle::nextID = 1
staticprivate

Definition at line 1111 of file G4INCLParticle.hh.

Referenced by Particle().

◆ outOfWell

G4bool G4INCL::Particle::outOfWell
private

Definition at line 1105 of file G4INCLParticle.hh.

Referenced by isOutOfWell(), setOutOfWell(), and swap().

◆ rpCorrelated

G4bool G4INCL::Particle::rpCorrelated
protected

Definition at line 1095 of file G4INCLParticle.hh.

Referenced by getReflectionMomentum(), rpCorrelate(), rpDecorrelate(), and swap().

◆ theA

G4int G4INCL::Particle::theA
protected

◆ theBiasCollisionVector

std::vector<G4int> G4INCL::Particle::theBiasCollisionVector
private

Time ordered vector of all biased vertices on the particle path.

Definition at line 1108 of file G4INCLParticle.hh.

Referenced by getBiasCollisionVector(), Particle(), setBiasCollisionVector(), and swap().

◆ theEnergy

G4double G4INCL::Particle::theEnergy
protected

◆ theFrozenEnergy

G4double G4INCL::Particle::theFrozenEnergy
protected

◆ theFrozenMomentum

G4INCL::ThreeVector G4INCL::Particle::theFrozenMomentum
protected

◆ theHelicity

G4double G4INCL::Particle::theHelicity
private

Definition at line 1103 of file G4INCLParticle.hh.

Referenced by getHelicity(), setHelicity(), and swap().

◆ theMass

G4double G4INCL::Particle::theMass
private

◆ theMomentum

G4INCL::ThreeVector G4INCL::Particle::theMomentum
protected

◆ theNKaon

G4int G4INCL::Particle::theNKaon
protected

The number of Kaons inside the nucleus (update during the cascade)

Definition at line 1100 of file G4INCLParticle.hh.

Referenced by G4INCL::Nucleus::emitInsideKaon(), getNumberOfKaon(), and setNumberOfKaon().

◆ theParticipantType

ParticipantType G4INCL::Particle::theParticipantType
protected

◆ theParticleBias

G4double G4INCL::Particle::theParticleBias
protected

Definition at line 1098 of file G4INCLParticle.hh.

Referenced by getParticleBias(), setParticleBias(), and swap().

◆ thePosition

G4INCL::ThreeVector G4INCL::Particle::thePosition
protected

◆ thePotentialEnergy

G4double G4INCL::Particle::thePotentialEnergy
protected

◆ thePropagationEnergy

G4double* G4INCL::Particle::thePropagationEnergy
protected

Definition at line 1084 of file G4INCLParticle.hh.

Referenced by freezePropagation(), Particle(), swap(), and thawPropagation().

◆ thePropagationMomentum

G4INCL::ThreeVector* G4INCL::Particle::thePropagationMomentum
protected

◆ theS

G4int G4INCL::Particle::theS
protected

◆ theType

G4INCL::ParticleType G4INCL::Particle::theType
protected

◆ theZ

G4int G4INCL::Particle::theZ
protected

◆ uncorrelatedMomentum

G4double G4INCL::Particle::uncorrelatedMomentum
protected

Definition at line 1096 of file G4INCLParticle.hh.

Referenced by getReflectionMomentum(), setUncorrelatedMomentum(), and swap().


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