Geant4  10.02.p03
G4INCL::Particle Class Reference

#include <G4INCLParticle.hh>

Inheritance diagram for G4INCL::Particle:
Collaboration diagram for G4INCL::Particle:

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. More...
 
Particleoperator= (const Particle &rhs)
 Assignment operator. More...
 
G4INCL::ParticleType getType () const
 
virtual G4INCL::ParticleSpecies getSpecies () const
 Get the particle species. More...
 
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? More...
 
G4bool isResonance () const
 Is it a resonance? More...
 
G4bool isDelta () const
 Is it a Delta? More...
 
G4int getA () const
 Returns the baryon number. More...
 
G4int getZ () const
 Returns the charge number. More...
 
G4double getBeta () const
 
ThreeVector boostVector () const
 
void boost (const ThreeVector &aBoostVector)
 
void lorentzContract (const ThreeVector &aBoostVector, const ThreeVector &refPos)
 Lorentz-contract the particle position around some center. More...
 
G4double getMass () const
 Get the cached particle mass. More...
 
G4double getINCLMass () const
 Get the INCL particle mass. More...
 
virtual G4double getTableMass () const
 Get the tabulated particle mass. More...
 
G4double getRealMass () const
 Get the real particle mass. More...
 
void setRealMass ()
 Set the mass of the Particle to its real mass. More...
 
void setTableMass ()
 Set the mass of the Particle to its table mass. More...
 
void setINCLMass ()
 Set the mass of the Particle to its table mass. More...
 
G4double getEmissionQValueCorrection (const G4int AParent, const G4int ZParent) const
 Computes correction on the emission Q-value. More...
 
G4double getTransferQValueCorrection (const G4int AFrom, const G4int ZFrom, const G4int ATo, const G4int ZTo) const
 Computes correction on the transfer Q-value. More...
 
G4double getInvariantMass () const
 Get the the particle invariant mass. More...
 
G4double getKineticEnergy () const
 Get the particle kinetic energy. More...
 
G4double getPotentialEnergy () const
 Get the particle potential energy. More...
 
void setPotentialEnergy (G4double v)
 Set the particle potential energy. More...
 
G4double getEnergy () const
 
void setMass (G4double mass)
 
void setEnergy (G4double energy)
 
const G4INCL::ThreeVectorgetMomentum () const
 
virtual G4INCL::ThreeVector getAngularMomentum () const
 
virtual void setMomentum (const G4INCL::ThreeVector &momentum)
 
const G4INCL::ThreeVectorgetPosition () 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. More...
 
void setNumberOfCollisions (G4int n)
 Set the number of collisions undergone by the particle. More...
 
void incrementNumberOfCollisions ()
 Increment the number of collisions undergone by the particle. More...
 
G4int getNumberOfDecays () const
 Return the number of decays undergone by the particle. More...
 
void setNumberOfDecays (G4int n)
 Set the number of decays undergone by the particle. More...
 
void incrementNumberOfDecays ()
 Increment the number of decays undergone by the particle. More...
 
void setOutOfWell ()
 Mark the particle as out of its potential well. More...
 
G4bool isOutOfWell () const
 Check if the particle is out of its potential well. More...
 
void setEmissionTime (G4double t)
 
G4double getEmissionTime ()
 
ThreeVector getTransversePosition () const
 Transverse component of the position w.r.t. the momentum. More...
 
ThreeVector getLongitudinalPosition () const
 Longitudinal component of the position w.r.t. the momentum. More...
 
const ThreeVectoradjustMomentumFromEnergy ()
 Rescale the momentum to match the total energy. More...
 
G4double adjustEnergyFromMomentum ()
 Recompute the energy to match the momentum. More...
 
G4bool isCluster () const
 
void setFrozenMomentum (const ThreeVector &momentum)
 Set the frozen particle momentum. More...
 
void setFrozenEnergy (const G4double energy)
 Set the frozen particle momentum. More...
 
ThreeVector getFrozenMomentum () const
 Get the frozen particle momentum. More...
 
G4double getFrozenEnergy () const
 Get the frozen particle momentum. More...
 
ThreeVector getPropagationVelocity () const
 Get the propagation velocity of the particle. More...
 
void freezePropagation ()
 Freeze particle propagation. More...
 
void thawPropagation ()
 Unfreeze particle propagation. More...
 
virtual void rotatePositionAndMomentum (const G4double angle, const ThreeVector &axis)
 Rotate the particle position and momentum. More...
 
virtual void rotatePosition (const G4double angle, const ThreeVector &axis)
 Rotate the particle position. More...
 
virtual void rotateMomentum (const G4double angle, const ThreeVector &axis)
 Rotate the particle momentum. More...
 
std::string print () const
 
std::string dump () const
 
long getID () const
 
ParticleList const * getParticles () const
 
G4double getReflectionMomentum () const
 Return the reflection momentum. More...
 
void setUncorrelatedMomentum (const G4double p)
 Set the uncorrelated momentum. More...
 
void rpCorrelate ()
 Make the particle follow a strict r-p correlation. More...
 
void rpDecorrelate ()
 Make the particle not follow a strict r-p correlation. More...
 
G4double getCosRPAngle () const
 Get the cosine of the angle between position and momentum. More...
 

Protected Member Functions

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

Protected Attributes

G4int theZ
 
G4int theA
 
ParticipantType theParticipantType
 
G4INCL::ParticleType theType
 
G4double theEnergy
 
G4doublethePropagationEnergy
 
G4double theFrozenEnergy
 
G4INCL::ThreeVector theMomentum
 
G4INCL::ThreeVectorthePropagationMomentum
 
G4INCL::ThreeVector theFrozenMomentum
 
G4INCL::ThreeVector thePosition
 
G4int nCollisions
 
G4int nDecays
 
G4double thePotentialEnergy
 
long ID
 
G4bool rpCorrelated
 
G4double uncorrelatedMomentum
 

Private Member Functions

 INCL_DECLARE_ALLOCATION_POOL (Particle)
 

Private Attributes

G4double theHelicity
 
G4double emissionTime
 
G4bool outOfWell
 
G4double theMass
 

Static Private Attributes

static G4ThreadLocal long nextID = 1
 

Detailed Description

Definition at line 73 of file G4INCLParticle.hh.

Constructor & Destructor Documentation

◆ Particle() [1/4]

G4INCL::Particle::Particle ( )

Definition at line 52 of file G4INCLParticle.cc.

53  : theZ(0), theA(0),
56  theEnergy(0.0),
59  theMomentum(ThreeVector(0.,0.,0.)),
62  thePosition(ThreeVector(0.,0.,0.)),
63  nCollisions(0),
64  nDecays(0),
65  thePotentialEnergy(0.0),
66  rpCorrelated(false),
68  theHelicity(0.0),
69  emissionTime(0.0),
70  outOfWell(false),
71  theMass(0.)
72  {
73  ID = nextID;
74  nextID++;
75  }
ParticipantType theParticipantType
G4double * thePropagationEnergy
G4INCL::ThreeVector theFrozenMomentum
G4double uncorrelatedMomentum
static G4ThreadLocal long nextID
G4INCL::ThreeVector thePosition
G4INCL::ParticleType theType
G4INCL::ThreeVector theMomentum
G4double theFrozenEnergy
G4INCL::ThreeVector * thePropagationMomentum
G4double thePotentialEnergy
Here is the caller graph for this function:

◆ Particle() [2/4]

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

Definition at line 77 of file G4INCLParticle.cc.

79  : theEnergy(energy),
82  theMomentum(momentum),
86  nCollisions(0), nDecays(0),
88  rpCorrelated(false),
90  theHelicity(0.0),
91  emissionTime(0.0), outOfWell(false)
92  {
94  ID = nextID;
95  nextID++;
96  if(theEnergy <= 0.0) {
97  INCL_WARN("Particle with energy " << theEnergy << " created." << '\n');
98  }
99  setType(t);
101  }
ParticipantType theParticipantType
void setMass(G4double mass)
G4double * thePropagationEnergy
G4INCL::ThreeVector theFrozenMomentum
#define INCL_WARN(x)
G4double mag() const
G4double uncorrelatedMomentum
double energy
Definition: plottest35.C:25
static G4ThreadLocal long nextID
G4INCL::ThreeVector thePosition
void setType(ParticleType t)
G4INCL::ThreeVector theMomentum
G4double theFrozenEnergy
G4INCL::ThreeVector * thePropagationMomentum
G4double getInvariantMass() const
Get the the particle invariant mass.
G4double thePotentialEnergy
Here is the call graph for this function:

◆ Particle() [3/4]

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

Definition at line 103 of file G4INCLParticle.cc.

106  theMomentum(momentum),
110  nCollisions(0), nDecays(0),
111  thePotentialEnergy(0.),
112  rpCorrelated(false),
114  theHelicity(0.0),
115  emissionTime(0.0), outOfWell(false)
116  {
118  ID = nextID;
119  nextID++;
120  setType(t);
121  if( isResonance() ) {
122  INCL_ERROR("Cannot create resonance without specifying its momentum four-vector." << '\n');
123  }
124  G4double energy = std::sqrt(theMomentum.mag2() + theMass*theMass);
125  theEnergy = energy;
127  }
ParticipantType theParticipantType
G4double * thePropagationEnergy
G4double mag2() const
G4INCL::ThreeVector theFrozenMomentum
#define INCL_ERROR(x)
G4double mag() const
G4double uncorrelatedMomentum
double energy
Definition: plottest35.C:25
static G4ThreadLocal long nextID
G4INCL::ThreeVector thePosition
void setType(ParticleType t)
G4INCL::ThreeVector theMomentum
double G4double
Definition: G4Types.hh:76
G4double theFrozenEnergy
G4INCL::ThreeVector * thePropagationMomentum
G4bool isResonance() const
Is it a resonance?
G4double thePotentialEnergy
Here is the call graph for this function:

◆ ~Particle()

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

Definition at line 78 of file G4INCLParticle.hh.

78 {}

◆ Particle() [4/4]

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

Copy constructor.

Does not copy the particle ID.

Definition at line 84 of file G4INCLParticle.hh.

84  :
85  theZ(rhs.theZ),
86  theA(rhs.theA),
87  theParticipantType(rhs.theParticipantType),
88  theType(rhs.theType),
89  theEnergy(rhs.theEnergy),
90  theFrozenEnergy(rhs.theFrozenEnergy),
91  theMomentum(rhs.theMomentum),
92  theFrozenMomentum(rhs.theFrozenMomentum),
93  thePosition(rhs.thePosition),
94  nCollisions(rhs.nCollisions),
95  nDecays(rhs.nDecays),
96  thePotentialEnergy(rhs.thePotentialEnergy),
97  rpCorrelated(rhs.rpCorrelated),
98  uncorrelatedMomentum(rhs.uncorrelatedMomentum),
99  theHelicity(rhs.theHelicity),
100  emissionTime(rhs.emissionTime),
101  outOfWell(rhs.outOfWell),
102  theMass(rhs.theMass)
103  {
104  if(rhs.thePropagationEnergy == &(rhs.theFrozenEnergy))
106  else
108  if(rhs.thePropagationMomentum == &(rhs.theFrozenMomentum))
110  else
112  // ID intentionally not copied
113  ID = nextID++;
114  }
ParticipantType theParticipantType
G4double * thePropagationEnergy
G4INCL::ThreeVector theFrozenMomentum
G4double uncorrelatedMomentum
static G4ThreadLocal long nextID
G4INCL::ThreeVector thePosition
G4INCL::ParticleType theType
G4INCL::ThreeVector theMomentum
G4double theFrozenEnergy
G4INCL::ThreeVector * thePropagationMomentum
G4double thePotentialEnergy

Member Function Documentation

◆ adjustEnergyFromMomentum()

G4double G4INCL::Particle::adjustEnergyFromMomentum ( )

Recompute the energy to match the momentum.

Definition at line 142 of file G4INCLParticle.cc.

142  {
143  theEnergy = std::sqrt(theMomentum.mag2() + theMass*theMass);
144  return theEnergy;
145  }
G4double mag2() const
G4INCL::ThreeVector theMomentum
Here is the call graph for this function:
Here is the caller graph for this function:

◆ adjustMomentumFromEnergy()

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

Rescale the momentum to match the total energy.

Definition at line 129 of file G4INCLParticle.cc.

129  {
130  const G4double p2 = theMomentum.mag2();
132  if( newp2<0.0 ) {
133  INCL_ERROR("Particle has E^2 < m^2." << '\n' << print());
134  newp2 = 0.0;
135  theEnergy = theMass;
136  }
137 
138  theMomentum *= std::sqrt(newp2/p2);
139  return theMomentum;
140  }
G4double mag2() const
#define INCL_ERROR(x)
std::string print() const
G4INCL::ThreeVector theMomentum
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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 306 of file G4INCLParticle.hh.

306  {
307  const G4double beta2 = aBoostVector.mag2();
308  const G4double gamma = 1.0 / std::sqrt(1.0 - beta2);
309  const G4double bp = theMomentum.dot(aBoostVector);
310  const G4double alpha = (gamma*gamma)/(1.0 + gamma);
311 
312  theMomentum = theMomentum + aBoostVector * (alpha * bp - gamma * theEnergy);
313  theEnergy = gamma * (theEnergy - bp);
314  }
G4double dot(const ThreeVector &v) const
static const G4double bp
G4INCL::ThreeVector theMomentum
double G4double
Definition: G4Types.hh:76
static const G4double alpha
Here is the call graph for this function:
Here is the caller graph for this function:

◆ boostVector()

ThreeVector G4INCL::Particle::boostVector ( ) const
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 296 of file G4INCLParticle.hh.

296  {
297  return theMomentum / theEnergy;
298  }
G4INCL::ThreeVector theMomentum
Here is the caller graph for this function:

◆ dump()

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

Definition at line 727 of file G4INCLParticle.hh.

727  {
728  std::stringstream ss;
729  ss << "(particle " << ID << " ";
731  ss << '\n'
732  << thePosition.dump()
733  << '\n'
734  << theMomentum.dump()
735  << '\n'
736  << theEnergy << ")" << '\n';
737  return ss.str();
738  };
std::string dump() const
G4INCL::ThreeVector thePosition
G4INCL::ParticleType theType
G4INCL::ThreeVector theMomentum
std::string getName(const ParticleType t)
Get the native INCL name of the particle.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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 667 of file G4INCLParticle.hh.

667  {
670  }
G4double * thePropagationEnergy
G4INCL::ThreeVector theFrozenMomentum
G4double theFrozenEnergy
G4INCL::ThreeVector * thePropagationMomentum

◆ getA()

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

Returns the baryon number.

Definition at line 280 of file G4INCLParticle.hh.

280 { return theA; }
Here is the caller graph for this function:

◆ 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 559 of file G4INCLParticle.hh.

560  {
562  };
ThreeVector vector(const ThreeVector &v) const
G4INCL::ThreeVector thePosition
G4INCL::ThreeVector theMomentum
Here is the call graph for this function:
Here is the caller graph for this function:

◆ getBeta()

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

Definition at line 285 of file G4INCLParticle.hh.

285  {
286  const G4double P = theMomentum.mag();
287  return P/theEnergy;
288  }
G4double mag() const
static double P[]
G4INCL::ThreeVector theMomentum
double G4double
Definition: G4Types.hh:76

◆ getCosRPAngle()

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

Get the cosine of the angle between position and momentum.

Definition at line 773 of file G4INCLParticle.hh.

773  {
775  if(norm>0.)
776  return thePosition.dot(*thePropagationMomentum) / std::sqrt(norm);
777  else
778  return 1.;
779  }
G4double mag2() const
Float_t norm
G4double dot(const ThreeVector &v) const
G4INCL::ThreeVector thePosition
double G4double
Definition: G4Types.hh:76
G4INCL::ThreeVector * thePropagationMomentum
Here is the caller graph for this function:

◆ getEmissionQValueCorrection()

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 444 of file G4INCLParticle.hh.

444  {
445  const G4int ADaughter = AParent - theA;
446  const G4int ZDaughter = ZParent - theZ;
447 
448  // Note the minus sign here
449  G4double theQValue;
450  if(isCluster())
451  theQValue = -ParticleTable::getTableQValue(theA, theZ, ADaughter, ZDaughter);
452  else {
453  const G4double massTableParent = ParticleTable::getTableMass(AParent,ZParent);
454  const G4double massTableDaughter = ParticleTable::getTableMass(ADaughter,ZDaughter);
455  const G4double massTableParticle = getTableMass();
456  theQValue = massTableParent - massTableDaughter - massTableParticle;
457  }
458 
459  const G4double massINCLParent = ParticleTable::getINCLMass(AParent,ZParent);
460  const G4double massINCLDaughter = ParticleTable::getINCLMass(ADaughter,ZDaughter);
461  const G4double massINCLParticle = getINCLMass();
462 
463  // The rhs corresponds to the INCL Q-value
464  return theQValue - (massINCLParent-massINCLDaughter-massINCLParticle);
465  }
G4double getINCLMass() const
Get the INCL particle mass.
G4double getTableQValue(const G4int A1, const G4int Z1, const G4int A2, const G4int Z2)
Get Q-value (in MeV/c^2)
int G4int
Definition: G4Types.hh:78
G4bool isCluster() const
virtual G4double getTableMass() const
Get the tabulated particle mass.
G4double getINCLMass(const G4int A, const G4int Z)
Get INCL nuclear mass (in MeV/c^2)
G4ThreadLocal NuclearMassFn getTableMass
Static pointer to the mass function for nuclei.
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ getEmissionTime()

G4double G4INCL::Particle::getEmissionTime ( )
inline

Definition at line 624 of file G4INCLParticle.hh.

624 { return emissionTime; };
Here is the caller graph for this function:

◆ getEnergy()

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

Get the energy of the particle in MeV.

Definition at line 529 of file G4INCLParticle.hh.

530  {
531  return theEnergy;
532  };
Here is the caller graph for this function:

◆ getFrozenEnergy()

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

Get the frozen particle momentum.

Definition at line 656 of file G4INCLParticle.hh.

656 { return theFrozenEnergy; }
G4double theFrozenEnergy

◆ getFrozenMomentum()

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

Get the frozen particle momentum.

Definition at line 653 of file G4INCLParticle.hh.

653 { return theFrozenMomentum; }
G4INCL::ThreeVector theFrozenMomentum

◆ getHelicity()

G4double G4INCL::Particle::getHelicity ( )
inline

Definition at line 585 of file G4INCLParticle.hh.

585 { return theHelicity; };

◆ getID()

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

Definition at line 740 of file G4INCLParticle.hh.

740 { return ID; };
Here is the caller graph for this function:

◆ getINCLMass()

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

Get the INCL particle mass.

Definition at line 338 of file G4INCLParticle.hh.

338  {
339  switch(theType) {
340  case Proton:
341  case Neutron:
342  case PiPlus:
343  case PiMinus:
344  case PiZero:
346  break;
347 
348  case DeltaPlusPlus:
349  case DeltaPlus:
350  case DeltaZero:
351  case DeltaMinus:
352  return theMass;
353  break;
354 
355  case Composite:
357  break;
358 
359  default:
360  INCL_ERROR("Particle::getINCLMass: Unknown particle type." << '\n');
361  return 0.0;
362  break;
363  }
364  }
#define INCL_ERROR(x)
G4INCL::ParticleType theType
G4double getINCLMass(const G4int A, const G4int Z)
Get INCL nuclear mass (in MeV/c^2)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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 507 of file G4INCLParticle.hh.

507  {
508  const G4double mass = std::pow(theEnergy, 2) - theMomentum.dot(theMomentum);
509  if(mass < 0.0) {
510  INCL_ERROR("E*E - p*p is negative." << '\n');
511  return 0.0;
512  } else {
513  return std::sqrt(mass);
514  }
515  };
#define INCL_ERROR(x)
G4double dot(const ThreeVector &v) const
G4INCL::ThreeVector theMomentum
double G4double
Definition: G4Types.hh:76
Here is the caller graph for this function:

◆ getKineticEnergy()

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

Get the particle kinetic energy.

Definition at line 518 of file G4INCLParticle.hh.

518 { return theEnergy - theMass; }
Here is the caller graph for this function:

◆ getLongitudinalPosition()

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

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

Definition at line 632 of file G4INCLParticle.hh.

632  {
634  }
G4double mag2() const
G4double dot(const ThreeVector &v) const
G4INCL::ThreeVector thePosition
G4INCL::ThreeVector * thePropagationMomentum
Here is the caller graph for this function:

◆ getMass()

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

Get the cached particle mass.

Definition at line 335 of file G4INCLParticle.hh.

335 { return theMass; }
Here is the caller graph for this function:

◆ getMomentum()

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

Get the momentum vector.

Definition at line 553 of file G4INCLParticle.hh.

554  {
555  return theMomentum;
556  };
G4INCL::ThreeVector theMomentum
Here is the caller graph for this function:

◆ getNumberOfCollisions()

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

Return the number of collisions undergone by the particle.

Definition at line 593 of file G4INCLParticle.hh.

593 { return nCollisions; }
Here is the caller graph for this function:

◆ getNumberOfDecays()

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

Return the number of decays undergone by the particle.

Definition at line 602 of file G4INCLParticle.hh.

602 { return nDecays; }

◆ getParticipantType()

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

Definition at line 235 of file G4INCLParticle.hh.

235  {
236  return theParticipantType;
237  }
ParticipantType theParticipantType
Here is the caller graph for this function:

◆ getParticles()

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

Return a NULL pointer

Definition at line 745 of file G4INCLParticle.hh.

745  {
746  INCL_WARN("Particle::getParticles() method was called on a Particle object" << '\n');
747  return 0;
748  }
#define INCL_WARN(x)

◆ getPosition()

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

Set the position vector.

Definition at line 575 of file G4INCLParticle.hh.

576  {
577  return thePosition;
578  };
G4INCL::ThreeVector thePosition
Here is the caller graph for this function:

◆ getPotentialEnergy()

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

Get the particle potential energy.

Definition at line 521 of file G4INCLParticle.hh.

521 { return thePotentialEnergy; }
G4double thePotentialEnergy
Here is the caller graph for this function:

◆ getPropagationVelocity()

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

Get the propagation velocity of the particle.

Definition at line 659 of file G4INCLParticle.hh.

659 { return (*thePropagationMomentum)/(*thePropagationEnergy); }
G4INCL::ThreeVector * thePropagationMomentum
Here is the caller graph for this function:

◆ getRealMass()

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

Get the real particle mass.

Definition at line 396 of file G4INCLParticle.hh.

396  {
397  switch(theType) {
398  case Proton:
399  case Neutron:
400  case PiPlus:
401  case PiMinus:
402  case PiZero:
404  break;
405 
406  case DeltaPlusPlus:
407  case DeltaPlus:
408  case DeltaZero:
409  case DeltaMinus:
410  return theMass;
411  break;
412 
413  case Composite:
415  break;
416 
417  default:
418  INCL_ERROR("Particle::getRealMass: Unknown particle type." << '\n');
419  return 0.0;
420  break;
421  }
422  }
#define INCL_ERROR(x)
G4double getRealMass(const G4INCL::ParticleType t)
Get particle mass (in MeV/c^2)
G4INCL::ParticleType theType
Here is the call graph for this function:
Here is the caller graph for this function:

◆ 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 756 of file G4INCLParticle.hh.

756  {
757  if(rpCorrelated)
758  return theMomentum.mag();
759  else
760  return uncorrelatedMomentum;
761  }
G4double mag() const
G4double uncorrelatedMomentum
G4INCL::ThreeVector theMomentum
Here is the caller graph for this function:

◆ getSpecies()

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

Get the particle species.

Reimplemented in G4INCL::Cluster.

Definition at line 171 of file G4INCLParticle.hh.

171  {
172  return ParticleSpecies(theType);
173  };
G4INCL::ParticleType theType
Here is the caller graph for this function:

◆ getTableMass()

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

Get the tabulated particle mass.

Reimplemented in G4INCL::Cluster.

Definition at line 367 of file G4INCLParticle.hh.

367  {
368  switch(theType) {
369  case Proton:
370  case Neutron:
371  case PiPlus:
372  case PiMinus:
373  case PiZero:
375  break;
376 
377  case DeltaPlusPlus:
378  case DeltaPlus:
379  case DeltaZero:
380  case DeltaMinus:
381  return theMass;
382  break;
383 
384  case Composite:
386  break;
387 
388  default:
389  INCL_ERROR("Particle::getTableMass: Unknown particle type." << '\n');
390  return 0.0;
391  break;
392  }
393  }
G4ThreadLocal ParticleMassFn getTableParticleMass
Static pointer to the mass function for particles.
#define INCL_ERROR(x)
G4INCL::ParticleType theType
G4ThreadLocal NuclearMassFn getTableMass
Static pointer to the mass function for nuclei.
Here is the caller graph for this function:

◆ getTransferQValueCorrection()

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 482 of file G4INCLParticle.hh.

482  {
483  const G4int AFromDaughter = AFrom - theA;
484  const G4int ZFromDaughter = ZFrom - theZ;
485  const G4int AToDaughter = ATo + theA;
486  const G4int ZToDaughter = ZTo + theZ;
487  const G4double theQValue = ParticleTable::getTableQValue(AToDaughter,ZToDaughter,AFromDaughter,ZFromDaughter,AFrom,ZFrom);
488 
489  const G4double massINCLTo = ParticleTable::getINCLMass(ATo,ZTo);
490  const G4double massINCLToDaughter = ParticleTable::getINCLMass(AToDaughter,ZToDaughter);
491  /* Note that here we have to use the table mass in the INCL Q-value. We
492  * cannot use theMass, because at this stage the particle is probably
493  * still off-shell; and we cannot use getINCLMass(), because it leads to
494  * violations of global energy conservation.
495  */
496  const G4double massINCLParticle = getTableMass();
497 
498  // The rhs corresponds to the INCL Q-value for particle absorption
499  return theQValue - (massINCLToDaughter-massINCLTo-massINCLParticle);
500  }
G4double getTableQValue(const G4int A1, const G4int Z1, const G4int A2, const G4int Z2)
Get Q-value (in MeV/c^2)
int G4int
Definition: G4Types.hh:78
virtual G4double getTableMass() const
Get the tabulated particle mass.
G4double getINCLMass(const G4int A, const G4int Z)
Get INCL nuclear mass (in MeV/c^2)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ getTransversePosition()

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

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

Definition at line 627 of file G4INCLParticle.hh.

627  {
629  }
G4INCL::ThreeVector thePosition
ThreeVector getLongitudinalPosition() const
Longitudinal component of the position w.r.t. the momentum.
Here is the caller graph for this function:

◆ getType()

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

Get the particle type.

See also
G4INCL::ParticleType

Definition at line 166 of file G4INCLParticle.hh.

166  {
167  return theType;
168  };
G4INCL::ParticleType theType
Here is the caller graph for this function:

◆ getZ()

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

Returns the charge number.

Definition at line 283 of file G4INCLParticle.hh.

283 { return theZ; }
Here is the caller graph for this function:

◆ INCL_DECLARE_ALLOCATION_POOL()

G4INCL::Particle::INCL_DECLARE_ALLOCATION_POOL ( Particle  )
private

◆ incrementNumberOfCollisions()

void G4INCL::Particle::incrementNumberOfCollisions ( )
inline

Increment the number of collisions undergone by the particle.

Definition at line 599 of file G4INCLParticle.hh.

599 { nCollisions++; }

◆ incrementNumberOfDecays()

void G4INCL::Particle::incrementNumberOfDecays ( )
inline

Increment the number of decays undergone by the particle.

Definition at line 608 of file G4INCLParticle.hh.

608 { nDecays++; }

◆ isCluster()

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

Definition at line 642 of file G4INCLParticle.hh.

642  {
643  return (theType == Composite);
644  }
G4INCL::ParticleType theType

◆ isDelta()

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

Is it a Delta?

Definition at line 274 of file G4INCLParticle.hh.

Here is the caller graph for this function:

◆ isNucleon()

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

Is this a nucleon?

Definition at line 228 of file G4INCLParticle.hh.

228  {
230  return true;
231  else
232  return false;
233  };
G4INCL::ParticleType theType
Here is the caller graph for this function:

◆ isOutOfWell()

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

Check if the particle is out of its potential well.

Definition at line 621 of file G4INCLParticle.hh.

621 { return outOfWell; }
Here is the caller graph for this function:

◆ isParticipant()

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

Definition at line 243 of file G4INCLParticle.hh.

243  {
245  }
ParticipantType theParticipantType
Here is the caller graph for this function:

◆ isPion()

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

Is this a pion?

Definition at line 268 of file G4INCLParticle.hh.

Here is the caller graph for this function:

◆ isProjectileSpectator()

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

Definition at line 251 of file G4INCLParticle.hh.

◆ isResonance()

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

Is it a resonance?

Definition at line 271 of file G4INCLParticle.hh.

271 { return isDelta(); }
G4bool isDelta() const
Is it a Delta?
Here is the caller graph for this function:

◆ isTargetSpectator()

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

Definition at line 247 of file G4INCLParticle.hh.

247  {
249  }
ParticipantType theParticipantType
Here is the caller graph for this function:

◆ 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 324 of file G4INCLParticle.hh.

324  {
325  const G4double beta2 = aBoostVector.mag2();
326  const G4double gamma = 1.0 / std::sqrt(1.0 - beta2);
327  const ThreeVector theRelativePosition = thePosition - refPos;
328  const ThreeVector transversePosition = theRelativePosition - aBoostVector * (theRelativePosition.dot(aBoostVector) / aBoostVector.mag2());
329  const ThreeVector longitudinalPosition = theRelativePosition - transversePosition;
330 
331  thePosition = refPos + transversePosition + longitudinalPosition / gamma;
332  }
G4INCL::ThreeVector thePosition
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ makeParticipant()

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

Reimplemented in G4INCL::Cluster.

Definition at line 255 of file G4INCLParticle.hh.

255  {
257  }
ParticipantType theParticipantType
Here is the caller graph for this function:

◆ makeProjectileSpectator()

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

Reimplemented in G4INCL::Cluster.

Definition at line 263 of file G4INCLParticle.hh.

Here is the caller graph for this function:

◆ makeTargetSpectator()

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

Reimplemented in G4INCL::Cluster.

Definition at line 259 of file G4INCLParticle.hh.

Here is the caller graph for this function:

◆ operator=()

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

Assignment operator.

Does not copy the particle ID.

Definition at line 156 of file G4INCLParticle.hh.

156  {
157  Particle temporaryParticle(rhs);
158  swap(temporaryParticle);
159  return *this;
160  }
void swap(Particle &rhs)
Helper method for the assignment operator.
Here is the caller graph for this function:

◆ print()

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

Definition at line 712 of file G4INCLParticle.hh.

712  {
713  std::stringstream ss;
714  ss << "Particle (ID = " << ID << ") type = ";
716  ss << '\n'
717  << " energy = " << theEnergy << '\n'
718  << " momentum = "
719  << theMomentum.print()
720  << '\n'
721  << " position = "
722  << thePosition.print()
723  << '\n';
724  return ss.str();
725  };
std::string print() const
G4INCL::ThreeVector thePosition
G4INCL::ParticleType theType
G4INCL::ThreeVector theMomentum
std::string getName(const ParticleType t)
Get the native INCL name of the particle.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ propagate()

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

Definition at line 588 of file G4INCLParticle.hh.

588  {
589  thePosition += ((*thePropagationMomentum)*(step/(*thePropagationEnergy)));
590  };
G4INCL::ThreeVector thePosition
Here is the caller graph for this function:

◆ 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 707 of file G4INCLParticle.hh.

707  {
708  theMomentum.rotate(angle, axis);
710  }
void rotate(const G4double angle, const ThreeVector &axis)
Rotate the vector by a given angle around a given axis.
G4INCL::ThreeVector theFrozenMomentum
static G4double angle[DIM]
G4INCL::ThreeVector theMomentum
Here is the caller graph for this function:

◆ 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 698 of file G4INCLParticle.hh.

698  {
699  thePosition.rotate(angle, axis);
700  }
void rotate(const G4double angle, const ThreeVector &axis)
Rotate the vector by a given angle around a given axis.
static G4double angle[DIM]
G4INCL::ThreeVector thePosition
Here is the caller graph for this function:

◆ 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 688 of file G4INCLParticle.hh.

688  {
689  rotatePosition(angle, axis);
690  rotateMomentum(angle, axis);
691  }
virtual void rotatePosition(const G4double angle, const ThreeVector &axis)
Rotate the particle position.
static G4double angle[DIM]
virtual void rotateMomentum(const G4double angle, const ThreeVector &axis)
Rotate the particle momentum.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ rpCorrelate()

void G4INCL::Particle::rpCorrelate ( )
inline

Make the particle follow a strict r-p correlation.

Definition at line 767 of file G4INCLParticle.hh.

767 { rpCorrelated = true; }
Here is the caller graph for this function:

◆ rpDecorrelate()

void G4INCL::Particle::rpDecorrelate ( )
inline

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

Definition at line 770 of file G4INCLParticle.hh.

770 { rpCorrelated = false; }

◆ setEmissionTime()

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

Definition at line 623 of file G4INCLParticle.hh.

623 { emissionTime = t; }
Here is the caller graph for this function:

◆ setEnergy()

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

Set the energy of the particle in MeV.

Definition at line 545 of file G4INCLParticle.hh.

546  {
547  this->theEnergy = energy;
548  };
double energy
Definition: plottest35.C:25
Here is the caller graph for this function:

◆ setFrozenEnergy()

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

Set the frozen particle momentum.

Definition at line 650 of file G4INCLParticle.hh.

650 { theFrozenEnergy = energy; }
double energy
Definition: plottest35.C:25
G4double theFrozenEnergy

◆ setFrozenMomentum()

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

Set the frozen particle momentum.

Definition at line 647 of file G4INCLParticle.hh.

647 { theFrozenMomentum = momentum; }
G4INCL::ThreeVector theFrozenMomentum

◆ setHelicity()

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

Definition at line 586 of file G4INCLParticle.hh.

586 { theHelicity = h; };
Here is the caller graph for this function:

◆ setINCLMass()

void G4INCL::Particle::setINCLMass ( )
inline

Set the mass of the Particle to its table mass.

Definition at line 431 of file G4INCLParticle.hh.

431 { setMass(getINCLMass()); }
void setMass(G4double mass)
G4double getINCLMass() const
Get the INCL particle mass.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ setMass()

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

Set the mass of the particle in MeV/c^2.

Definition at line 537 of file G4INCLParticle.hh.

538  {
539  this->theMass = mass;
540  }
Here is the caller graph for this function:

◆ setMomentum()

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

Set the momentum vector.

Definition at line 567 of file G4INCLParticle.hh.

568  {
569  this->theMomentum = momentum;
570  };
G4INCL::ThreeVector theMomentum
Here is the caller graph for this function:

◆ setNumberOfCollisions()

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

Set the number of collisions undergone by the particle.

Definition at line 596 of file G4INCLParticle.hh.

596 { nCollisions = n; }
Char_t n[5]

◆ setNumberOfDecays()

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

Set the number of decays undergone by the particle.

Definition at line 605 of file G4INCLParticle.hh.

605 { nDecays = n; }
Char_t n[5]

◆ 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 618 of file G4INCLParticle.hh.

618 { outOfWell = true; }

◆ setParticipantType()

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

Definition at line 239 of file G4INCLParticle.hh.

239  {
240  theParticipantType = p;
241  }
ParticipantType theParticipantType

◆ setPosition()

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

Reimplemented in G4INCL::Cluster.

Definition at line 580 of file G4INCLParticle.hh.

581  {
582  this->thePosition = position;
583  };
#define position
Definition: xmlparse.cc:622
G4INCL::ThreeVector thePosition
Here is the caller graph for this function:

◆ setPotentialEnergy()

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

Set the particle potential energy.

Definition at line 524 of file G4INCLParticle.hh.

Here is the caller graph for this function:

◆ setRealMass()

void G4INCL::Particle::setRealMass ( )
inline

Set the mass of the Particle to its real mass.

Definition at line 425 of file G4INCLParticle.hh.

425 { setMass(getRealMass()); }
void setMass(G4double mass)
G4double getRealMass() const
Get the real particle mass.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ setTableMass()

void G4INCL::Particle::setTableMass ( )
inline

Set the mass of the Particle to its table mass.

Definition at line 428 of file G4INCLParticle.hh.

428 { setMass(getTableMass()); }
void setMass(G4double mass)
virtual G4double getTableMass() const
Get the tabulated particle mass.
Here is the caller graph for this function:

◆ setType()

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

Definition at line 175 of file G4INCLParticle.hh.

175  {
176  theType = t;
177  switch(theType)
178  {
179  case DeltaPlusPlus:
180  theA = 1;
181  theZ = 2;
182  break;
183  case Proton:
184  case DeltaPlus:
185  theA = 1;
186  theZ = 1;
187  break;
188  case Neutron:
189  case DeltaZero:
190  theA = 1;
191  theZ = 0;
192  break;
193  case DeltaMinus:
194  theA = 1;
195  theZ = -1;
196  break;
197  case PiPlus:
198  theA = 0;
199  theZ = 1;
200  break;
201  case PiZero:
202  theA = 0;
203  theZ = 0;
204  break;
205  case PiMinus:
206  theA = 0;
207  theZ = -1;
208  break;
209  case Composite:
210  // INCL_ERROR("Trying to set particle type to Composite! Construct a Cluster object instead" << '\n');
211  theA = 0;
212  theZ = 0;
213  break;
214  case UnknownParticle:
215  theA = 0;
216  theZ = 0;
217  INCL_ERROR("Trying to set particle type to Unknown!" << '\n');
218  break;
219  }
220 
221  if( !isResonance() && t!=Composite )
222  setINCLMass();
223  }
#define INCL_ERROR(x)
void setINCLMass()
Set the mass of the Particle to its table mass.
G4INCL::ParticleType theType
G4bool isResonance() const
Is it a resonance?
Here is the caller graph for this function:

◆ setUncorrelatedMomentum()

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

Set the uncorrelated momentum.

Definition at line 764 of file G4INCLParticle.hh.

764 { uncorrelatedMomentum = p; }
G4double uncorrelatedMomentum
Here is the caller graph for this function:

◆ swap()

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

Helper method for the assignment operator.

Definition at line 118 of file G4INCLParticle.hh.

118  {
119  std::swap(theZ, rhs.theZ);
120  std::swap(theA, rhs.theA);
121  std::swap(theParticipantType, rhs.theParticipantType);
122  std::swap(theType, rhs.theType);
123  if(rhs.thePropagationEnergy == &(rhs.theFrozenEnergy))
125  else
127  std::swap(theEnergy, rhs.theEnergy);
128  std::swap(theFrozenEnergy, rhs.theFrozenEnergy);
129  if(rhs.thePropagationMomentum == &(rhs.theFrozenMomentum))
131  else
133  std::swap(theMomentum, rhs.theMomentum);
134  std::swap(theFrozenMomentum, rhs.theFrozenMomentum);
135  std::swap(thePosition, rhs.thePosition);
136  std::swap(nCollisions, rhs.nCollisions);
137  std::swap(nDecays, rhs.nDecays);
138  std::swap(thePotentialEnergy, rhs.thePotentialEnergy);
139  // ID intentionally not swapped
140 
141  std::swap(theHelicity, rhs.theHelicity);
142  std::swap(emissionTime, rhs.emissionTime);
143  std::swap(outOfWell, rhs.outOfWell);
144 
145  std::swap(theMass, rhs.theMass);
146  std::swap(rpCorrelated, rhs.rpCorrelated);
147  std::swap(uncorrelatedMomentum, rhs.uncorrelatedMomentum);
148  }
ParticipantType theParticipantType
G4double * thePropagationEnergy
G4INCL::ThreeVector theFrozenMomentum
G4double uncorrelatedMomentum
G4INCL::ThreeVector thePosition
G4INCL::ParticleType theType
G4INCL::ThreeVector theMomentum
G4double theFrozenEnergy
G4INCL::ThreeVector * thePropagationMomentum
G4double thePotentialEnergy
Here is the caller graph for this function:

◆ 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 678 of file G4INCLParticle.hh.

678  {
681  }
G4double * thePropagationEnergy
G4INCL::ThreeVector theMomentum
G4INCL::ThreeVector * thePropagationMomentum

Member Data Documentation

◆ emissionTime

G4double G4INCL::Particle::emissionTime
private

Definition at line 802 of file G4INCLParticle.hh.

◆ ID

long G4INCL::Particle::ID
protected

Definition at line 795 of file G4INCLParticle.hh.

◆ nCollisions

G4int G4INCL::Particle::nCollisions
protected

Definition at line 792 of file G4INCLParticle.hh.

◆ nDecays

G4int G4INCL::Particle::nDecays
protected

Definition at line 793 of file G4INCLParticle.hh.

◆ nextID

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

Definition at line 806 of file G4INCLParticle.hh.

◆ outOfWell

G4bool G4INCL::Particle::outOfWell
private

Definition at line 803 of file G4INCLParticle.hh.

◆ rpCorrelated

G4bool G4INCL::Particle::rpCorrelated
protected

Definition at line 797 of file G4INCLParticle.hh.

◆ theA

G4int G4INCL::Particle::theA
protected

Definition at line 782 of file G4INCLParticle.hh.

◆ theEnergy

G4double G4INCL::Particle::theEnergy
protected

Definition at line 785 of file G4INCLParticle.hh.

◆ theFrozenEnergy

G4double G4INCL::Particle::theFrozenEnergy
protected

Definition at line 787 of file G4INCLParticle.hh.

◆ theFrozenMomentum

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

Definition at line 790 of file G4INCLParticle.hh.

◆ theHelicity

G4double G4INCL::Particle::theHelicity
private

Definition at line 801 of file G4INCLParticle.hh.

◆ theMass

G4double G4INCL::Particle::theMass
private

Definition at line 805 of file G4INCLParticle.hh.

◆ theMomentum

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

Definition at line 788 of file G4INCLParticle.hh.

◆ theParticipantType

ParticipantType G4INCL::Particle::theParticipantType
protected

Definition at line 783 of file G4INCLParticle.hh.

◆ thePosition

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

Definition at line 791 of file G4INCLParticle.hh.

◆ thePotentialEnergy

G4double G4INCL::Particle::thePotentialEnergy
protected

Definition at line 794 of file G4INCLParticle.hh.

◆ thePropagationEnergy

G4double* G4INCL::Particle::thePropagationEnergy
protected

Definition at line 786 of file G4INCLParticle.hh.

◆ thePropagationMomentum

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

Definition at line 789 of file G4INCLParticle.hh.

◆ theType

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

Definition at line 784 of file G4INCLParticle.hh.

◆ theZ

G4int G4INCL::Particle::theZ
protected

Definition at line 782 of file G4INCLParticle.hh.

◆ uncorrelatedMomentum

G4double G4INCL::Particle::uncorrelatedMomentum
protected

Definition at line 798 of file G4INCLParticle.hh.


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