Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions | Protected Attributes
G4INCL::Particle Class Reference

#include <G4INCLParticle.hh>

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

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 isInList (ParticleList const &l) const
 Check if the particle belongs to a given list. 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 rotate (const G4double angle, const ThreeVector &axis)
 Rotate the particle position and momentum. More...
 
std::string print () const
 
std::string dump () const
 
long getID () const
 
ParticleList const * getParticles () const
 
G4double getReflectionMomentum () const
 Return the reflection momentum. More...
 
void setUncorrelatedMomentum (const G4double p)
 Set the uncorrelated momentum. More...
 
void rpCorrelate ()
 Make the particle follow a strict r-p correlation. More...
 
void rpDecorrelate ()
 Make the particle not follow a strict r-p correlation. More...
 
G4double getCosRPAngle () const
 Get the cosine of the angle between position and momentum. More...
 

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
 

Detailed Description

Definition at line 94 of file G4INCLParticle.hh.

Constructor & Destructor Documentation

G4INCL::Particle::Particle ( )

Definition at line 51 of file G4INCLParticle.cc.

References ID.

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

52  : theZ(0), theA(0),
55  theEnergy(0.0),
58  theMomentum(ThreeVector(0.,0.,0.)),
61  thePosition(ThreeVector(0.,0.,0.)),
62  nCollisions(0),
63  nDecays(0),
64  thePotentialEnergy(0.0),
65  rpCorrelated(false),
67  theHelicity(0.0),
68  emissionTime(0.0),
69  outOfWell(false),
70  theMass(0.)
71  {
72  ID = nextID;
73  nextID++;
74  }
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
G4INCL::Particle::Particle ( ParticleType  t,
G4double  energy,
ThreeVector const &  momentum,
ThreeVector const &  position 
)

Definition at line 76 of file G4INCLParticle.cc.

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

78  : theEnergy(energy),
81  theMomentum(momentum),
85  nCollisions(0), nDecays(0),
87  rpCorrelated(false),
89  theHelicity(0.0),
90  emissionTime(0.0), outOfWell(false)
91  {
93  ID = nextID;
94  nextID++;
95  if(theEnergy <= 0.0) {
96  INCL_WARN("Particle with energy " << theEnergy << " created." << std::endl);
97  }
98  setType(t);
100  }
ParticipantType theParticipantType
void setMass(G4double mass)
G4double * thePropagationEnergy
G4INCL::ThreeVector theFrozenMomentum
#define INCL_WARN(x)
G4double uncorrelatedMomentum
G4double getInvariantMass() const
Get the the particle invariant mass.
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
G4INCL::ThreeVector thePosition
void setType(ParticleType t)
G4INCL::ThreeVector theMomentum
G4double mag() const
G4double theFrozenEnergy
G4INCL::ThreeVector * thePropagationMomentum
G4double thePotentialEnergy
G4INCL::Particle::Particle ( ParticleType  t,
ThreeVector const &  momentum,
ThreeVector const &  position 
)

Definition at line 102 of file G4INCLParticle.cc.

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

105  theMomentum(momentum),
109  nCollisions(0), nDecays(0),
110  thePotentialEnergy(0.),
111  rpCorrelated(false),
113  theHelicity(0.0),
114  emissionTime(0.0), outOfWell(false)
115  {
117  ID = nextID;
118  nextID++;
119  setType(t);
120  if( isResonance() ) {
121  INCL_ERROR("Cannot create resonance without specifying its momentum four-vector." << std::endl);
122  }
123  G4double energy = std::sqrt(theMomentum.mag2() + theMass*theMass);
124  theEnergy = energy;
126  }
ParticipantType theParticipantType
G4bool isResonance() const
Is it a resonance?
G4double * thePropagationEnergy
G4INCL::ThreeVector theFrozenMomentum
#define INCL_ERROR(x)
G4double uncorrelatedMomentum
G4double mag2() const
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
G4INCL::ThreeVector thePosition
void setType(ParticleType t)
G4INCL::ThreeVector theMomentum
G4double mag() const
double G4double
Definition: G4Types.hh:76
G4double theFrozenEnergy
G4INCL::ThreeVector * thePropagationMomentum
G4double thePotentialEnergy
virtual G4INCL::Particle::~Particle ( )
inlinevirtual

Definition at line 99 of file G4INCLParticle.hh.

99 {}
G4INCL::Particle::Particle ( const Particle rhs)
inline

Copy constructor.

Does not copy the particle ID.

Definition at line 105 of file G4INCLParticle.hh.

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

105  :
106  theZ(rhs.theZ),
107  theA(rhs.theA),
108  theParticipantType(rhs.theParticipantType),
109  theType(rhs.theType),
110  theEnergy(rhs.theEnergy),
111  theFrozenEnergy(rhs.theFrozenEnergy),
112  theMomentum(rhs.theMomentum),
113  theFrozenMomentum(rhs.theFrozenMomentum),
114  thePosition(rhs.thePosition),
115  nCollisions(rhs.nCollisions),
116  nDecays(rhs.nDecays),
117  thePotentialEnergy(rhs.thePotentialEnergy),
118  rpCorrelated(rhs.rpCorrelated),
119  uncorrelatedMomentum(rhs.uncorrelatedMomentum),
120  theHelicity(rhs.theHelicity),
121  emissionTime(rhs.emissionTime),
122  outOfWell(rhs.outOfWell),
123  theMass(rhs.theMass)
124  {
125  if(rhs.thePropagationEnergy == &(rhs.theFrozenEnergy))
127  else
129  if(rhs.thePropagationMomentum == &(rhs.theFrozenMomentum))
131  else
133  // ID intentionally not copied
134  ID = nextID++;
135  }
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

Member Function Documentation

G4double G4INCL::Particle::adjustEnergyFromMomentum ( )

Recompute the energy to match the momentum.

Definition at line 141 of file G4INCLParticle.cc.

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

Referenced by G4INCL::Nucleus::computeRecoilKinematics(), G4INCL::PhaseSpaceDecay::decay(), G4INCL::Nucleus::decayOutgoingDeltas(), G4INCL::DeltaDecayChannel::getFinalState(), and G4INCL::RecombinationChannel::getFinalState().

141  {
142  theEnergy = std::sqrt(theMomentum.mag2() + theMass*theMass);
143  return theEnergy;
144  }
G4double mag2() const
G4INCL::ThreeVector theMomentum
const ThreeVector & G4INCL::Particle::adjustMomentumFromEnergy ( )

Rescale the momentum to match the total energy.

Definition at line 128 of file G4INCLParticle.cc.

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

Referenced by G4INCL::Cluster::Cluster(), G4INCL::Nucleus::emitInsidePions(), G4INCL::CrossSections::interactionDistanceNN(), G4INCL::CrossSections::interactionDistancePiN(), G4INCL::StandardPropagationModel::shootParticle(), and G4INCL::KinematicsUtils::transformToLocalEnergyFrame().

128  {
129  const G4double p2 = theMomentum.mag2();
130  G4double newp2 = theEnergy*theEnergy - theMass*theMass;
131  if( newp2<0.0 ) {
132  INCL_ERROR("Particle has E^2 < m^2." << std::endl << print());
133  newp2 = 0.0;
134  theEnergy = theMass;
135  }
136 
137  theMomentum *= std::sqrt(newp2/p2);
138  return theMomentum;
139  }
#define INCL_ERROR(x)
std::string print() const
G4double mag2() const
G4INCL::ThreeVector theMomentum
double G4double
Definition: G4Types.hh:76
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 327 of file G4INCLParticle.hh.

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

Referenced by G4INCL::Cluster::boost(), G4INCL::PhaseSpaceDecay::decay(), G4INCL::Nucleus::decayOutgoingDeltas(), and G4INCL::InteractionAvatar::preInteraction().

327  {
328  const G4double beta2 = aBoostVector.mag2();
329  const G4double gamma = 1.0 / std::sqrt(1.0 - beta2);
330  const G4double bp = theMomentum.dot(aBoostVector);
331  const G4double alpha = (gamma*gamma)/(1.0 + gamma);
332 
333  theMomentum = theMomentum + aBoostVector * (alpha * bp - gamma * theEnergy);
334  theEnergy = gamma * (theEnergy - bp);
335  }
G4double dot(const ThreeVector &v) const
G4INCL::ThreeVector theMomentum
double G4double
Definition: G4Types.hh:76
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 317 of file G4INCLParticle.hh.

References theEnergy, and theMomentum.

Referenced by G4INCL::PhaseSpaceDecay::decay(), G4INCL::StandardPropagationModel::shootComposite(), and G4INCL::StandardPropagationModel::shootParticle().

317  {
318  return theMomentum / theEnergy;
319  }
G4INCL::ThreeVector theMomentum
std::string G4INCL::Particle::dump ( ) const
inline

Definition at line 735 of file G4INCLParticle.hh.

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

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

735  {
736  std::stringstream ss;
737  ss << "(particle " << ID << " ";
739  ss << std::endl
740  << thePosition.dump()
741  << std::endl
742  << theMomentum.dump()
743  << std::endl
744  << theEnergy << ")" << std::endl;
745  return ss.str();
746  };
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.
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 693 of file G4INCLParticle.hh.

References theFrozenEnergy, theFrozenMomentum, thePropagationEnergy, and thePropagationMomentum.

693  {
696  }
G4double * thePropagationEnergy
G4INCL::ThreeVector theFrozenMomentum
G4double theFrozenEnergy
G4INCL::ThreeVector * thePropagationMomentum
G4int G4INCL::Particle::getA ( ) const
inline
virtual G4INCL::ThreeVector G4INCL::Particle::getAngularMomentum ( ) const
inlinevirtual

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

Reimplemented in G4INCL::Cluster.

Definition at line 580 of file G4INCLParticle.hh.

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

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

581  {
583  };
ThreeVector vector(const ThreeVector &v) const
G4INCL::ThreeVector thePosition
G4INCL::ThreeVector theMomentum
G4double G4INCL::Particle::getBeta ( ) const
inline

Definition at line 306 of file G4INCLParticle.hh.

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

306  {
307  const G4double P = theMomentum.mag();
308  return P/theEnergy;
309  }
G4INCL::ThreeVector theMomentum
G4double mag() const
double G4double
Definition: G4Types.hh:76
G4double G4INCL::Particle::getCosRPAngle ( ) const
inline

Get the cosine of the angle between position and momentum.

Definition at line 781 of file G4INCLParticle.hh.

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

781  {
783  if(norm>0.)
784  return thePosition.dot(*thePropagationMomentum) / std::sqrt(norm);
785  else
786  return 1.;
787  }
G4double dot(const ThreeVector &v) const
G4double mag2() const
G4INCL::ThreeVector thePosition
double G4double
Definition: G4Types.hh:76
G4INCL::ThreeVector * thePropagationMomentum
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 465 of file G4INCLParticle.hh.

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

Referenced by G4INCL::Nucleus::emitInsidePions(), G4INCL::ParticleEntryChannel::getFinalState(), G4INCL::TransmissionChannel::getFinalState(), and G4INCL::SurfaceAvatar::getTransmissionProbability().

465  {
466  const G4int ADaughter = AParent - theA;
467  const G4int ZDaughter = ZParent - theZ;
468 
469  // Note the minus sign here
470  G4double theQValue;
471  if(isCluster())
472  theQValue = -ParticleTable::getTableQValue(theA, theZ, ADaughter, ZDaughter);
473  else {
474  const G4double massTableParent = ParticleTable::getTableMass(AParent,ZParent);
475  const G4double massTableDaughter = ParticleTable::getTableMass(ADaughter,ZDaughter);
476  const G4double massTableParticle = getTableMass();
477  theQValue = massTableParent - massTableDaughter - massTableParticle;
478  }
479 
480  const G4double massINCLParent = ParticleTable::getINCLMass(AParent,ZParent);
481  const G4double massINCLDaughter = ParticleTable::getINCLMass(ADaughter,ZDaughter);
482  const G4double massINCLParticle = getINCLMass();
483 
484  // The rhs corresponds to the INCL Q-value
485  return theQValue - (massINCLParent-massINCLDaughter-massINCLParticle);
486  }
G4double getTableQValue(const G4int A1, const G4int Z1, const G4int A2, const G4int Z2)
Get Q-value (in MeV/c^2)
G4double getINCLMass() const
Get the INCL particle mass.
int G4int
Definition: G4Types.hh:78
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.
G4bool isCluster() const
virtual G4double getTableMass() const
Get the tabulated particle mass.
double G4double
Definition: G4Types.hh:76
G4double G4INCL::Particle::getEmissionTime ( )
inline

Definition at line 645 of file G4INCLParticle.hh.

645 { return emissionTime; };
G4double G4INCL::Particle::getEnergy ( ) const
inline
G4double G4INCL::Particle::getFrozenEnergy ( ) const
inline

Get the frozen particle momentum.

Definition at line 682 of file G4INCLParticle.hh.

References theFrozenEnergy.

682 { return theFrozenEnergy; }
G4double theFrozenEnergy
ThreeVector G4INCL::Particle::getFrozenMomentum ( ) const
inline

Get the frozen particle momentum.

Definition at line 679 of file G4INCLParticle.hh.

References theFrozenMomentum.

679 { return theFrozenMomentum; }
G4INCL::ThreeVector theFrozenMomentum
G4double G4INCL::Particle::getHelicity ( )
inline

Definition at line 606 of file G4INCLParticle.hh.

606 { return theHelicity; };
long G4INCL::Particle::getID ( ) const
inline
G4double G4INCL::Particle::getINCLMass ( ) const
inline

Get the INCL particle mass.

Definition at line 359 of file G4INCLParticle.hh.

References G4INCL::Composite, G4INCL::DeltaMinus, G4INCL::DeltaPlus, G4INCL::DeltaPlusPlus, G4INCL::DeltaZero, G4INCL::ParticleTable::getINCLMass(), INCL_ERROR, G4INCL::Neutron, G4INCL::PiMinus, G4INCL::PiPlus, G4INCL::PiZero, G4INCL::Proton, theA, theType, and theZ.

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

359  {
360  switch(theType) {
361  case Proton:
362  case Neutron:
363  case PiPlus:
364  case PiMinus:
365  case PiZero:
367  break;
368 
369  case DeltaPlusPlus:
370  case DeltaPlus:
371  case DeltaZero:
372  case DeltaMinus:
373  return theMass;
374  break;
375 
376  case Composite:
378  break;
379 
380  default:
381  INCL_ERROR("Particle::getINCLMass: Unknown particle type." << std::endl);
382  return 0.0;
383  break;
384  }
385  }
#define INCL_ERROR(x)
G4INCL::ParticleType theType
G4double getINCLMass(const G4int A, const G4int Z)
Get INCL nuclear mass (in MeV/c^2)
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 528 of file G4INCLParticle.hh.

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

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

528  {
529  const G4double mass = std::pow(theEnergy, 2) - theMomentum.dot(theMomentum);
530  if(mass < 0.0) {
531  INCL_ERROR("E*E - p*p is negative." << std::endl);
532  return 0.0;
533  } else {
534  return std::sqrt(mass);
535  }
536  };
G4double dot(const ThreeVector &v) const
#define INCL_ERROR(x)
G4INCL::ThreeVector theMomentum
double G4double
Definition: G4Types.hh:76
G4double G4INCL::Particle::getKineticEnergy ( ) const
inline
ThreeVector G4INCL::Particle::getLongitudinalPosition ( ) const
inline

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

Definition at line 653 of file G4INCLParticle.hh.

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

Referenced by getTransversePosition().

653  {
655  }
G4double dot(const ThreeVector &v) const
G4double mag2() const
G4INCL::ThreeVector thePosition
G4INCL::ThreeVector * thePropagationMomentum
G4double G4INCL::Particle::getMass ( ) const
inline
const G4INCL::ThreeVector& G4INCL::Particle::getMomentum ( ) const
inline
G4int G4INCL::Particle::getNumberOfCollisions ( ) const
inline

Return the number of collisions undergone by the particle.

Definition at line 614 of file G4INCLParticle.hh.

References nCollisions.

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

614 { return nCollisions; }
G4int G4INCL::Particle::getNumberOfDecays ( ) const
inline

Return the number of decays undergone by the particle.

Definition at line 623 of file G4INCLParticle.hh.

References nDecays.

623 { return nDecays; }
ParticipantType G4INCL::Particle::getParticipantType ( ) const
inline

Definition at line 256 of file G4INCLParticle.hh.

References theParticipantType.

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

256  {
257  return theParticipantType;
258  }
ParticipantType theParticipantType
ParticleList const* G4INCL::Particle::getParticles ( ) const
inline

Return a NULL pointer

Definition at line 753 of file G4INCLParticle.hh.

References INCL_WARN.

753  {
754  INCL_WARN("Particle::getParticles() method was called on a Particle object" << std::endl);
755  return 0;
756  }
#define INCL_WARN(x)
const G4INCL::ThreeVector& G4INCL::Particle::getPosition ( ) const
inline
G4double G4INCL::Particle::getPotentialEnergy ( ) const
inline
ThreeVector G4INCL::Particle::getPropagationVelocity ( ) const
inline

Get the propagation velocity of the particle.

Definition at line 685 of file G4INCLParticle.hh.

References thePropagationMomentum.

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

685 { return (*thePropagationMomentum)/(*thePropagationEnergy); }
G4INCL::ThreeVector * thePropagationMomentum
G4double G4INCL::Particle::getRealMass ( ) const
inline

Get the real particle mass.

Definition at line 417 of file G4INCLParticle.hh.

References G4INCL::Composite, G4INCL::DeltaMinus, G4INCL::DeltaPlus, G4INCL::DeltaPlusPlus, G4INCL::DeltaZero, G4INCL::ParticleTable::getRealMass(), INCL_ERROR, G4INCL::Neutron, G4INCL::PiMinus, G4INCL::PiPlus, G4INCL::PiZero, G4INCL::Proton, theA, theType, and theZ.

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

417  {
418  switch(theType) {
419  case Proton:
420  case Neutron:
421  case PiPlus:
422  case PiMinus:
423  case PiZero:
425  break;
426 
427  case DeltaPlusPlus:
428  case DeltaPlus:
429  case DeltaZero:
430  case DeltaMinus:
431  return theMass;
432  break;
433 
434  case Composite:
436  break;
437 
438  default:
439  INCL_ERROR("Particle::getRealMass: Unknown particle type." << std::endl);
440  return 0.0;
441  break;
442  }
443  }
#define INCL_ERROR(x)
G4double getRealMass(const G4INCL::ParticleType t)
Get particle mass (in MeV/c^2)
G4INCL::ParticleType theType
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 764 of file G4INCLParticle.hh.

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

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

764  {
765  if(rpCorrelated)
766  return theMomentum.mag();
767  else
768  return uncorrelatedMomentum;
769  }
G4double uncorrelatedMomentum
G4INCL::ThreeVector theMomentum
G4double mag() const
virtual G4INCL::ParticleSpecies G4INCL::Particle::getSpecies ( ) const
inlinevirtual

Get the particle species.

Reimplemented in G4INCL::Cluster.

Definition at line 192 of file G4INCLParticle.hh.

References theType.

Referenced by G4INCL::CoulombDistortion::maxImpactParameter(), and G4INCL::StandardPropagationModel::shootParticle().

192  {
193  return ParticleSpecies(theType);
194  };
G4INCL::ParticleType theType
virtual G4double G4INCL::Particle::getTableMass ( ) const
inlinevirtual

Get the tabulated particle mass.

Reimplemented in G4INCL::Cluster.

Definition at line 388 of file G4INCLParticle.hh.

References G4INCL::Composite, G4INCL::DeltaMinus, G4INCL::DeltaPlus, G4INCL::DeltaPlusPlus, G4INCL::DeltaZero, G4INCL::ParticleTable::getTableMass, G4INCL::ParticleTable::getTableParticleMass, INCL_ERROR, G4INCL::Neutron, G4INCL::PiMinus, G4INCL::PiPlus, G4INCL::PiZero, G4INCL::Proton, theA, theType, and theZ.

Referenced by G4INCL::Nucleus::decayOutgoingDeltas(), getEmissionQValueCorrection(), G4INCL::ParticleEntryChannel::getFinalState(), G4INCL::TransmissionChannel::getFinalState(), getTransferQValueCorrection(), and setTableMass().

388  {
389  switch(theType) {
390  case Proton:
391  case Neutron:
392  case PiPlus:
393  case PiMinus:
394  case PiZero:
396  break;
397 
398  case DeltaPlusPlus:
399  case DeltaPlus:
400  case DeltaZero:
401  case DeltaMinus:
402  return theMass;
403  break;
404 
405  case Composite:
407  break;
408 
409  default:
410  INCL_ERROR("Particle::getTableMass: Unknown particle type." << std::endl);
411  return 0.0;
412  break;
413  }
414  }
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.
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 503 of file G4INCLParticle.hh.

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

503  {
504  const G4int AFromDaughter = AFrom - theA;
505  const G4int ZFromDaughter = ZFrom - theZ;
506  const G4int AToDaughter = ATo + theA;
507  const G4int ZToDaughter = ZTo + theZ;
508  const G4double theQValue = ParticleTable::getTableQValue(AToDaughter,ZToDaughter,AFromDaughter,ZFromDaughter,AFrom,ZFrom);
509 
510  const G4double massINCLTo = ParticleTable::getINCLMass(ATo,ZTo);
511  const G4double massINCLToDaughter = ParticleTable::getINCLMass(AToDaughter,ZToDaughter);
512  /* Note that here we have to use the table mass in the INCL Q-value. We
513  * cannot use theMass, because at this stage the particle is probably
514  * still off-shell; and we cannot use getINCLMass(), because it leads to
515  * violations of global energy conservation.
516  */
517  const G4double massINCLParticle = getTableMass();
518 
519  // The rhs corresponds to the INCL Q-value for particle absorption
520  return theQValue - (massINCLToDaughter-massINCLTo-massINCLParticle);
521  }
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
G4double getINCLMass(const G4int A, const G4int Z)
Get INCL nuclear mass (in MeV/c^2)
virtual G4double getTableMass() const
Get the tabulated particle mass.
double G4double
Definition: G4Types.hh:76
ThreeVector G4INCL::Particle::getTransversePosition ( ) const
inline

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

Definition at line 648 of file G4INCLParticle.hh.

References getLongitudinalPosition(), and thePosition.

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

648  {
650  }
G4INCL::ThreeVector thePosition
ThreeVector getLongitudinalPosition() const
Longitudinal component of the position w.r.t. the momentum.
G4INCL::ParticleType G4INCL::Particle::getType ( ) const
inline
G4int G4INCL::Particle::getZ ( ) const
inline
void G4INCL::Particle::incrementNumberOfCollisions ( )
inline

Increment the number of collisions undergone by the particle.

Definition at line 620 of file G4INCLParticle.hh.

References nCollisions.

620 { nCollisions++; }
void G4INCL::Particle::incrementNumberOfDecays ( )
inline

Increment the number of decays undergone by the particle.

Definition at line 629 of file G4INCLParticle.hh.

References nDecays.

629 { nDecays++; }
G4bool G4INCL::Particle::isCluster ( ) const
inline

Definition at line 668 of file G4INCLParticle.hh.

References G4INCL::Composite, and theType.

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

668  {
669  return (theType == Composite);
670  }
G4INCL::ParticleType theType
G4bool G4INCL::Particle::isDelta ( ) const
inline
G4bool G4INCL::Particle::isInList ( ParticleList const &  l) const
inline

Check if the particle belongs to a given list.

Definition at line 664 of file G4INCLParticle.hh.

664  {
665  return (std::find(l.begin(), l.end(), this)!=l.end());
666  }
G4bool G4INCL::Particle::isNucleon ( ) const
inline
G4bool G4INCL::Particle::isOutOfWell ( ) const
inline

Check if the particle is out of its potential well.

Definition at line 642 of file G4INCLParticle.hh.

Referenced by G4INCL::NuclearPotential::INuclearPotential::computePionPotentialEnergy().

642 { return outOfWell; }
G4bool G4INCL::Particle::isParticipant ( ) const
inline
G4bool G4INCL::Particle::isPion ( ) const
inline
G4bool G4INCL::Particle::isProjectileSpectator ( ) const
inline
G4bool G4INCL::Particle::isResonance ( ) const
inline

Is it a resonance?

Definition at line 292 of file G4INCLParticle.hh.

References isDelta().

Referenced by G4INCL::StandardPropagationModel::generateBinaryCollisionAvatar(), G4INCL::SurfaceAvatar::getChannel(), Particle(), G4INCL::CDPP::processOneParticle(), and setType().

292 { return isDelta(); }
G4bool isDelta() const
Is it a Delta?
G4bool G4INCL::Particle::isTargetSpectator ( ) const
inline
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 345 of file G4INCLParticle.hh.

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

345  {
346  const G4double beta2 = aBoostVector.mag2();
347  const G4double gamma = 1.0 / std::sqrt(1.0 - beta2);
348  const ThreeVector theRelativePosition = thePosition - refPos;
349  const ThreeVector transversePosition = theRelativePosition - aBoostVector * (theRelativePosition.dot(aBoostVector) / aBoostVector.mag2());
350  const ThreeVector longitudinalPosition = theRelativePosition - transversePosition;
351 
352  thePosition = refPos + transversePosition + longitudinalPosition / gamma;
353  }
G4INCL::ThreeVector thePosition
double G4double
Definition: G4Types.hh:76
virtual void G4INCL::Particle::makeParticipant ( )
inlinevirtual

Reimplemented in G4INCL::Cluster.

Definition at line 276 of file G4INCLParticle.hh.

References G4INCL::Participant, and theParticipantType.

Referenced by G4INCL::Store::loadParticles(), and G4INCL::Cluster::makeParticipant().

276  {
278  }
ParticipantType theParticipantType
virtual void G4INCL::Particle::makeProjectileSpectator ( )
inlinevirtual
virtual void G4INCL::Particle::makeTargetSpectator ( )
inlinevirtual

Reimplemented in G4INCL::Cluster.

Definition at line 280 of file G4INCLParticle.hh.

References G4INCL::TargetSpectator, and theParticipantType.

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

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

Assignment operator.

Does not copy the particle ID.

Definition at line 177 of file G4INCLParticle.hh.

References swap().

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

177  {
178  Particle temporaryParticle(rhs);
179  swap(temporaryParticle);
180  return *this;
181  }
void swap(Particle &rhs)
Helper method for the assignment operator.
std::string G4INCL::Particle::print ( ) const
inline

Definition at line 720 of file G4INCLParticle.hh.

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

Referenced by adjustMomentumFromEnergy(), G4INCL::Nucleus::emitInsidePions(), G4INCL::BinaryCollisionAvatar::getChannel(), G4INCL::ParticleEntryChannel::getFinalState(), G4INCL::KinematicsUtils::getLocalEnergy(), G4INCL::StandardPropagationModel::getReflectionTime(), and G4INCL::ProjectileRemnant::removeParticle().

720  {
721  std::stringstream ss;
722  ss << "Particle (ID = " << ID << ") type = ";
724  ss << std::endl
725  << " energy = " << theEnergy << std::endl
726  << " momentum = "
727  << theMomentum.print()
728  << std::endl
729  << " position = "
730  << thePosition.print()
731  << std::endl;
732  return ss.str();
733  };
G4INCL::ThreeVector thePosition
G4INCL::ParticleType theType
G4INCL::ThreeVector theMomentum
std::string print() const
std::string getName(const ParticleType t)
Get the native INCL name of the particle.
void G4INCL::Particle::propagate ( G4double  step)
inline

Definition at line 609 of file G4INCLParticle.hh.

References thePosition.

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

609  {
610  thePosition += ((*thePropagationMomentum)*(step/(*thePropagationEnergy)));
611  };
G4INCL::ThreeVector thePosition
virtual void G4INCL::Particle::rotate ( const G4double  angle,
const ThreeVector axis 
)
inlinevirtual

Rotate the particle position and momentum.

Parameters
anglethe rotation angle
axisa unit vector representing the rotation axis

Reimplemented in G4INCL::Cluster.

Definition at line 714 of file G4INCLParticle.hh.

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

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

714  {
715  thePosition.rotate(angle, axis);
716  theMomentum.rotate(angle, axis);
717  theFrozenMomentum.rotate(angle, axis);
718  }
void rotate(const G4double angle, const ThreeVector &axis)
Rotate the vector by a given angle around a given axis.
G4INCL::ThreeVector theFrozenMomentum
G4INCL::ThreeVector thePosition
G4INCL::ThreeVector theMomentum
void G4INCL::Particle::rpCorrelate ( )
inline

Make the particle follow a strict r-p correlation.

Definition at line 775 of file G4INCLParticle.hh.

References rpCorrelated.

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

775 { rpCorrelated = true; }
void G4INCL::Particle::rpDecorrelate ( )
inline

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

Definition at line 778 of file G4INCLParticle.hh.

References rpCorrelated.

778 { rpCorrelated = false; }
void G4INCL::Particle::setEmissionTime ( G4double  t)
inline
void G4INCL::Particle::setEnergy ( G4double  energy)
inline
void G4INCL::Particle::setFrozenEnergy ( const G4double  energy)
inline

Set the frozen particle momentum.

Definition at line 676 of file G4INCLParticle.hh.

References energy(), and theFrozenEnergy.

676 { theFrozenEnergy = energy; }
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
G4double theFrozenEnergy
void G4INCL::Particle::setFrozenMomentum ( const ThreeVector momentum)
inline

Set the frozen particle momentum.

Definition at line 673 of file G4INCLParticle.hh.

References theFrozenMomentum.

673 { theFrozenMomentum = momentum; }
G4INCL::ThreeVector theFrozenMomentum
void G4INCL::Particle::setHelicity ( G4double  h)
inline

Definition at line 607 of file G4INCLParticle.hh.

Referenced by G4INCL::DeltaProductionChannel::getFinalState(), and G4INCL::DeltaDecayChannel::getFinalState().

607 { theHelicity = h; };
void G4INCL::Particle::setINCLMass ( )
inline

Set the mass of the Particle to its table mass.

Definition at line 452 of file G4INCLParticle.hh.

References getINCLMass(), and setMass().

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

452 { setMass(getINCLMass()); }
void setMass(G4double mass)
G4double getINCLMass() const
Get the INCL particle mass.
void G4INCL::Particle::setMass ( G4double  mass)
inline
virtual void G4INCL::Particle::setMomentum ( const G4INCL::ThreeVector momentum)
inlinevirtual
void G4INCL::Particle::setNumberOfCollisions ( G4int  n)
inline

Set the number of collisions undergone by the particle.

Definition at line 617 of file G4INCLParticle.hh.

References n, and nCollisions.

617 { nCollisions = n; }
const G4int n
void G4INCL::Particle::setNumberOfDecays ( G4int  n)
inline

Set the number of decays undergone by the particle.

Definition at line 626 of file G4INCLParticle.hh.

References n, and nDecays.

626 { nDecays = n; }
const G4int n
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 639 of file G4INCLParticle.hh.

639 { outOfWell = true; }
void G4INCL::Particle::setParticipantType ( ParticipantType const  p)
inline

Definition at line 260 of file G4INCLParticle.hh.

References theParticipantType.

260  {
262  }
ParticipantType theParticipantType
const char * p
Definition: xmltok.h:285
virtual void G4INCL::Particle::setPosition ( const G4INCL::ThreeVector position)
inlinevirtual
void G4INCL::Particle::setPotentialEnergy ( G4double  v)
inline

Set the particle potential energy.

Definition at line 545 of file G4INCLParticle.hh.

References thePotentialEnergy, and test::v.

Referenced by G4INCL::Nucleus::emitInsidePions(), G4INCL::Store::loadParticles(), and G4INCL::Nucleus::updatePotentialEnergy().

545 { thePotentialEnergy = v; }
G4double thePotentialEnergy
void G4INCL::Particle::setRealMass ( )
inline

Set the mass of the Particle to its real mass.

Definition at line 446 of file G4INCLParticle.hh.

References getRealMass(), and setMass().

Referenced by G4INCL::ClusterDecay::decay().

446 { setMass(getRealMass()); }
void setMass(G4double mass)
G4double getRealMass() const
Get the real particle mass.
void G4INCL::Particle::setTableMass ( )
inline

Set the mass of the Particle to its table mass.

Definition at line 449 of file G4INCLParticle.hh.

References getTableMass(), and setMass().

Referenced by G4INCL::Nucleus::decayOutgoingDeltas(), G4INCL::Nucleus::emitInsidePions(), G4INCL::ProjectileRemnant::ProjectileRemnant(), and G4INCL::ProjectileRemnant::reset().

449 { setMass(getTableMass()); }
void setMass(G4double mass)
virtual G4double getTableMass() const
Get the tabulated particle mass.
void G4INCL::Particle::setType ( ParticleType  t)
inline

Definition at line 196 of file G4INCLParticle.hh.

References G4INCL::Composite, G4INCL::DeltaMinus, G4INCL::DeltaPlus, G4INCL::DeltaPlusPlus, G4INCL::DeltaZero, INCL_ERROR, isResonance(), G4INCL::Neutron, G4INCL::PiMinus, G4INCL::PiPlus, G4INCL::PiZero, G4INCL::Proton, setINCLMass(), theA, theType, theZ, and G4INCL::UnknownParticle.

Referenced by G4INCL::Cluster::Cluster(), G4INCL::ClusterDecay::decay(), G4INCL::DeltaProductionChannel::getFinalState(), G4INCL::PionNucleonChannel::getFinalState(), G4INCL::ElasticChannel::getFinalState(), G4INCL::DeltaDecayChannel::getFinalState(), G4INCL::RecombinationChannel::getFinalState(), and Particle().

196  {
197  theType = t;
198  switch(theType)
199  {
200  case DeltaPlusPlus:
201  theA = 1;
202  theZ = 2;
203  break;
204  case Proton:
205  case DeltaPlus:
206  theA = 1;
207  theZ = 1;
208  break;
209  case Neutron:
210  case DeltaZero:
211  theA = 1;
212  theZ = 0;
213  break;
214  case DeltaMinus:
215  theA = 1;
216  theZ = -1;
217  break;
218  case PiPlus:
219  theA = 0;
220  theZ = 1;
221  break;
222  case PiZero:
223  theA = 0;
224  theZ = 0;
225  break;
226  case PiMinus:
227  theA = 0;
228  theZ = -1;
229  break;
230  case Composite:
231  // INCL_ERROR("Trying to set particle type to Composite! Construct a Cluster object instead" << std::endl);
232  theA = 0;
233  theZ = 0;
234  break;
235  case UnknownParticle:
236  theA = 0;
237  theZ = 0;
238  INCL_ERROR("Trying to set particle type to Unknown!" << std::endl);
239  break;
240  }
241 
242  if( !isResonance() && t!=Composite )
243  setINCLMass();
244  }
G4bool isResonance() const
Is it a resonance?
#define INCL_ERROR(x)
void setINCLMass()
Set the mass of the Particle to its table mass.
G4INCL::ParticleType theType
void G4INCL::Particle::setUncorrelatedMomentum ( const G4double  p)
inline

Set the uncorrelated momentum.

Definition at line 772 of file G4INCLParticle.hh.

References uncorrelatedMomentum.

772 { uncorrelatedMomentum = p; }
const char * p
Definition: xmltok.h:285
G4double uncorrelatedMomentum
void G4INCL::Particle::swap ( Particle rhs)
inlineprotected

Helper method for the assignment operator.

Definition at line 139 of file G4INCLParticle.hh.

References nCollisions, nDecays, rpCorrelated, CLHEP::swap(), theA, theEnergy, theFrozenEnergy, theFrozenMomentum, theMomentum, theParticipantType, thePosition, thePotentialEnergy, thePropagationEnergy, thePropagationMomentum, theType, theZ, and uncorrelatedMomentum.

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

139  {
140  std::swap(theZ, rhs.theZ);
141  std::swap(theA, rhs.theA);
142  std::swap(theParticipantType, rhs.theParticipantType);
143  std::swap(theType, rhs.theType);
144  if(rhs.thePropagationEnergy == &(rhs.theFrozenEnergy))
146  else
148  std::swap(theEnergy, rhs.theEnergy);
149  std::swap(theFrozenEnergy, rhs.theFrozenEnergy);
150  if(rhs.thePropagationMomentum == &(rhs.theFrozenMomentum))
152  else
154  std::swap(theMomentum, rhs.theMomentum);
155  std::swap(theFrozenMomentum, rhs.theFrozenMomentum);
156  std::swap(thePosition, rhs.thePosition);
157  std::swap(nCollisions, rhs.nCollisions);
158  std::swap(nDecays, rhs.nDecays);
159  std::swap(thePotentialEnergy, rhs.thePotentialEnergy);
160  // ID intentionally not swapped
161 
162  std::swap(theHelicity, rhs.theHelicity);
163  std::swap(emissionTime, rhs.emissionTime);
164  std::swap(outOfWell, rhs.outOfWell);
165 
166  std::swap(theMass, rhs.theMass);
167  std::swap(rpCorrelated, rhs.rpCorrelated);
168  std::swap(uncorrelatedMomentum, rhs.uncorrelatedMomentum);
169  }
ParticipantType theParticipantType
G4double * thePropagationEnergy
G4INCL::ThreeVector theFrozenMomentum
G4double uncorrelatedMomentum
void swap(shared_ptr< P > &, shared_ptr< P > &)
Definition: memory.h:1247
G4INCL::ThreeVector thePosition
G4INCL::ParticleType theType
G4INCL::ThreeVector theMomentum
G4double theFrozenEnergy
G4INCL::ThreeVector * thePropagationMomentum
G4double thePotentialEnergy
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 704 of file G4INCLParticle.hh.

References theEnergy, theMomentum, thePropagationEnergy, and thePropagationMomentum.

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

704  {
707  }
G4double * thePropagationEnergy
G4INCL::ThreeVector theMomentum
G4INCL::ThreeVector * thePropagationMomentum

Field Documentation

long G4INCL::Particle::ID
protected

Definition at line 803 of file G4INCLParticle.hh.

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

G4int G4INCL::Particle::nCollisions
protected
G4int G4INCL::Particle::nDecays
protected
G4bool G4INCL::Particle::rpCorrelated
protected

Definition at line 805 of file G4INCLParticle.hh.

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

G4int G4INCL::Particle::theA
protected
G4double G4INCL::Particle::theEnergy
protected
G4double G4INCL::Particle::theFrozenEnergy
protected
G4INCL::ThreeVector G4INCL::Particle::theFrozenMomentum
protected
G4INCL::ThreeVector G4INCL::Particle::theMomentum
protected
ParticipantType G4INCL::Particle::theParticipantType
protected
G4INCL::ThreeVector G4INCL::Particle::thePosition
protected
G4double G4INCL::Particle::thePotentialEnergy
protected
G4double* G4INCL::Particle::thePropagationEnergy
protected

Definition at line 794 of file G4INCLParticle.hh.

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

G4INCL::ThreeVector* G4INCL::Particle::thePropagationMomentum
protected
G4INCL::ParticleType G4INCL::Particle::theType
protected
G4int G4INCL::Particle::theZ
protected
G4double G4INCL::Particle::uncorrelatedMomentum
protected

Definition at line 806 of file G4INCLParticle.hh.

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


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