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 momentum, ThreeVector position)
 Particle (ParticleType t, ThreeVector momentum, ThreeVector position)
virtual ~Particle ()
 Particle (const Particle &rhs)
 Copy constructor.
Particleoperator= (const Particle &rhs)
 Assignment operator.
G4INCL::ParticleType getType () const
virtual G4INCL::ParticleSpecies getSpecies () const
 Get the particle species.
void setType (ParticleType t)
G4bool isNucleon () const
ParticipantType getParticipantType () const
void setParticipantType (ParticipantType const p)
G4bool isParticipant () const
G4bool isTargetSpectator () const
G4bool isProjectileSpectator () const
virtual void makeParticipant ()
virtual void makeTargetSpectator ()
virtual void makeProjectileSpectator ()
G4bool isPion () const
 Is this a pion?
G4bool isResonance () const
 Is it a resonance?
G4bool isDelta () const
 Is it a Delta?
G4int getA () const
 Returns the baryon number.
G4int getZ () const
 Returns the charge number.
G4double getBeta () const
ThreeVector boostVector () const
void boost (const ThreeVector &aBoostVector)
void lorentzContract (const ThreeVector &aBoostVector, const ThreeVector &refPos)
 Lorentz-contract the particle position around some center.
G4double getMass () const
 Get the cached particle mass.
G4double getINCLMass () const
 Get the INCL particle mass.
virtual G4double getTableMass () const
 Get the tabulated particle mass.
G4double getRealMass () const
 Get the real particle mass.
void setRealMass ()
 Set the mass of the Particle to its real mass.
void setTableMass ()
 Set the mass of the Particle to its table mass.
void setINCLMass ()
 Set the mass of the Particle to its table mass.
G4double getEmissionQValueCorrection (const G4int AParent, const G4int ZParent) const
 Computes correction on the emission Q-value.
G4double getTransferQValueCorrection (const G4int AFrom, const G4int ZFrom, const G4int ATo, const G4int ZTo) const
 Computes correction on the transfer Q-value.
G4double getInvariantMass () const
 Get the the particle invariant mass.
G4double getKineticEnergy () const
 Get the particle kinetic energy.
G4double getPotentialEnergy () const
 Get the particle potential energy.
void setPotentialEnergy (G4double v)
 Set the particle potential energy.
G4double getEnergy () const
void setMass (G4double mass)
void setEnergy (G4double energy)
const G4INCL::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.
void setNumberOfCollisions (G4int n)
 Set the number of collisions undergone by the particle.
void incrementNumberOfCollisions ()
 Increment the number of collisions undergone by the particle.
G4int getNumberOfDecays () const
 Return the number of decays undergone by the particle.
void setNumberOfDecays (G4int n)
 Set the number of decays undergone by the particle.
void incrementNumberOfDecays ()
 Increment the number of decays undergone by the particle.
void setOutOfWell ()
 Mark the particle as out of its potential well.
G4bool isOutOfWell () const
 Check if the particle is out of its potential well.
void setEmissionTime (G4double t)
G4double getEmissionTime ()
ThreeVector getTransversePosition () const
 Transverse component of the position w.r.t. the momentum.
ThreeVector getLongitudinalPosition () const
 Longitudinal component of the position w.r.t. the momentum.
const ThreeVectoradjustMomentumFromEnergy ()
 Rescale the momentum to match the total energy.
G4double adjustEnergyFromMomentum ()
 Recompute the energy to match the momentum.
G4bool isInList (ParticleList const &l) const
 Check if the particle belongs to a given list.
G4bool isCluster () const
void setFrozenMomentum (const ThreeVector &momentum)
 Set the frozen particle momentum.
void setFrozenEnergy (const G4double energy)
 Set the frozen particle momentum.
ThreeVector getFrozenMomentum () const
 Get the frozen particle momentum.
G4double getFrozenEnergy () const
 Get the frozen particle momentum.
ThreeVector getPropagationVelocity () const
 Get the propagation velocity of the particle.
void freezePropagation ()
 Freeze particle propagation.
void thawPropagation ()
 Unfreeze particle propagation.
virtual void rotate (const G4double angle, const ThreeVector &axis)
 Rotate the particle position and momentum.
std::string print () const
std::string dump () const
long getID () const
ParticleList const * getParticles () const

Protected Member Functions

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

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

Detailed Description

Definition at line 64 of file G4INCLParticle.hh.


Constructor & Destructor Documentation

G4INCL::Particle::Particle (  ) 

Definition at line 51 of file G4INCLParticle.cc.

References ID.

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

00052     : theZ(0), theA(0),
00053     theParticipantType(TargetSpectator),
00054     theType(UnknownParticle),
00055     theEnergy(0.0),
00056     thePropagationEnergy(&theEnergy),
00057     theFrozenEnergy(theEnergy),
00058     theMomentum(ThreeVector(0.,0.,0.)),
00059     thePropagationMomentum(&theMomentum),
00060     theFrozenMomentum(theMomentum),
00061     thePosition(ThreeVector(0.,0.,0.)),
00062     nCollisions(0),
00063     nDecays(0),
00064     thePotentialEnergy(0.0),
00065     theHelicity(0.0),
00066     emissionTime(0.0),
00067     outOfWell(false),
00068     theMass(0.)
00069   {
00070     ID = nextID;
00071     nextID++;
00072   }

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

Definition at line 74 of file G4INCLParticle.cc.

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

00076     : theEnergy(energy),
00077     thePropagationEnergy(&theEnergy),
00078     theFrozenEnergy(theEnergy),
00079     theMomentum(momentum),
00080     thePropagationMomentum(&theMomentum),
00081     theFrozenMomentum(theMomentum),
00082     thePosition(position),
00083     nCollisions(0), nDecays(0),
00084       thePotentialEnergy(0.), theHelicity(0.0),
00085       emissionTime(0.0), outOfWell(false)
00086   {
00087     theParticipantType = TargetSpectator;
00088     ID = nextID;
00089     nextID++;
00090     if(theEnergy <= 0.0) {
00091       WARN("Particle with energy " << theEnergy << " created." << std::endl);
00092     }
00093     setType(t);
00094     setMass(getInvariantMass());
00095   }

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

Definition at line 97 of file G4INCLParticle.cc.

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

00099     : thePropagationEnergy(&theEnergy),
00100     theMomentum(momentum),
00101     thePropagationMomentum(&theMomentum),
00102     theFrozenMomentum(theMomentum),
00103     thePosition(position),
00104     nCollisions(0), nDecays(0),
00105       thePotentialEnergy(0.), theHelicity(0.0),
00106       emissionTime(0.0), outOfWell(false)
00107   {
00108     theParticipantType = TargetSpectator;
00109     ID = nextID;
00110     nextID++;
00111     setType(t);
00112     if( isResonance() ) {
00113       ERROR("Cannot create resonance without specifying its momentum four-vector." << std::endl);
00114     }
00115     G4double energy = std::sqrt(theMomentum.mag2() + theMass*theMass);
00116     theEnergy = energy;
00117     theFrozenEnergy = theEnergy;
00118   }

virtual G4INCL::Particle::~Particle (  )  [inline, virtual]

Definition at line 69 of file G4INCLParticle.hh.

00069 {}

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

Copy constructor.

Does not copy the particle ID.

Definition at line 75 of file G4INCLParticle.hh.

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

00075                                   :
00076       theZ(rhs.theZ),
00077       theA(rhs.theA),
00078       theParticipantType(rhs.theParticipantType),
00079       theType(rhs.theType),
00080       theEnergy(rhs.theEnergy),
00081       theFrozenEnergy(rhs.theFrozenEnergy),
00082       theMomentum(rhs.theMomentum),
00083       theFrozenMomentum(rhs.theFrozenMomentum),
00084       thePosition(rhs.thePosition),
00085       nCollisions(rhs.nCollisions),
00086       nDecays(rhs.nDecays),
00087       thePotentialEnergy(rhs.thePotentialEnergy),
00088       theHelicity(rhs.theHelicity),
00089       emissionTime(rhs.emissionTime),
00090       outOfWell(rhs.outOfWell),
00091       theMass(rhs.theMass)
00092       {
00093         if(rhs.thePropagationEnergy == &(rhs.theFrozenEnergy))
00094           thePropagationEnergy = &theFrozenEnergy;
00095         else
00096           thePropagationEnergy = &theEnergy;
00097         if(rhs.thePropagationMomentum == &(rhs.theFrozenMomentum))
00098           thePropagationMomentum = &theFrozenMomentum;
00099         else
00100           thePropagationMomentum = &theMomentum;
00101         // ID intentionally not copied
00102         ID = nextID++;
00103       }


Member Function Documentation

G4double G4INCL::Particle::adjustEnergyFromMomentum (  ) 

Recompute the energy to match the momentum.

Definition at line 133 of file G4INCLParticle.cc.

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

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

00133                                               {
00134     theEnergy = std::sqrt(theMomentum.mag2() + theMass*theMass);
00135     return theEnergy;
00136   }

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

Rescale the momentum to match the total energy.

Definition at line 120 of file G4INCLParticle.cc.

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

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

00120                                                         {
00121     const G4double p2 = theMomentum.mag2();
00122     G4double newp2 = theEnergy*theEnergy - theMass*theMass;
00123     if( newp2<0.0 ) {
00124       ERROR("Particle has E^2 < m^2." << std::endl << print());
00125       newp2 = 0.0;
00126       theEnergy = theMass;
00127     }
00128 
00129     theMomentum *= std::sqrt(newp2/p2);
00130     return theMomentum;
00131   }

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());

Reimplemented in G4INCL::Cluster.

Definition at line 293 of file G4INCLParticle.hh.

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

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

00293                                                 {
00294       const G4double beta2 = aBoostVector.mag2();
00295       const G4double gamma = 1.0 / std::sqrt(1.0 - beta2);
00296       const G4double bp = theMomentum.dot(aBoostVector);
00297       const G4double alpha = (gamma*gamma)/(1.0 + gamma);
00298 
00299       theMomentum = theMomentum + aBoostVector * (alpha * bp - gamma * theEnergy);
00300       theEnergy = gamma * (theEnergy - bp);
00301     }

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

References theEnergy, and theMomentum.

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

00283                                     {
00284       return theMomentum / theEnergy;
00285     }

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

Definition at line 701 of file G4INCLParticle.hh.

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

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

00701                            {
00702       std::stringstream ss;
00703       ss << "(particle " << ID << " ";
00704       ss << ParticleTable::getName(theType);
00705       ss << std::endl
00706         << thePosition.dump()
00707         << std::endl
00708         << theMomentum.dump()
00709         << std::endl
00710         << theEnergy << ")" << std::endl;
00711       return ss.str();
00712     };

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

References theFrozenEnergy, theFrozenMomentum, thePropagationEnergy, and thePropagationMomentum.

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

Returns the baryon number.

Definition at line 267 of file G4INCLParticle.hh.

References theA.

Referenced by G4INCL::Cluster::addParticle(), G4INCL::ClusteringModelIntercomparison::clusterCanEscape(), G4INCL::ClusterDecay::decay(), G4INCL::Nucleus::fillEventInfo(), G4INCL::Nucleus::finalizeProjectileRemnant(), G4INCL::ClusterUtils::getA(), G4INCL::ClusteringModelIntercomparison::getCluster(), G4INCL::Nucleus::getConservationBalance(), G4INCL::ParticleEntryChannel::getFinalState(), G4INCL::ClusterUtils::getPhaseSpace(), G4INCL::SurfaceAvatar::getTransmissionProbability(), G4INCL::NuclearDensity::getTransmissionRadius(), G4INCL::Nucleus::insertParticle(), G4INCL::Nucleus::isEventTransparent(), G4INCL::ClusterDecay::isStable(), G4INCL::ProjectileRemnant::removeParticle(), G4INCL::StandardPropagationModel::shootComposite(), and G4INCL::StandardPropagationModel::shootParticle().

00267 { return theA; }

virtual G4INCL::ThreeVector G4INCL::Particle::getAngularMomentum (  )  const [inline, virtual]

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

Reimplemented in G4INCL::Cluster.

Definition at line 546 of file G4INCLParticle.hh.

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

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

00547     {
00548       return thePosition.vector(theMomentum);
00549     };

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

Definition at line 272 of file G4INCLParticle.hh.

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

00272                              {
00273       const G4double P = theMomentum.mag();
00274       return P/theEnergy;
00275     }

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:
AParent the mass number of the emitting nucleus
ZParent the charge number of the emitting nucleus
Returns:
the correction

Definition at line 431 of file G4INCLParticle.hh.

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

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

00431                                                                                          {
00432       const G4int ADaughter = AParent - theA;
00433       const G4int ZDaughter = ZParent - theZ;
00434 
00435       // Note the minus sign here
00436       G4double theQValue;
00437       if(isCluster())
00438         theQValue = -ParticleTable::getTableQValue(theA, theZ, ADaughter, ZDaughter);
00439       else {
00440         const G4double massTableParent = ParticleTable::getTableMass(AParent,ZParent);
00441         const G4double massTableDaughter = ParticleTable::getTableMass(ADaughter,ZDaughter);
00442         const G4double massTableParticle = getTableMass();
00443         theQValue = massTableParent - massTableDaughter - massTableParticle;
00444       }
00445 
00446       const G4double massINCLParent = ParticleTable::getINCLMass(AParent,ZParent);
00447       const G4double massINCLDaughter = ParticleTable::getINCLMass(ADaughter,ZDaughter);
00448       const G4double massINCLParticle = getINCLMass();
00449 
00450       // The rhs corresponds to the INCL Q-value
00451       return theQValue - (massINCLParent-massINCLDaughter-massINCLParticle);
00452     }

G4double G4INCL::Particle::getEmissionTime (  )  [inline]

Definition at line 611 of file G4INCLParticle.hh.

00611 { return emissionTime; };

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

Get the energy of the particle in MeV.

Definition at line 516 of file G4INCLParticle.hh.

References theEnergy.

Referenced by G4INCL::Cluster::addParticle(), G4INCL::DeltaDecayChannel::computeDecayTime(), G4INCL::ClusteringModelIntercomparison::getCluster(), G4INCL::TransmissionChannel::getFinalState(), G4INCL::ReflectionChannel::getFinalState(), G4INCL::PionNucleonChannel::getFinalState(), G4INCL::ParticleEntryChannel::getFinalState(), G4INCL::KinematicsUtils::makeBoostVector(), G4INCL::KinematicsUtils::momentumInCM(), G4INCL::InteractionAvatar::preInteraction(), G4INCL::InteractionAvatar::preInteractionBlocking(), G4INCL::ProjectileRemnant::removeParticle(), G4INCL::StandardPropagationModel::shootParticle(), G4INCL::KinematicsUtils::squareTotalEnergyInCM(), and G4INCL::KinematicsUtils::transformToLocalEnergyFrame().

00517     {
00518       return theEnergy;
00519     };

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

Get the frozen particle momentum.

Definition at line 648 of file G4INCLParticle.hh.

References theFrozenEnergy.

00648 { return theFrozenEnergy; }

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

Get the frozen particle momentum.

Definition at line 645 of file G4INCLParticle.hh.

References theFrozenMomentum.

00645 { return theFrozenMomentum; }

G4double G4INCL::Particle::getHelicity (  )  [inline]

Definition at line 572 of file G4INCLParticle.hh.

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

00572 { return theHelicity; };

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

Definition at line 714 of file G4INCLParticle.hh.

References ID.

Referenced by G4INCL::Store::add(), G4INCL::PauliStandard::getBlockingProbability(), G4INCL::SurfaceAvatar::getChannel(), G4INCL::ClusteringModelIntercomparison::getCluster(), G4INCL::ReflectionChannel::getFinalState(), and G4INCL::ParticleEntryChannel::getFinalState().

00714 { return ID; };

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

Get the INCL particle mass.

Definition at line 325 of file G4INCLParticle.hh.

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

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

00325                                         {
00326       switch(theType) {
00327         case Proton:
00328         case Neutron:
00329         case PiPlus:
00330         case PiMinus:
00331         case PiZero:
00332           return ParticleTable::getINCLMass(theType);
00333           break;
00334 
00335         case DeltaPlusPlus:
00336         case DeltaPlus:
00337         case DeltaZero:
00338         case DeltaMinus:
00339           return theMass;
00340           break;
00341 
00342         case Composite:
00343           return ParticleTable::getINCLMass(theA,theZ);
00344           break;
00345 
00346         default:
00347           ERROR("Particle::getINCLMass: Unknown particle type." << std::endl);
00348           return 0.0;
00349           break;
00350       }
00351     }

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

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

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

00494                                       {
00495       const G4double mass = std::pow(theEnergy, 2) - theMomentum.dot(theMomentum);
00496       if(mass < 0.0) {
00497         ERROR("E*E - p*p is negative." << std::endl);
00498         return 0.0;
00499       } else {
00500         return std::sqrt(mass);
00501       }
00502     };

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

Get the particle kinetic energy.

Definition at line 505 of file G4INCLParticle.hh.

References theEnergy.

Referenced by G4INCL::NuclearPotential::NuclearPotentialEnergyIsospinSmooth::computePotentialEnergy(), G4INCL::NuclearPotential::NuclearPotentialEnergyIsospin::computePotentialEnergy(), G4INCL::Nucleus::fillEventInfo(), G4INCL::SurfaceAvatar::getChannel(), G4INCL::Nucleus::getConservationBalance(), G4INCL::ParticleEntryChannel::getFinalState(), G4INCL::KinematicsUtils::getLocalEnergy(), G4INCL::SurfaceAvatar::getTransmissionProbability(), and G4INCL::CoulombDistortion::maxImpactParameter().

00505 { return theEnergy - theMass; }

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

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

Definition at line 619 of file G4INCLParticle.hh.

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

Referenced by getTransversePosition().

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

Get the cached particle mass.

Definition at line 322 of file G4INCLParticle.hh.

Referenced by G4INCL::DeltaDecayChannel::computeDecayTime(), G4INCL::InteractionAvatar::enforceEnergyConservation(), G4INCL::Cluster::freezeInternalMotion(), G4INCL::NuclearPotential::INuclearPotential::getFermiMomentum(), G4INCL::TransmissionChannel::getFinalState(), G4INCL::RecombinationChannel::getFinalState(), G4INCL::DeltaDecayChannel::getFinalState(), G4INCL::KinematicsUtils::getLocalEnergy(), G4INCL::SurfaceAvatar::getTransmissionProbability(), G4INCL::CrossSections::interactionDistanceNN(), G4INCL::CrossSections::interactionDistancePiN(), G4INCL::Cluster::internalBoostToCM(), G4INCL::KinematicsUtils::momentumInCM(), G4INCL::KinematicsUtils::momentumInLab(), G4INCL::InteractionAvatar::preInteractionBlocking(), G4INCL::Cluster::print(), G4INCL::ProjectileRemnant::ProjectileRemnant(), G4INCL::CrossSections::recombination(), and G4INCL::StandardPropagationModel::shootParticle().

00322 { return theMass; }

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

Get the momentum vector.

Definition at line 540 of file G4INCLParticle.hh.

References theMomentum.

Referenced by G4INCL::Cluster::addParticle(), G4INCL::ClusteringModelIntercomparison::clusterCanEscape(), G4INCL::Nucleus::decayOutgoingDeltas(), G4INCL::Nucleus::fillEventInfo(), G4INCL::PauliStandard::getBlockingProbability(), G4INCL::ClusteringModelIntercomparison::getCluster(), G4INCL::Nucleus::getConservationBalance(), G4INCL::ReflectionChannel::getFinalState(), G4INCL::PionNucleonChannel::getFinalState(), G4INCL::ParticleEntryChannel::getFinalState(), G4INCL::ElasticChannel::getFinalState(), G4INCL::DeltaProductionChannel::getFinalState(), G4INCL::ClusterUtils::getPhaseSpace(), G4INCL::Nucleus::getSurfaceRadius(), G4INCL::KinematicsUtils::makeBoostVector(), G4INCL::KinematicsUtils::momentumInCM(), G4INCL::InteractionAvatar::preInteraction(), G4INCL::InteractionAvatar::preInteractionBlocking(), G4INCL::ProjectileRemnant::removeParticle(), G4INCL::ParticleSampler::sampleParticles(), and G4INCL::StandardPropagationModel::shootParticle().

00541     {
00542       return theMomentum;
00543     };

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

Return the number of collisions undergone by the particle.

Definition at line 580 of file G4INCLParticle.hh.

References nCollisions.

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

00580 { return nCollisions; }

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

Return the number of decays undergone by the particle.

Definition at line 589 of file G4INCLParticle.hh.

References nDecays.

00589 { return nDecays; }

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

Definition at line 222 of file G4INCLParticle.hh.

References theParticipantType.

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

00222                                                {
00223       return theParticipantType;
00224     }

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

Return a NULL pointer

Reimplemented in G4INCL::Cluster.

Definition at line 719 of file G4INCLParticle.hh.

References WARN.

00719                                              {
00720       WARN("Particle::getParticles() method was called on a Particle object" << std::endl);
00721       return 0;
00722     }

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

Set the position vector.

Definition at line 562 of file G4INCLParticle.hh.

References thePosition.

Referenced by G4INCL::Cluster::addParticle(), G4INCL::InteractionAvatar::bringParticleInside(), G4INCL::CoulombNone::bringToSurface(), G4INCL::ClusteringModelIntercomparison::clusterCanEscape(), G4INCL::StandardPropagationModel::generateBinaryCollisionAvatar(), G4INCL::PauliStandard::getBlockingProbability(), G4INCL::BinaryCollisionAvatar::getChannel(), G4INCL::ClusteringModelIntercomparison::getCluster(), G4INCL::ReflectionChannel::getFinalState(), G4INCL::DeltaDecayChannel::getFinalState(), G4INCL::KinematicsUtils::getLocalEnergy(), G4INCL::ClusterUtils::getNewPositionVector(), G4INCL::ClusterUtils::getPhaseSpace(), G4INCL::StandardPropagationModel::getReflectionTime(), G4INCL::StandardPropagationModel::getTime(), G4INCL::InteractionAvatar::preInteractionBlocking(), and G4INCL::ParticleSampler::sampleParticles().

00563     {
00564       return thePosition;
00565     };

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

Get the particle potential energy.

Definition at line 508 of file G4INCLParticle.hh.

References thePotentialEnergy.

Referenced by G4INCL::Cluster::addParticle(), G4INCL::ClusteringModelIntercomparison::getCluster(), G4INCL::TransmissionChannel::getFinalState(), G4INCL::ReflectionChannel::getFinalState(), G4INCL::PionNucleonChannel::getFinalState(), G4INCL::KinematicsUtils::getLocalEnergy(), G4INCL::SurfaceAvatar::getTransmissionProbability(), and G4INCL::InteractionAvatar::preInteractionBlocking().

00508 { return thePotentialEnergy; }

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

Get the propagation velocity of the particle.

Definition at line 651 of file G4INCLParticle.hh.

References thePropagationMomentum.

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

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

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

Get the real particle mass.

Definition at line 383 of file G4INCLParticle.hh.

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

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

00383                                         {
00384       switch(theType) {
00385         case Proton:
00386         case Neutron:
00387         case PiPlus:
00388         case PiMinus:
00389         case PiZero:
00390           return ParticleTable::getRealMass(theType);
00391           break;
00392 
00393         case DeltaPlusPlus:
00394         case DeltaPlus:
00395         case DeltaZero:
00396         case DeltaMinus:
00397           return theMass;
00398           break;
00399 
00400         case Composite:
00401           return ParticleTable::getRealMass(theA,theZ);
00402           break;
00403 
00404         default:
00405           ERROR("Particle::getRealMass: Unknown particle type." << std::endl);
00406           return 0.0;
00407           break;
00408       }
00409     }

virtual G4INCL::ParticleSpecies G4INCL::Particle::getSpecies (  )  const [inline, virtual]

Get the particle species.

Reimplemented in G4INCL::Cluster.

Definition at line 158 of file G4INCLParticle.hh.

References theType.

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

00158                                                      {
00159       return ParticleSpecies(theType);
00160     };

virtual G4double G4INCL::Particle::getTableMass (  )  const [inline, virtual]

Get the tabulated particle mass.

Reimplemented in G4INCL::Cluster.

Definition at line 354 of file G4INCLParticle.hh.

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

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

00354                                                  {
00355       switch(theType) {
00356         case Proton:
00357         case Neutron:
00358         case PiPlus:
00359         case PiMinus:
00360         case PiZero:
00361           return ParticleTable::getTableParticleMass(theType);
00362           break;
00363 
00364         case DeltaPlusPlus:
00365         case DeltaPlus:
00366         case DeltaZero:
00367         case DeltaMinus:
00368           return theMass;
00369           break;
00370 
00371         case Composite:
00372           return ParticleTable::getTableMass(theA,theZ);
00373           break;
00374 
00375         default:
00376           ERROR("Particle::getTableMass: Unknown particle type." << std::endl);
00377           return 0.0;
00378           break;
00379       }
00380     }

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:
AFrom the mass number of the donating nucleus
ZFrom the charge number of the donating nucleus
ATo the mass number of the receiving nucleus
ZTo the charge number of the receiving nucleus
Returns:
the correction

Definition at line 469 of file G4INCLParticle.hh.

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

00469                                                                                                                        {
00470       const G4int AFromDaughter = AFrom - theA;
00471       const G4int ZFromDaughter = ZFrom - theZ;
00472       const G4int AToDaughter = ATo + theA;
00473       const G4int ZToDaughter = ZTo + theZ;
00474       const G4double theQValue = ParticleTable::getTableQValue(AToDaughter,ZToDaughter,AFromDaughter,ZFromDaughter,AFrom,ZFrom);
00475 
00476       const G4double massINCLTo = ParticleTable::getINCLMass(ATo,ZTo);
00477       const G4double massINCLToDaughter = ParticleTable::getINCLMass(AToDaughter,ZToDaughter);
00478       /* Note that here we have to use the table mass in the INCL Q-value. We
00479        * cannot use theMass, because at this stage the particle is probably
00480        * still off-shell; and we cannot use getINCLMass(), because it leads to
00481        * violations of global energy conservation.
00482        */
00483       const G4double massINCLParticle = getTableMass();
00484 
00485       // The rhs corresponds to the INCL Q-value for particle absorption
00486       return theQValue - (massINCLToDaughter-massINCLTo-massINCLParticle);
00487     }

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

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

Definition at line 614 of file G4INCLParticle.hh.

References getLongitudinalPosition(), and thePosition.

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

00614                                               {
00615       return thePosition - getLongitudinalPosition();
00616     }

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

Get the particle type.

See also:
G4INCL::ParticleType

Definition at line 153 of file G4INCLParticle.hh.

References theType.

Referenced by G4INCL::NuclearPotential::INuclearPotential::computePionPotentialEnergy(), G4INCL::NuclearPotential::NuclearPotentialIsospin::computePotentialEnergy(), G4INCL::NuclearPotential::NuclearPotentialConstant::computePotentialEnergy(), G4INCL::CrossSections::deltaProduction(), G4INCL::PauliStandard::getBlockingProbability(), G4INCL::NuclearPotential::INuclearPotential::getFermiEnergy(), G4INCL::NuclearPotential::INuclearPotential::getFermiMomentum(), G4INCL::RecombinationChannel::getFinalState(), G4INCL::ElasticChannel::getFinalState(), G4INCL::DeltaProductionChannel::getFinalState(), G4INCL::DeltaDecayChannel::getFinalState(), G4INCL::KinematicsUtils::getLocalEnergy(), G4INCL::NuclearPotential::INuclearPotential::getSeparationEnergy(), G4INCL::NuclearDensity::getTransmissionRadius(), G4INCL::Nucleus::insertParticle(), G4INCL::ParticleConfig::isPair(), G4INCL::CrossSections::pionNucleon(), G4INCL::InteractionAvatar::preInteractionBlocking(), and G4INCL::CrossSections::recombination().

00153                                        {
00154       return theType;
00155     };

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

Returns the charge number.

Definition at line 270 of file G4INCLParticle.hh.

References theZ.

Referenced by G4INCL::Cluster::addParticle(), G4INCL::CoulombNonRelativistic::bringToSurface(), G4INCL::ClusterDecay::decay(), G4INCL::CoulombNonRelativistic::distortOut(), G4INCL::Nucleus::fillEventInfo(), G4INCL::Nucleus::finalizeProjectileRemnant(), G4INCL::ClusteringModelIntercomparison::getCluster(), G4INCL::Nucleus::getConservationBalance(), G4INCL::ParticleEntryChannel::getFinalState(), G4INCL::Nucleus::getTransmissionBarrier(), G4INCL::SurfaceAvatar::getTransmissionProbability(), G4INCL::NuclearDensity::getTransmissionRadius(), G4INCL::ClusterUtils::getZ(), G4INCL::Nucleus::insertParticle(), G4INCL::Nucleus::isEventTransparent(), G4INCL::ClusterDecay::isStable(), G4INCL::ProjectileRemnant::removeParticle(), G4INCL::StandardPropagationModel::shootComposite(), and G4INCL::StandardPropagationModel::shootParticle().

00270 { return theZ; }

void G4INCL::Particle::incrementNumberOfCollisions (  )  [inline]

Increment the number of collisions undergone by the particle.

Definition at line 586 of file G4INCLParticle.hh.

References nCollisions.

00586 { nCollisions++; }

void G4INCL::Particle::incrementNumberOfDecays (  )  [inline]

Increment the number of decays undergone by the particle.

Definition at line 595 of file G4INCLParticle.hh.

References nDecays.

00595 { nDecays++; }

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

Definition at line 634 of file G4INCLParticle.hh.

References G4INCL::Composite, and theType.

Referenced by getEmissionQValueCorrection().

00634                              {
00635       return (theType == Composite);
00636     }

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

Is it a Delta?

Definition at line 261 of file G4INCLParticle.hh.

References G4INCL::DeltaMinus, G4INCL::DeltaPlus, G4INCL::DeltaPlusPlus, G4INCL::DeltaZero, and theType.

Referenced by G4INCL::DecayAvatar::getChannel(), G4INCL::BinaryCollisionAvatar::getChannel(), G4INCL::NuclearPotential::INuclearPotential::getFermiMomentum(), G4INCL::DeltaProductionChannel::getFinalState(), isResonance(), G4INCL::CrossSections::recombination(), G4INCL::RecombinationChannel::RecombinationChannel(), and G4INCL::CrossSections::total().

00261                                   {
00262       return (theType==DeltaPlusPlus || theType==DeltaPlus ||
00263           theType==DeltaZero || theType==DeltaMinus);
00264     }

G4bool G4INCL::Particle::isInList ( ParticleList const &  l  )  const [inline]

Check if the particle belongs to a given list.

Definition at line 630 of file G4INCLParticle.hh.

00630                                                  {
00631       return (std::find(l.begin(), l.end(), this)!=l.end());
00632     }

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

Is this a nucleon?

Definition at line 215 of file G4INCLParticle.hh.

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

Referenced by G4INCL::NuclearPotential::NuclearPotentialEnergyIsospinSmooth::computePotentialEnergy(), G4INCL::NuclearPotential::NuclearPotentialEnergyIsospin::computePotentialEnergy(), G4INCL::StandardPropagationModel::generateBinaryCollisionAvatar(), G4INCL::SurfaceAvatar::getChannel(), G4INCL::BinaryCollisionAvatar::getChannel(), G4INCL::PionNucleonChannel::getFinalState(), G4INCL::ParticleEntryChannel::getFinalState(), G4INCL::Nucleus::insertParticle(), G4INCL::CrossSections::pionNucleon(), G4INCL::StandardPropagationModel::shootParticle(), and G4INCL::CrossSections::total().

00215                              {
00216       if(theType == G4INCL::Proton || theType == G4INCL::Neutron)
00217         return true;
00218       else
00219         return false;
00220     };

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

Check if the particle is out of its potential well.

Definition at line 608 of file G4INCLParticle.hh.

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

00608 { return outOfWell; }

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

Definition at line 230 of file G4INCLParticle.hh.

References G4INCL::Participant, and theParticipantType.

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

00230                                  {
00231       return (theParticipantType==Participant);
00232     }

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

Is this a pion?

Definition at line 255 of file G4INCLParticle.hh.

References G4INCL::PiMinus, G4INCL::PiPlus, G4INCL::PiZero, and theType.

Referenced by G4INCL::CrossSections::elastic(), G4INCL::StandardPropagationModel::generateBinaryCollisionAvatar(), G4INCL::BinaryCollisionAvatar::getChannel(), G4INCL::Nucleus::getSurfaceRadius(), G4INCL::CrossSections::pionNucleon(), G4INCL::InteractionAvatar::preInteractionLocalEnergy(), and G4INCL::CrossSections::total().

00255 { return (theType == PiPlus || theType == PiZero || theType == PiMinus); }

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

Definition at line 238 of file G4INCLParticle.hh.

References G4INCL::ProjectileSpectator, and theParticipantType.

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

00238                                          {
00239       return (theParticipantType==ProjectileSpectator);
00240     }

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

Is it a resonance?

Definition at line 258 of file G4INCLParticle.hh.

References isDelta().

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

00258 { return isDelta(); }

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

Definition at line 234 of file G4INCLParticle.hh.

References G4INCL::TargetSpectator, and theParticipantType.

Referenced by G4INCL::SurfaceAvatar::getChannel(), G4INCL::Nucleus::insertParticle(), and G4INCL::SurfaceAvatar::postInteraction().

00234                                      {
00235       return (theParticipantType==TargetSpectator);
00236     }

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:
aBoostVector the boost vector (velocity) [c]
refPos the reference position

Definition at line 311 of file G4INCLParticle.hh.

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

00311                                                                                      {
00312       const G4double beta2 = aBoostVector.mag2();
00313       const G4double gamma = 1.0 / std::sqrt(1.0 - beta2);
00314       const ThreeVector theRelativePosition = thePosition - refPos;
00315       const ThreeVector transversePosition = theRelativePosition - aBoostVector * (theRelativePosition.dot(aBoostVector) / aBoostVector.mag2());
00316       const ThreeVector longitudinalPosition = theRelativePosition - transversePosition;
00317 
00318       thePosition = refPos + transversePosition + longitudinalPosition / gamma;
00319     }

virtual void G4INCL::Particle::makeParticipant (  )  [inline, virtual]

Reimplemented in G4INCL::Cluster.

Definition at line 242 of file G4INCLParticle.hh.

References G4INCL::Participant, and theParticipantType.

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

00242                                    {
00243       theParticipantType = Participant;
00244     }

virtual void G4INCL::Particle::makeProjectileSpectator (  )  [inline, virtual]

Reimplemented in G4INCL::Cluster.

Definition at line 250 of file G4INCLParticle.hh.

References G4INCL::ProjectileSpectator, and theParticipantType.

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

00250                                            {
00251       theParticipantType = ProjectileSpectator;
00252     }

virtual void G4INCL::Particle::makeTargetSpectator (  )  [inline, virtual]

Reimplemented in G4INCL::Cluster.

Definition at line 246 of file G4INCLParticle.hh.

References G4INCL::TargetSpectator, and theParticipantType.

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

00246                                        {
00247       theParticipantType = TargetSpectator;
00248     }

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

Assignment operator.

Does not copy the particle ID.

Definition at line 143 of file G4INCLParticle.hh.

References swap().

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

00143                                              {
00144       Particle temporaryParticle(rhs);
00145       swap(temporaryParticle);
00146       return *this;
00147     }

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

Reimplemented in G4INCL::Cluster.

Definition at line 686 of file G4INCLParticle.hh.

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

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

00686                             {
00687       std::stringstream ss;
00688       ss << "Particle (ID = " << ID << ") type = ";
00689       ss << ParticleTable::getName(theType);
00690       ss << std::endl
00691         << "   energy = " << theEnergy << std::endl
00692         << "   momentum = "
00693         << theMomentum.print()
00694         << std::endl
00695         << "   position = "
00696         << thePosition.print()
00697         << std::endl;
00698       return ss.str();
00699     };

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

Definition at line 575 of file G4INCLParticle.hh.

References thePosition.

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

00575                                   {
00576       thePosition += ((*thePropagationMomentum)*(step/(*thePropagationEnergy)));
00577     };

virtual void G4INCL::Particle::rotate ( const G4double  angle,
const ThreeVector axis 
) [inline, virtual]

Rotate the particle position and momentum.

Parameters:
angle the rotation angle
axis a unit vector representing the rotation axis

Reimplemented in G4INCL::Cluster.

Definition at line 680 of file G4INCLParticle.hh.

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

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

00680                                                                        {
00681       thePosition.rotate(angle, axis);
00682       theMomentum.rotate(angle, axis);
00683       theFrozenMomentum.rotate(angle, axis);
00684     }

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

Definition at line 610 of file G4INCLParticle.hh.

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

00610 { emissionTime = t; }

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

Set the energy of the particle in MeV.

Definition at line 532 of file G4INCLParticle.hh.

References theEnergy.

Referenced by G4INCL::PionNucleonChannel::getFinalState(), G4INCL::DeltaProductionChannel::getFinalState(), G4INCL::CrossSections::interactionDistanceNN(), G4INCL::CrossSections::interactionDistancePiN(), G4INCL::InteractionAvatar::restoreParticles(), G4INCL::StandardPropagationModel::shootParticle(), G4INCL::KinematicsUtils::transformToLocalEnergyFrame(), and G4INCL::Nucleus::useFusionKinematics().

00533     {
00534       this->theEnergy = energy;
00535     };

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

Set the frozen particle momentum.

Definition at line 642 of file G4INCLParticle.hh.

References theFrozenEnergy.

00642 { theFrozenEnergy = energy; }

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

Set the frozen particle momentum.

Definition at line 639 of file G4INCLParticle.hh.

References theFrozenMomentum.

00639 { theFrozenMomentum = momentum; }

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

Definition at line 573 of file G4INCLParticle.hh.

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

00573 { theHelicity = h; };

void G4INCL::Particle::setINCLMass (  )  [inline]

Set the mass of the Particle to its table mass.

Definition at line 418 of file G4INCLParticle.hh.

References getINCLMass(), and setMass().

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

00418 { setMass(getINCLMass()); }

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

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

Definition at line 524 of file G4INCLParticle.hh.

Referenced by G4INCL::Nucleus::computeRecoilKinematics(), G4INCL::Nucleus::finalizeProjectileRemnant(), G4INCL::PionNucleonChannel::getFinalState(), G4INCL::DeltaProductionChannel::getFinalState(), Particle(), G4INCL::InteractionAvatar::restoreParticles(), setINCLMass(), setRealMass(), setTableMass(), and G4INCL::Nucleus::useFusionKinematics().

00525     {
00526       this->theMass = mass;
00527     }

virtual void G4INCL::Particle::setMomentum ( const G4INCL::ThreeVector momentum  )  [inline, virtual]

Set the momentum vector.

Definition at line 554 of file G4INCLParticle.hh.

References theMomentum.

Referenced by G4INCL::Nucleus::decayOutgoingDeltas(), G4INCL::ReflectionChannel::getFinalState(), G4INCL::RecombinationChannel::getFinalState(), G4INCL::PionNucleonChannel::getFinalState(), G4INCL::ElasticChannel::getFinalState(), G4INCL::DeltaProductionChannel::getFinalState(), G4INCL::DeltaDecayChannel::getFinalState(), G4INCL::InteractionAvatar::restoreParticles(), and G4INCL::Nucleus::useFusionKinematics().

00555     {
00556       this->theMomentum = momentum;
00557     };

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

Set the number of collisions undergone by the particle.

Definition at line 583 of file G4INCLParticle.hh.

References nCollisions.

00583 { nCollisions = n; }

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

Set the number of decays undergone by the particle.

Definition at line 592 of file G4INCLParticle.hh.

References nDecays.

00592 { nDecays = 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 605 of file G4INCLParticle.hh.

00605 { outOfWell = true; }

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

Definition at line 226 of file G4INCLParticle.hh.

References theParticipantType.

00226                                                      {
00227       theParticipantType = p;
00228     }

virtual void G4INCL::Particle::setPosition ( const G4INCL::ThreeVector position  )  [inline, virtual]

Reimplemented in G4INCL::Cluster.

Definition at line 567 of file G4INCLParticle.hh.

References position, and thePosition.

Referenced by G4INCL::InteractionAvatar::bringParticleInside(), G4INCL::CoulombNone::bringToSurface(), G4INCL::ClusteringModelIntercomparison::getCluster(), G4INCL::ReflectionChannel::getFinalState(), G4INCL::InteractionAvatar::restoreParticles(), G4INCL::ParticleSampler::sampleParticles(), G4INCL::Cluster::setPosition(), and G4INCL::StandardPropagationModel::shootParticle().

00568     {
00569       this->thePosition = position;
00570     };

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

Set the particle potential energy.

Definition at line 511 of file G4INCLParticle.hh.

References thePotentialEnergy.

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

00511 { thePotentialEnergy = v; }

void G4INCL::Particle::setRealMass (  )  [inline]

Set the mass of the Particle to its real mass.

Definition at line 412 of file G4INCLParticle.hh.

References getRealMass(), and setMass().

00412 { setMass(getRealMass()); }

void G4INCL::Particle::setTableMass (  )  [inline]

Set the mass of the Particle to its table mass.

Definition at line 415 of file G4INCLParticle.hh.

References getTableMass(), and setMass().

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

00415 { setMass(getTableMass()); }

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

Definition at line 162 of file G4INCLParticle.hh.

References G4INCL::Composite, G4INCL::DeltaMinus, G4INCL::DeltaPlus, G4INCL::DeltaPlusPlus, G4INCL::DeltaZero, 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::RecombinationChannel::getFinalState(), G4INCL::PionNucleonChannel::getFinalState(), G4INCL::ElasticChannel::getFinalState(), G4INCL::DeltaProductionChannel::getFinalState(), G4INCL::DeltaDecayChannel::getFinalState(), Particle(), and G4INCL::InteractionAvatar::restoreParticles().

00162                                  {
00163       theType = t;
00164       switch(theType)
00165       {
00166         case DeltaPlusPlus:
00167           theA = 1;
00168           theZ = 2;
00169           break;
00170         case Proton:
00171         case DeltaPlus:
00172           theA = 1;
00173           theZ = 1;
00174           break;
00175         case Neutron:
00176         case DeltaZero:
00177           theA = 1;
00178           theZ = 0;
00179           break;
00180         case DeltaMinus:
00181           theA = 1;
00182           theZ = -1;
00183           break;
00184         case PiPlus:
00185           theA = 0;
00186           theZ = 1;
00187           break;
00188         case PiZero:
00189           theA = 0;
00190           theZ = 0;
00191           break;
00192         case PiMinus:
00193           theA = 0;
00194           theZ = -1;
00195           break;
00196         case Composite:
00197          // ERROR("Trying to set particle type to Composite! Construct a Cluster object instead" << std::endl);
00198           theA = 0;
00199           theZ = 0;
00200           break;
00201         case UnknownParticle:
00202           theA = 0;
00203           theZ = 0;
00204           ERROR("Trying to set particle type to Unknown!" << std::endl);
00205           break;
00206       }
00207 
00208       if( !isResonance() && t!=Composite )
00209         setINCLMass();
00210     }

void G4INCL::Particle::swap ( Particle rhs  )  [inline, protected]

Helper method for the assignment operator.

Definition at line 107 of file G4INCLParticle.hh.

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

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

00107                              {
00108       std::swap(theZ, rhs.theZ);
00109       std::swap(theA, rhs.theA);
00110       std::swap(theParticipantType, rhs.theParticipantType);
00111       std::swap(theType, rhs.theType);
00112       if(rhs.thePropagationEnergy == &(rhs.theFrozenEnergy))
00113         thePropagationEnergy = &theFrozenEnergy;
00114       else
00115         thePropagationEnergy = &theEnergy;
00116       std::swap(theEnergy, rhs.theEnergy);
00117       std::swap(theFrozenEnergy, rhs.theFrozenEnergy);
00118       if(rhs.thePropagationMomentum == &(rhs.theFrozenMomentum))
00119         thePropagationMomentum = &theFrozenMomentum;
00120       else
00121         thePropagationMomentum = &theMomentum;
00122       std::swap(theMomentum, rhs.theMomentum);
00123       std::swap(theFrozenMomentum, rhs.theFrozenMomentum);
00124       std::swap(thePosition, rhs.thePosition);
00125       std::swap(nCollisions, rhs.nCollisions);
00126       std::swap(nDecays, rhs.nDecays);
00127       std::swap(thePotentialEnergy, rhs.thePotentialEnergy);
00128       // ID intentionally not swapped
00129 
00130       std::swap(theHelicity, rhs.theHelicity);
00131       std::swap(emissionTime, rhs.emissionTime);
00132       std::swap(outOfWell, rhs.outOfWell);
00133 
00134       std::swap(theMass, rhs.theMass);
00135     }

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

References theEnergy, theMomentum, thePropagationEnergy, and thePropagationMomentum.

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


Field Documentation

long G4INCL::Particle::ID [protected]

Definition at line 738 of file G4INCLParticle.hh.

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

G4int G4INCL::Particle::nCollisions [protected]

Definition at line 735 of file G4INCLParticle.hh.

Referenced by G4INCL::Cluster::addParticle(), getNumberOfCollisions(), incrementNumberOfCollisions(), G4INCL::ProjectileRemnant::reset(), setNumberOfCollisions(), and swap().

G4int G4INCL::Particle::nDecays [protected]

Definition at line 736 of file G4INCLParticle.hh.

Referenced by getNumberOfDecays(), incrementNumberOfDecays(), setNumberOfDecays(), and swap().

G4int G4INCL::Particle::theA [protected]

Definition at line 725 of file G4INCLParticle.hh.

Referenced by G4INCL::ProjectileRemnant::addMostDynamicalSpectators(), G4INCL::Cluster::addParticle(), G4INCL::Nucleus::applyFinalState(), G4INCL::Cluster::Cluster(), G4INCL::ProjectileRemnant::computeExcitationEnergy(), G4INCL::Nucleus::computeRecoilKinematics(), G4INCL::Nucleus::decayInsideDeltas(), G4INCL::Nucleus::decayMe(), G4INCL::Nucleus::emitInsidePions(), getA(), G4INCL::Nucleus::getCoulombRadius(), getEmissionQValueCorrection(), getINCLMass(), getRealMass(), G4INCL::Cluster::getSpecies(), getTableMass(), getTransferQValueCorrection(), G4INCL::Cluster::initializeParticles(), G4INCL::Nucleus::insertParticle(), G4INCL::Cluster::internalBoostToCM(), G4INCL::Nucleus::Nucleus(), G4INCL::Cluster::print(), G4INCL::ProjectileRemnant::removeParticle(), G4INCL::ProjectileRemnant::reset(), G4INCL::Cluster::setA(), setType(), and swap().

G4double G4INCL::Particle::theEnergy [protected]

Definition at line 728 of file G4INCLParticle.hh.

Referenced by G4INCL::ProjectileRemnant::addMostDynamicalSpectators(), G4INCL::Cluster::addParticle(), adjustEnergyFromMomentum(), adjustMomentumFromEnergy(), boost(), boostVector(), dump(), getBeta(), getEnergy(), getInvariantMass(), getKineticEnergy(), G4INCL::Cluster::internalBoostToCM(), Particle(), print(), G4INCL::Cluster::print(), G4INCL::ProjectileRemnant::removeParticle(), G4INCL::ProjectileRemnant::reset(), setEnergy(), swap(), thawPropagation(), and G4INCL::Nucleus::useFusionKinematics().

G4double G4INCL::Particle::theFrozenEnergy [protected]

Definition at line 730 of file G4INCLParticle.hh.

Referenced by freezePropagation(), getFrozenEnergy(), Particle(), setFrozenEnergy(), and swap().

G4INCL::ThreeVector G4INCL::Particle::theFrozenMomentum [protected]

Definition at line 733 of file G4INCLParticle.hh.

Referenced by freezePropagation(), getFrozenMomentum(), Particle(), rotate(), setFrozenMomentum(), and swap().

G4INCL::ThreeVector G4INCL::Particle::theMomentum [protected]

Definition at line 731 of file G4INCLParticle.hh.

Referenced by G4INCL::ProjectileRemnant::addMostDynamicalSpectators(), G4INCL::Cluster::addParticle(), adjustEnergyFromMomentum(), adjustMomentumFromEnergy(), boost(), boostVector(), G4INCL::Nucleus::computeRecoilKinematics(), dump(), G4INCL::Cluster::freezeInternalMotion(), getAngularMomentum(), getBeta(), getInvariantMass(), getMomentum(), G4INCL::Cluster::internalBoostToCM(), Particle(), print(), G4INCL::Cluster::print(), G4INCL::ProjectileRemnant::removeParticle(), G4INCL::ProjectileRemnant::reset(), rotate(), setMomentum(), swap(), thawPropagation(), and G4INCL::Nucleus::useFusionKinematics().

ParticipantType G4INCL::Particle::theParticipantType [protected]

Definition at line 726 of file G4INCLParticle.hh.

Referenced by getParticipantType(), isParticipant(), isProjectileSpectator(), isTargetSpectator(), makeParticipant(), makeProjectileSpectator(), makeTargetSpectator(), Particle(), setParticipantType(), and swap().

G4INCL::ThreeVector G4INCL::Particle::thePosition [protected]

Definition at line 734 of file G4INCLParticle.hh.

Referenced by G4INCL::Cluster::addParticle(), G4INCL::Cluster::boost(), G4INCL::Cluster::Cluster(), G4INCL::Nucleus::computeRecoilKinematics(), dump(), getAngularMomentum(), getLongitudinalPosition(), getPosition(), getTransversePosition(), G4INCL::Nucleus::initializeParticles(), G4INCL::Cluster::initializeParticles(), G4INCL::Cluster::internalBoostToCM(), lorentzContract(), print(), G4INCL::Cluster::print(), propagate(), G4INCL::ProjectileRemnant::reset(), rotate(), setPosition(), G4INCL::Cluster::setPosition(), and swap().

G4double G4INCL::Particle::thePotentialEnergy [protected]

Definition at line 737 of file G4INCLParticle.hh.

Referenced by G4INCL::Cluster::addParticle(), getPotentialEnergy(), G4INCL::ProjectileRemnant::reset(), setPotentialEnergy(), and swap().

G4double* G4INCL::Particle::thePropagationEnergy [protected]

Definition at line 729 of file G4INCLParticle.hh.

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

G4INCL::ThreeVector* G4INCL::Particle::thePropagationMomentum [protected]

Definition at line 732 of file G4INCLParticle.hh.

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

G4INCL::ParticleType G4INCL::Particle::theType [protected]

Definition at line 727 of file G4INCLParticle.hh.

Referenced by dump(), getINCLMass(), getRealMass(), getSpecies(), getTableMass(), getType(), isCluster(), isDelta(), isNucleon(), isPion(), print(), G4INCL::Cluster::print(), setType(), and swap().

G4int G4INCL::Particle::theZ [protected]

Definition at line 725 of file G4INCLParticle.hh.

Referenced by G4INCL::ProjectileRemnant::addMostDynamicalSpectators(), G4INCL::Cluster::addParticle(), G4INCL::Nucleus::applyFinalState(), G4INCL::Cluster::Cluster(), G4INCL::Nucleus::computeRecoilKinematics(), G4INCL::Nucleus::decayInsideDeltas(), G4INCL::Nucleus::decayMe(), G4INCL::Nucleus::emitInsidePions(), G4INCL::Nucleus::getCoulombRadius(), getEmissionQValueCorrection(), getINCLMass(), getRealMass(), G4INCL::Cluster::getSpecies(), getTableMass(), getTransferQValueCorrection(), G4INCL::Nucleus::getTransmissionBarrier(), getZ(), G4INCL::Cluster::initializeParticles(), G4INCL::Nucleus::insertParticle(), G4INCL::Nucleus::Nucleus(), G4INCL::Cluster::print(), G4INCL::ProjectileRemnant::removeParticle(), G4INCL::ProjectileRemnant::reset(), setType(), G4INCL::Cluster::setZ(), and swap().


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:54:08 2013 for Geant4 by  doxygen 1.4.7