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

#include <G4INCLCluster.hh>

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

Public Member Functions

 Cluster (const G4int Z, const G4int A, const G4bool createParticleSampler=true)
 Standard Cluster constructor. More...
 
template<class Iterator >
 Cluster (Iterator begin, Iterator end)
 
virtual ~Cluster ()
 
 Cluster (const Cluster &rhs)
 Copy constructor. More...
 
Clusteroperator= (const Cluster &rhs)
 Assignment operator. More...
 
void swap (Cluster &rhs)
 Helper method for the assignment operator. More...
 
ParticleSpecies getSpecies () const
 Get the particle species. More...
 
void deleteParticles ()
 
void clearParticles ()
 
void setZ (const G4int Z)
 Set the charge number of the cluster. More...
 
void setA (const G4int A)
 Set the mass number of the cluster. More...
 
G4double getExcitationEnergy () const
 Get the excitation energy of the cluster. More...
 
void setExcitationEnergy (const G4double e)
 Set the excitation energy of the cluster. More...
 
virtual G4double getTableMass () const
 Get the real particle mass. More...
 
ParticleList const & getParticles () const
 
void removeParticle (Particle *const p)
 Remove a particle from the cluster components. More...
 
void addParticle (Particle *const p)
 
void addParticles (ParticleList const &pL)
 Add a list of particles to the cluster. More...
 
ParticleList getParticleList () const
 Returns the list of particles that make up the cluster. More...
 
std::string print () const
 
virtual void initializeParticles ()
 Initialise the NuclearDensity pointer and sample the particles. More...
 
void internalBoostToCM ()
 Boost to the CM of the component particles. More...
 
void putParticlesOffShell ()
 Put the cluster components off shell. More...
 
void setPosition (const ThreeVector &position)
 Set the position of the cluster. More...
 
void boost (const ThreeVector &aBoostVector)
 Boost the cluster with the indicated velocity. More...
 
void freezeInternalMotion ()
 Freeze the internal motion of the particles. More...
 
virtual void rotate (const G4double angle, const ThreeVector &axis)
 Rotate position and momentum of all the particles. More...
 
virtual void makeProjectileSpectator ()
 Make all the components projectile spectators, too. More...
 
virtual void makeTargetSpectator ()
 Make all the components target spectators, too. More...
 
virtual void makeParticipant ()
 Make all the components participants, too. More...
 
ThreeVector const & getSpin () const
 Get the spin of the nucleus. More...
 
void setSpin (const ThreeVector &j)
 Set the spin of the nucleus. More...
 
G4INCL::ThreeVector getAngularMomentum () const
 Get the total angular momentum (orbital + spin) More...
 
- Public Member Functions inherited from G4INCL::Particle
 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
 
void setType (ParticleType t)
 
G4bool isNucleon () const
 
ParticipantType getParticipantType () const
 
void setParticipantType (ParticipantType const p)
 
G4bool isParticipant () const
 
G4bool isTargetSpectator () const
 
G4bool isProjectileSpectator () const
 
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...
 
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 void setMomentum (const G4INCL::ThreeVector &momentum)
 
const G4INCL::ThreeVectorgetPosition () const
 
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...
 
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 Attributes

ParticleList particles
 
G4double theExcitationEnergy
 
ThreeVector theSpin
 
ParticleSamplertheParticleSampler
 
- Protected Attributes inherited from G4INCL::Particle
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
 

Additional Inherited Members

- Protected Member Functions inherited from G4INCL::Particle
void swap (Particle &rhs)
 Helper method for the assignment operator. More...
 

Detailed Description

Cluster is a particle (inherits from the Particle class) that is actually a collection of elementary particles.

Definition at line 50 of file G4INCLCluster.hh.

Constructor & Destructor Documentation

G4INCL::Cluster::Cluster ( const G4int  Z,
const G4int  A,
const G4bool  createParticleSampler = true 
)
inline

Standard Cluster constructor.

This constructor should mainly be used when constructing Nucleus or when constructing Clusters to be used as composite projectiles.

Definition at line 58 of file G4INCLCluster.hh.

References G4INCL::Composite, G4INCL::Particle::setINCLMass(), G4INCL::Particle::setType(), G4INCL::Particle::theA, theParticleSampler, and G4INCL::Particle::theZ.

58  :
59  Particle(),
61  theSpin(0.,0.,0.),
62  theParticleSampler(NULL)
63  {
65  theZ = Z;
66  theA = A;
67  setINCLMass();
68  if(createParticleSampler)
69  theParticleSampler = new ParticleSampler(A,Z);
70  }
void setINCLMass()
Set the mass of the Particle to its table mass.
ThreeVector theSpin
ParticleSampler * theParticleSampler
G4double theExcitationEnergy
void setType(ParticleType t)
template<class Iterator >
G4INCL::Cluster::Cluster ( Iterator  begin,
Iterator  end 
)
inline

A cluster can be directly built from a list of particles.

Definition at line 76 of file G4INCLCluster.hh.

References addParticle(), G4INCL::Particle::adjustMomentumFromEnergy(), G4INCL::Composite, G4INCL::Particle::setINCLMass(), G4INCL::Particle::setType(), G4INCL::Particle::theA, and G4INCL::Particle::thePosition.

76  :
77  Particle(),
79  theSpin(0.,0.,0.),
80  theParticleSampler(NULL)
81  {
83  for(Iterator i = begin; i != end; ++i) {
84  addParticle(*i);
85  }
86  thePosition /= theA;
87  setINCLMass();
89  }
void setINCLMass()
Set the mass of the Particle to its table mass.
ThreeVector theSpin
G4INCL::ThreeVector thePosition
ParticleSampler * theParticleSampler
G4double theExcitationEnergy
void setType(ParticleType t)
void addParticle(Particle *const p)
const ThreeVector & adjustMomentumFromEnergy()
Rescale the momentum to match the total energy.
virtual G4INCL::Cluster::~Cluster ( )
inlinevirtual

Definition at line 91 of file G4INCLCluster.hh.

References theParticleSampler.

91  {
92  delete theParticleSampler;
93  }
ParticleSampler * theParticleSampler
G4INCL::Cluster::Cluster ( const Cluster rhs)
inline

Copy constructor.

Definition at line 96 of file G4INCLCluster.hh.

References G4INCL::Particle::Particle(), particles, G4INCL::Particle::theA, theParticleSampler, and G4INCL::Particle::theZ.

96  :
97  Particle(rhs),
98  theExcitationEnergy(rhs.theExcitationEnergy),
99  theSpin(rhs.theSpin)
100  {
101  for(ParticleIter p=rhs.particles.begin(), e=rhs.particles.end(); p!=e; ++p) {
102  particles.push_back(new Particle(**p));
103  }
104  if(rhs.theParticleSampler)
105  theParticleSampler = new ParticleSampler(rhs.theA,rhs.theZ);
106  else
107  theParticleSampler = NULL;
108  }
const char * p
Definition: xmltok.h:285
ThreeVector theSpin
ParticleList particles
ParticleSampler * theParticleSampler
G4double theExcitationEnergy
ParticleList::const_iterator ParticleIter

Member Function Documentation

void G4INCL::Cluster::addParticle ( Particle *const  p)
inline

Add one particle to the cluster. This updates the cluster mass, energy, size, etc.

Definition at line 172 of file G4INCLCluster.hh.

References G4INCL::Particle::getA(), G4INCL::Particle::getEnergy(), G4INCL::Particle::getMomentum(), G4INCL::Particle::getNumberOfCollisions(), G4INCL::Particle::getPosition(), G4INCL::Particle::getPotentialEnergy(), G4INCL::Particle::getZ(), G4INCL::Particle::nCollisions, particles, G4INCL::Particle::theA, G4INCL::Particle::theEnergy, G4INCL::Particle::theMomentum, G4INCL::Particle::thePosition, G4INCL::Particle::thePotentialEnergy, and G4INCL::Particle::theZ.

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

172  {
173 // assert(p->isNucleon());
174  particles.push_back(p);
175  theEnergy += p->getEnergy();
176  thePotentialEnergy += p->getPotentialEnergy();
177  theMomentum += p->getMomentum();
178  thePosition += p->getPosition();
179  theA += p->getA();
180  theZ += p->getZ();
181  nCollisions += p->getNumberOfCollisions();
182  }
const char * p
Definition: xmltok.h:285
ParticleList particles
G4INCL::ThreeVector thePosition
G4INCL::ThreeVector theMomentum
G4double thePotentialEnergy
void G4INCL::Cluster::addParticles ( ParticleList const &  pL)
inline

Add a list of particles to the cluster.

Definition at line 185 of file G4INCLCluster.hh.

References addParticle().

Referenced by initializeParticles().

185  {
186  for(ParticleIter p=pL.begin(), e=pL.end(); p!=e; ++p)
187  addParticle(*p);
188  }
const char * p
Definition: xmltok.h:285
void addParticle(Particle *const p)
ParticleList::const_iterator ParticleIter
void G4INCL::Cluster::boost ( const ThreeVector aBoostVector)
inline

Boost the cluster with the indicated velocity.

The Cluster is boosted as a whole, just like any Particle object; moreover, the internal components (particles list) are also boosted, according to Alain Boudard's off-shell recipe.

Parameters
aBoostVectorthe velocity to boost to [c]

Definition at line 315 of file G4INCLCluster.hh.

References G4INCL::Particle::boost(), G4INCL::ThreeVector::getX(), G4INCL::ThreeVector::getY(), G4INCL::ThreeVector::getZ(), INCL_DEBUG, particles, print(), and G4INCL::Particle::thePosition.

Referenced by G4INCL::ProjectileRemnant::ProjectileRemnant().

315  {
316  Particle::boost(aBoostVector);
317  for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
318  (*p)->boost(aBoostVector);
319  // Apply Lorentz contraction to the particle position
320  (*p)->lorentzContract(aBoostVector,thePosition);
321  (*p)->rpCorrelate();
322  }
323 
324  INCL_DEBUG("Cluster was boosted with (bx,by,bz)=("
325  << aBoostVector.getX() << ", " << aBoostVector.getY() << ", " << aBoostVector.getZ() << "):"
326  << std::endl << print());
327 
328  }
const char * p
Definition: xmltok.h:285
void boost(const ThreeVector &aBoostVector)
std::string print() const
ParticleList particles
G4INCL::ThreeVector thePosition
#define INCL_DEBUG(x)
ParticleList::const_iterator ParticleIter
void G4INCL::Cluster::clearParticles ( )
inline

Definition at line 140 of file G4INCLCluster.hh.

References particles.

Referenced by deleteParticles().

140 { particles.clear(); }
ParticleList particles
void G4INCL::Cluster::deleteParticles ( )
inline

Definition at line 133 of file G4INCLCluster.hh.

References clearParticles(), and particles.

Referenced by G4INCL::Store::clearOutgoing(), G4INCL::Nucleus::decayOutgoingClusters(), G4INCL::ProjectileRemnant::reset(), and G4INCL::ProjectileRemnant::~ProjectileRemnant().

133  {
134  for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
135  delete (*p);
136  }
137  clearParticles();
138  }
const char * p
Definition: xmltok.h:285
ParticleList particles
ParticleList::const_iterator ParticleIter
void G4INCL::Cluster::freezeInternalMotion ( )
inline

Freeze the internal motion of the particles.

Each particle is assigned a frozen momentum four-vector determined by the collective cluster velocity. This is used for propagation, but not for dynamics. Normal propagation is restored by calling the Particle::thawPropagation() method, which should be done in InteractionAvatar::postInteraction.

Definition at line 338 of file G4INCLCluster.hh.

References G4INCL::Particle::getMass(), G4INCL::ThreeVector::mag2(), particles, and G4INCL::Particle::theMomentum.

Referenced by G4INCL::ProjectileRemnant::ProjectileRemnant().

338  {
339  const ThreeVector &normMomentum = theMomentum / getMass();
340  for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
341  const G4double pMass = (*p)->getMass();
342  const ThreeVector frozenMomentum = normMomentum * pMass;
343  const G4double frozenEnergy = std::sqrt(frozenMomentum.mag2()+pMass*pMass);
344  (*p)->setFrozenMomentum(frozenMomentum);
345  (*p)->setFrozenEnergy(frozenEnergy);
346  (*p)->freezePropagation();
347  }
348  }
G4double getMass() const
Get the cached particle mass.
const char * p
Definition: xmltok.h:285
ParticleList particles
G4INCL::ThreeVector theMomentum
double G4double
Definition: G4Types.hh:76
ParticleList::const_iterator ParticleIter
G4INCL::ThreeVector G4INCL::Cluster::getAngularMomentum ( ) const
inlinevirtual

Get the total angular momentum (orbital + spin)

Reimplemented from G4INCL::Particle.

Definition at line 395 of file G4INCLCluster.hh.

References G4INCL::Particle::getAngularMomentum(), and getSpin().

Referenced by G4INCL::Nucleus::computeRecoilKinematics(), and G4INCL::StandardPropagationModel::shootComposite().

395  {
397  }
virtual G4INCL::ThreeVector getAngularMomentum() const
ThreeVector const & getSpin() const
Get the spin of the nucleus.
G4double G4INCL::Cluster::getExcitationEnergy ( ) const
inline

Get the excitation energy of the cluster.

Definition at line 149 of file G4INCLCluster.hh.

References theExcitationEnergy.

Referenced by G4INCL::Nucleus::fillEventInfo(), and G4INCL::Nucleus::getConservationBalance().

149 { return theExcitationEnergy; }
G4double theExcitationEnergy
ParticleList G4INCL::Cluster::getParticleList ( ) const
inline

Returns the list of particles that make up the cluster.

Definition at line 191 of file G4INCLCluster.hh.

References particles.

191 { return particles; }
ParticleList particles
ParticleList const& G4INCL::Cluster::getParticles ( ) const
inline

Get the list of particles in the cluster.

Definition at line 163 of file G4INCLCluster.hh.

References particles.

Referenced by G4INCL::Nucleus::applyFinalState(), G4INCL::CoulombNone::bringToSurface(), and G4INCL::SurfaceAvatar::postInteraction().

163 { return particles; }
ParticleList particles
ParticleSpecies G4INCL::Cluster::getSpecies ( ) const
inlinevirtual

Get the particle species.

Reimplemented from G4INCL::Particle.

Definition at line 129 of file G4INCLCluster.hh.

References G4INCL::Particle::theA, and G4INCL::Particle::theZ.

129  {
130  return ParticleSpecies(theA, theZ);
131  }
ThreeVector const& G4INCL::Cluster::getSpin ( ) const
inline

Get the spin of the nucleus.

Definition at line 389 of file G4INCLCluster.hh.

References theSpin.

Referenced by G4INCL::Nucleus::fillEventInfo(), and getAngularMomentum().

389 { return theSpin; }
ThreeVector theSpin
virtual G4double G4INCL::Cluster::getTableMass ( ) const
inlinevirtual

Get the real particle mass.

Overloads the Particle method.

Reimplemented from G4INCL::Particle.

Definition at line 158 of file G4INCLCluster.hh.

References G4INCL::Particle::getRealMass().

Referenced by G4INCL::Nucleus::useFusionKinematics().

158 { return getRealMass(); }
G4double getRealMass() const
Get the real particle mass.
void G4INCL::Cluster::initializeParticles ( )
virtual

Initialise the NuclearDensity pointer and sample the particles.

Reimplemented in G4INCL::Nucleus.

Definition at line 42 of file G4INCLCluster.cc.

References addParticles(), INCL_DEBUG, print(), G4INCL::ParticleSampler::sampleParticles(), G4INCL::Particle::theA, theParticleSampler, G4INCL::Particle::thePosition, and G4INCL::Particle::theZ.

Referenced by G4INCL::Nucleus::initializeParticles(), and G4INCL::ProjectileRemnant::ProjectileRemnant().

42  {
43 // assert(theA>=2);
44  const ThreeVector oldPosition = thePosition;
46 #if !defined(NDEBUG) && !defined(INCLXX_IN_GEANT4_MODE)
47  const G4int theMassNumber = theA;
48  const G4int theChargeNumber = theZ;
49 #endif
50  theA = 0;
51  theZ = 0;
52  addParticles(theParticles); // add the particles to the `particles' list
53  thePosition = oldPosition;
54 // assert(theMassNumber==theA && theChargeNumber==theZ);
55  INCL_DEBUG("Cluster initialized:" << std::endl << print());
56  }
ParticleList sampleParticles(ThreeVector const &position)
int G4int
Definition: G4Types.hh:78
UnorderedVector< Particle * > ParticleList
std::string print() const
G4INCL::ThreeVector thePosition
ParticleSampler * theParticleSampler
void addParticles(ParticleList const &pL)
Add a list of particles to the cluster.
#define INCL_DEBUG(x)
void G4INCL::Cluster::internalBoostToCM ( )
inline

Boost to the CM of the component particles.

The position of all particles in the particles list is shifted so that their centre of mass is in the origin and their total momentum is zero.

Definition at line 225 of file G4INCLCluster.hh.

References G4INCL::Particle::getMass(), INCL_DEBUG, particles, print(), G4INCL::ThreeVector::setX(), G4INCL::ThreeVector::setY(), G4INCL::ThreeVector::setZ(), G4INCL::Particle::theA, G4INCL::Particle::theEnergy, G4INCL::Particle::theMomentum, and G4INCL::Particle::thePosition.

Referenced by G4INCL::ProjectileRemnant::ProjectileRemnant().

225  {
226 
227  // First compute the current CM position and total momentum
228  ThreeVector theCMPosition, theTotalMomentum;
229  G4double theTotalEnergy = 0.0;
230  for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
231  theCMPosition += (*p)->getPosition();
232  theTotalMomentum += (*p)->getMomentum();
233  theTotalEnergy += (*p)->getEnergy();
234  }
235  theCMPosition /= theA;
236 // assert((unsigned int)theA==particles.size());
237 
238  // Now determine the CM velocity of the particles
239  // commented out because currently unused, see below
240  // ThreeVector betaCM = theTotalMomentum / theTotalEnergy;
241 
242  // The new particle positions and momenta are scaled by a factor of
243  // \f$\sqrt{A/(A-1)}\f$, so that the resulting density distributions in
244  // the CM have the same variance as the one we started with.
245  const G4double rescaling = std::sqrt(((G4double)theA)/((G4double)(theA-1)));
246 
247  // Loop again to boost and reposition
248  for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
249  // \bug{We should do the following, but the Fortran version actually
250  // does not!
251  // (*p)->boost(betaCM);
252  // Here is what the Fortran version does:}
253  (*p)->setMomentum(((*p)->getMomentum()-theTotalMomentum/theA)*rescaling);
254 
255  // Set the CM position of the particles
256  (*p)->setPosition(((*p)->getPosition()-theCMPosition)*rescaling);
257  }
258 
259  // Set the global cluster kinematic variables
260  thePosition.setX(0.0);
261  thePosition.setY(0.0);
262  thePosition.setZ(0.0);
263  theMomentum.setX(0.0);
264  theMomentum.setY(0.0);
265  theMomentum.setZ(0.0);
266  theEnergy = getMass();
267 
268  INCL_DEBUG("Cluster boosted to internal CM:" << std::endl << print());
269 
270  }
G4double getMass() const
Get the cached particle mass.
const char * p
Definition: xmltok.h:285
void setY(G4double ay)
Set the y coordinate.
std::string print() const
ParticleList particles
G4INCL::ThreeVector thePosition
G4INCL::ThreeVector theMomentum
void setX(G4double ax)
Set the x coordinate.
double G4double
Definition: G4Types.hh:76
void setZ(G4double az)
Set the z coordinate.
#define INCL_DEBUG(x)
ParticleList::const_iterator ParticleIter
virtual void G4INCL::Cluster::makeParticipant ( )
inlinevirtual

Make all the components participants, too.

Reimplemented from G4INCL::Particle.

Definition at line 381 of file G4INCLCluster.hh.

References G4INCL::Particle::makeParticipant(), and particles.

381  {
383  for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
384  (*p)->makeParticipant();
385  }
386  }
const char * p
Definition: xmltok.h:285
ParticleList particles
virtual void makeParticipant()
ParticleList::const_iterator ParticleIter
virtual void G4INCL::Cluster::makeProjectileSpectator ( )
inlinevirtual

Make all the components projectile spectators, too.

Reimplemented from G4INCL::Particle.

Definition at line 365 of file G4INCLCluster.hh.

References G4INCL::Particle::makeProjectileSpectator(), and particles.

Referenced by G4INCL::ProjectileRemnant::ProjectileRemnant().

365  {
367  for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
368  (*p)->makeProjectileSpectator();
369  }
370  }
const char * p
Definition: xmltok.h:285
virtual void makeProjectileSpectator()
ParticleList particles
ParticleList::const_iterator ParticleIter
virtual void G4INCL::Cluster::makeTargetSpectator ( )
inlinevirtual

Make all the components target spectators, too.

Reimplemented from G4INCL::Particle.

Definition at line 373 of file G4INCLCluster.hh.

References G4INCL::Particle::makeTargetSpectator(), and particles.

373  {
375  for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
376  (*p)->makeTargetSpectator();
377  }
378  }
virtual void makeTargetSpectator()
const char * p
Definition: xmltok.h:285
ParticleList particles
ParticleList::const_iterator ParticleIter
Cluster& G4INCL::Cluster::operator= ( const Cluster rhs)
inline

Assignment operator.

Definition at line 111 of file G4INCLCluster.hh.

References G4INCL::Particle::operator=(), and swap().

111  {
112  Cluster temporaryCluster(rhs);
113  Particle::operator=(temporaryCluster);
114  swap(temporaryCluster);
115  return *this;
116  }
void swap(Cluster &rhs)
Helper method for the assignment operator.
Particle & operator=(const Particle &rhs)
Assignment operator.
Cluster(const G4int Z, const G4int A, const G4bool createParticleSampler=true)
Standard Cluster constructor.
std::string G4INCL::Cluster::print ( ) const
inline

Definition at line 193 of file G4INCLCluster.hh.

References G4INCL::Particle::getMass(), G4INCL::ParticleTable::getName(), G4INCL::Particle::ID, particles, G4INCL::ThreeVector::print(), G4INCL::Particle::theA, G4INCL::Particle::theEnergy, G4INCL::Particle::theMomentum, G4INCL::Particle::thePosition, G4INCL::Particle::theType, and G4INCL::Particle::theZ.

Referenced by boost(), G4INCL::SurfaceAvatar::getChannel(), initializeParticles(), internalBoostToCM(), putParticlesOffShell(), G4INCL::ProjectileRemnant::removeParticle(), and G4INCL::ProjectileRemnant::reset().

193  {
194  std::stringstream ss;
195  ss << "Cluster (ID = " << ID << ") type = ";
197  ss << std::endl
198  << " A = " << theA << std::endl
199  << " Z = " << theZ << std::endl
200  << " mass = " << getMass() << std::endl
201  << " energy = " << theEnergy << std::endl
202  << " momentum = "
203  << theMomentum.print()
204  << std::endl
205  << " position = "
206  << thePosition.print()
207  << std::endl
208  << "Contains the following particles:"
209  << std::endl;
210  for(ParticleIter i=particles.begin(), e=particles.end(); i!=e; ++i)
211  ss << (*i)->print();
212  ss << std::endl;
213  return ss.str();
214  }
G4double getMass() const
Get the cached particle mass.
ParticleList particles
G4INCL::ThreeVector thePosition
G4INCL::ParticleType theType
G4INCL::ThreeVector theMomentum
std::string print() const
ParticleList::const_iterator ParticleIter
std::string getName(const ParticleType t)
Get the native INCL name of the particle.
void G4INCL::Cluster::putParticlesOffShell ( )
inline

Put the cluster components off shell.

The Cluster components are put off shell in such a way that their total energy equals the cluster mass.

Definition at line 277 of file G4INCLCluster.hh.

References energy(), INCL_DEBUG, G4INCL::ThreeVector::mag2(), particles, and print().

Referenced by G4INCL::ProjectileRemnant::ProjectileRemnant().

277  {
278  // Compute the dynamical potential
279  const G4double theDynamicalPotential = computeDynamicalPotential();
280  INCL_DEBUG("The dynamical potential is " << theDynamicalPotential << " MeV" << std::endl);
281 
282  for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
283  const G4double energy = (*p)->getEnergy() - theDynamicalPotential;
284  const ThreeVector &momentum = (*p)->getMomentum();
285  // Here particles are put off-shell so that we can satisfy the energy-
286  // and momentum-conservation laws
287  (*p)->setEnergy(energy);
288  (*p)->setMass(std::sqrt(energy*energy - momentum.mag2()));
289  }
290  INCL_DEBUG("Cluster components are now off shell:" << std::endl
291  << print());
292  }
const char * p
Definition: xmltok.h:285
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
std::string print() const
ParticleList particles
double G4double
Definition: G4Types.hh:76
#define INCL_DEBUG(x)
ParticleList::const_iterator ParticleIter
void G4INCL::Cluster::removeParticle ( Particle *const  p)
inline

Remove a particle from the cluster components.

Definition at line 166 of file G4INCLCluster.hh.

References particles, and G4INCL::UnorderedVector< T >::remove().

Referenced by G4INCL::ProjectileRemnant::removeParticle().

166 { particles.remove(p); }
const char * p
Definition: xmltok.h:285
void remove(const T &t)
ParticleList particles
virtual void G4INCL::Cluster::rotate ( const G4double  angle,
const ThreeVector axis 
)
inlinevirtual

Rotate position and momentum of all the particles.

This includes the cluster components. Overloads Particle::rotate().

Parameters
anglethe rotation angle
axisa unit vector representing the rotation axis

Reimplemented from G4INCL::Particle.

Definition at line 357 of file G4INCLCluster.hh.

References particles, and G4INCL::Particle::rotate().

357  {
358  Particle::rotate(angle, axis);
359  for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
360  (*p)->rotate(angle, axis);
361  }
362  }
const char * p
Definition: xmltok.h:285
virtual void rotate(const G4double angle, const ThreeVector &axis)
Rotate the particle position and momentum.
ParticleList particles
ParticleList::const_iterator ParticleIter
void G4INCL::Cluster::setA ( const G4int  A)
inline

Set the mass number of the cluster.

Definition at line 146 of file G4INCLCluster.hh.

References G4INCL::Particle::theA.

146 { theA = A; }
void G4INCL::Cluster::setExcitationEnergy ( const G4double  e)
inline

Set the excitation energy of the cluster.

Definition at line 152 of file G4INCLCluster.hh.

References theExcitationEnergy.

Referenced by G4INCL::Nucleus::finalizeProjectileRemnant().

152 { theExcitationEnergy=e; }
G4double theExcitationEnergy
void G4INCL::Cluster::setPosition ( const ThreeVector position)
inlinevirtual

Set the position of the cluster.

This overloads the Particle method to take into account that the positions of the cluster members must be updated as well.

Reimplemented from G4INCL::Particle.

Definition at line 299 of file G4INCLCluster.hh.

References particles, G4INCL::Particle::setPosition(), and G4INCL::Particle::thePosition.

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

299  {
300  ThreeVector shift(position-thePosition);
302  for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
303  (*p)->setPosition((*p)->getPosition()+shift);
304  }
305  }
const char * p
Definition: xmltok.h:285
ParticleList particles
virtual void setPosition(const G4INCL::ThreeVector &position)
G4INCL::ThreeVector thePosition
ParticleList::const_iterator ParticleIter
void G4INCL::Cluster::setSpin ( const ThreeVector j)
inline

Set the spin of the nucleus.

Definition at line 392 of file G4INCLCluster.hh.

References theSpin.

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

392 { theSpin = j; }
ThreeVector theSpin
void G4INCL::Cluster::setZ ( const G4int  Z)
inline

Set the charge number of the cluster.

Definition at line 143 of file G4INCLCluster.hh.

References G4INCL::Particle::theZ.

143 { theZ = Z; }
void G4INCL::Cluster::swap ( Cluster rhs)
inline

Helper method for the assignment operator.

Definition at line 119 of file G4INCLCluster.hh.

References particles, G4INCL::Particle::swap(), CLHEP::swap(), theExcitationEnergy, theParticleSampler, and theSpin.

Referenced by operator=().

119  {
120  Particle::swap(rhs);
121  std::swap(theExcitationEnergy, rhs.theExcitationEnergy);
122  std::swap(theSpin, rhs.theSpin);
123  // std::swap is overloaded by std::list and guaranteed to operate in
124  // constant time
125  std::swap(particles, rhs.particles);
126  std::swap(theParticleSampler, rhs.theParticleSampler);
127  }
ThreeVector theSpin
ParticleList particles
void swap(shared_ptr< P > &, shared_ptr< P > &)
Definition: memory.h:1247
ParticleSampler * theParticleSampler
G4double theExcitationEnergy
void swap(Particle &rhs)
Helper method for the assignment operator.

Field Documentation

ParticleList G4INCL::Cluster::particles
protected
G4double G4INCL::Cluster::theExcitationEnergy
protected
ParticleSampler* G4INCL::Cluster::theParticleSampler
protected
ThreeVector G4INCL::Cluster::theSpin
protected

Definition at line 420 of file G4INCLCluster.hh.

Referenced by G4INCL::Nucleus::computeRecoilKinematics(), getSpin(), setSpin(), and swap().


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