Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4INCL::Cluster Class Reference

#include <G4INCLCluster.hh>

Inheritance diagram for G4INCL::Cluster:

Public Member Functions

 Cluster (const G4int Z, const G4int A, const G4int S, const G4bool createParticleSampler=true)
 Standard Cluster constructor.
template<class Iterator>
 Cluster (Iterator begin, Iterator end)
virtual ~Cluster ()
 Cluster (const Cluster &rhs)
 Copy constructor.
Clusteroperator= (const Cluster &rhs)
 Assignment operator.
void swap (Cluster &rhs)
 Helper method for the assignment operator.
ParticleSpecies getSpecies () const
 Get the particle species.
void deleteParticles ()
void clearParticles ()
void setZ (const G4int Z)
 Set the charge number of the cluster.
void setA (const G4int A)
 Set the mass number of the cluster.
void setS (const G4int S)
 Set the strangess number of the cluster.
G4double getExcitationEnergy () const
 Get the excitation energy of the cluster.
void setExcitationEnergy (const G4double e)
 Set the excitation energy of the cluster.
virtual G4double getTableMass () const
 Get the real particle mass.
ParticleList const & getParticles () const
void removeParticle (Particle *const p)
 Remove a particle from the cluster components.
void addParticle (Particle *const p)
void updateClusterParameters ()
 Set total cluster mass, energy, size, etc. from the particles.
void addParticles (ParticleList const &pL)
 Add a list of particles to the cluster.
ParticleList getParticleList () const
 Returns the list of particles that make up the cluster.
std::string print () const
virtual void initializeParticles ()
 Initialise the NuclearDensity pointer and sample the particles.
void internalBoostToCM ()
 Boost to the CM of the component particles.
void putParticlesOffShell ()
 Put the cluster components off shell.
void setPosition (const ThreeVector &position)
 Set the position of the cluster.
void boost (const ThreeVector &aBoostVector)
 Boost the cluster with the indicated velocity.
void freezeInternalMotion ()
 Freeze the internal motion of the particles.
virtual void rotatePosition (const G4double angle, const ThreeVector &axis)
 Rotate position of all the particles.
virtual void rotateMomentum (const G4double angle, const ThreeVector &axis)
 Rotate momentum of all the particles.
virtual void makeProjectileSpectator ()
 Make all the components projectile spectators, too.
virtual void makeTargetSpectator ()
 Make all the components target spectators, too.
virtual void makeParticipant ()
 Make all the components participants, too.
ThreeVector const & getSpin () const
 Get the spin of the nucleus.
void setSpin (const ThreeVector &j)
 Set the spin of the nucleus.
G4INCL::ThreeVector getAngularMomentum () const
 Get the total angular momentum (orbital + spin).
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.
Particleoperator= (const Particle &rhs)
 Assignment operator.
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?
G4bool isEta () const
 Is this an eta?
G4bool isOmega () const
 Is this an omega?
G4bool isEtaPrime () const
 Is this an etaprime?
G4bool isPhoton () const
 Is this a photon?
G4bool isResonance () const
 Is it a resonance?
G4bool isDelta () const
 Is it a Delta?
G4bool isSigma () const
 Is this a Sigma?
G4bool isKaon () const
 Is this a Kaon?
G4bool isAntiKaon () const
 Is this an antiKaon?
G4bool isLambda () const
 Is this a Lambda?
G4bool isNucleonorLambda () const
 Is this a Nucleon or a Lambda?
G4bool isHyperon () const
 Is this an Hyperon?
G4bool isMeson () const
 Is this a Meson?
G4bool isBaryon () const
 Is this a Baryon?
G4bool isStrange () const
 Is this a Strange?
G4bool isXi () const
 Is this a Xi?
G4bool isAntiNucleon () const
 Is this an antinucleon?
G4bool isAntiSigma () const
 Is this an antiSigma?
G4bool isAntiXi () const
 Is this an antiXi?
G4bool isAntiLambda () const
 Is this an antiLambda?
G4bool isAntiHyperon () const
 Is this an antiHyperon?
G4bool isAntiBaryon () const
 Is this an antiBaryon?
G4bool isAntiNucleonorAntiLambda () const
 Is this an antiNucleon or an antiLambda?
G4int getA () const
 Returns the baryon number.
G4int getZ () const
 Returns the charge number.
G4int getS () const
 Returns the strangeness number.
G4int getSrcPair () const
 Returns the strangeness 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.
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 getEmissionQValueCorrection (const G4int AParent, const G4int ZParent, const G4int SParent) const
 Computes correction on the emission Q-value for hypernuclei.
G4double getTransferQValueCorrection (const G4int AFrom, const G4int ZFrom, const G4int SFrom, const G4int ATo, const G4int ZTo, const G4int STo) const
 Computes correction on the transfer Q-value for hypernuclei.
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 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.
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 setNumberOfSrcPair (int n)
 Set the number of srcpairs.
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 setSrcPartner ()
 Set and reset src partner.
void resetSrcPartner ()
G4bool isSrcPartner () const
 Check if the particle is a src partner.
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 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 rotatePositionAndMomentum (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
G4double getReflectionMomentum () const
 Return the reflection momentum.
void setUncorrelatedMomentum (const G4double p)
 Set the uncorrelated momentum.
void rpCorrelate ()
 Make the particle follow a strict r-p correlation.
void rpDecorrelate ()
 Make the particle not follow a strict r-p correlation.
G4double getCosRPAngle () const
 Get the cosine of the angle between position and momentum.
G4double getParticleBias () const
 Get the particle bias.
void setParticleBias (G4double ParticleBias)
 Set the particle bias.
std::vector< G4intgetBiasCollisionVector () const
 Get the vector list of biased vertices on the particle path.
void setBiasCollisionVector (std::vector< G4int > BiasCollisionVector)
 Set the vector list of biased vertices on the particle path.
G4int getNumberOfKaon () const
 Number of Kaon inside de nucleus.
void setNumberOfKaon (const G4int NK)
G4int getParentResonancePDGCode () const
void setParentResonancePDGCode (const G4int parentPDGCode)
G4int getParentResonanceID () const
void setParentResonanceID (const G4int parentID)

Protected Attributes

ParticleList particles
G4double theExcitationEnergy
ThreeVector theSpin
ParticleSamplertheParticleSampler
Protected Attributes inherited from G4INCL::Particle
G4int theZ
G4int theA
G4int theS
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
G4int nSrcPair
G4double thePotentialEnergy
long ID
G4bool rpCorrelated
G4double uncorrelatedMomentum
G4double theParticleBias
G4int theNKaon
 The number of Kaons inside the nucleus (update during the cascade).
G4int theParentResonancePDGCode
G4int theParentResonanceID

Additional Inherited Members

Static Public Member Functions inherited from G4INCL::Particle
static G4double getTotalBias ()
 General bias vector function.
static void setINCLBiasVector (std::vector< G4double > NewVector)
static void FillINCLBiasVector (G4double newBias)
static G4double getBiasFromVector (std::vector< G4int > VectorBias)
static std::vector< G4intMergeVectorBias (Particle const *const p1, Particle const *const p2)
static std::vector< G4intMergeVectorBias (std::vector< G4int > p1, Particle const *const p2)
Static Public Attributes inherited from G4INCL::Particle
static std::vector< G4doubleINCLBiasVector
 Time ordered vector of all bias applied.
static G4ThreadLocal G4int nextBiasedCollisionID = 0
Protected Member Functions inherited from G4INCL::Particle
void swap (Particle &rhs)
 Helper method for the assignment operator.

Detailed Description

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

Definition at line 52 of file G4INCLCluster.hh.

Constructor & Destructor Documentation

◆ Cluster() [1/3]

G4INCL::Cluster::Cluster ( const G4int Z,
const G4int A,
const G4int S,
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 60 of file G4INCLCluster.hh.

60 :
61 Particle(),
63 theSpin(0.,0.,0.),
65 {
66 if(A >= 0){
68 theZ = Z;
69 theA = A;
70 theS = S;
72 if(createParticleSampler)
73 theParticleSampler = new ParticleSampler(A,Z,S);
74 }
75 else {
77 theZ = Z;
78 theA = A;
79 theS = S;
81 if(createParticleSampler)
82 theParticleSampler = new ParticleSampler(A,Z,S);
83 }
84 }
G4double S(G4double temp)
const G4double A[17]
ParticleSampler * theParticleSampler
ThreeVector theSpin
G4double theExcitationEnergy
void setINCLMass()
Set the mass of the Particle to its table mass.
void setType(ParticleType t)

Referenced by G4INCL::Nucleus::applyFinalState(), Cluster(), G4INCL::Nucleus::decayOutgoingClusters(), G4INCL::Nucleus::Nucleus(), operator=(), G4INCL::ProjectileRemnant::ProjectileRemnant(), and swap().

◆ Cluster() [2/3]

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 90 of file G4INCLCluster.hh.

90 :
91 Particle(),
93 theSpin(0.,0.,0.),
95 {
97 for(Iterator i = begin; i != end; ++i) {
98 addParticle(*i);
99 }
100 if (theA < 0){
102 thePosition /= (-theA);
103 }
104 else
105 thePosition /= theA;
106 setINCLMass();
108 }
void addParticle(Particle *const p)
const ThreeVector & adjustMomentumFromEnergy()
Rescale the momentum to match the total energy.
G4INCL::ThreeVector thePosition

◆ ~Cluster()

virtual G4INCL::Cluster::~Cluster ( )
inlinevirtual

Definition at line 110 of file G4INCLCluster.hh.

110 {
111 delete theParticleSampler;
112 }

◆ Cluster() [3/3]

G4INCL::Cluster::Cluster ( const Cluster & rhs)
inline

Copy constructor.

Definition at line 115 of file G4INCLCluster.hh.

115 :
116 Particle(rhs),
117 theExcitationEnergy(rhs.theExcitationEnergy),
118 theSpin(rhs.theSpin)
119 {
120 for(ParticleIter p=rhs.particles.begin(), e=rhs.particles.end(); p!=e; ++p) {
121 particles.push_back(new Particle(**p));
122 }
123 if(rhs.theParticleSampler)
124 theParticleSampler = new ParticleSampler(rhs.theA,rhs.theZ,rhs.theS);
125 else
126 theParticleSampler = NULL;
127 }
ParticleList particles
ParticleList::const_iterator ParticleIter

Member Function Documentation

◆ addParticle()

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 194 of file G4INCLCluster.hh.

194 {
195 particles.push_back(p);
196 theEnergy += p->getEnergy();
197 thePotentialEnergy += p->getPotentialEnergy();
198 theMomentum += p->getMomentum();
199 thePosition += p->getPosition();
200 theA += p->getA();
201 theZ += p->getZ();
202 theS += p->getS();
203 nCollisions += p->getNumberOfCollisions();
204 }
G4INCL::ThreeVector theMomentum

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

◆ addParticles()

void G4INCL::Cluster::addParticles ( ParticleList const & pL)
inline

Add a list of particles to the cluster.

Definition at line 229 of file G4INCLCluster.hh.

229 {
230 particles = pL;
232 }
void updateClusterParameters()
Set total cluster mass, energy, size, etc. from the particles.

◆ boost()

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 373 of file G4INCLCluster.hh.

373 {
374 Particle::boost(aBoostVector);
375 for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
376 (*p)->boost(aBoostVector);
377 // Apply Lorentz contraction to the particle position
378 (*p)->lorentzContract(aBoostVector,thePosition);
379 (*p)->rpCorrelate();
380 }
381
382 INCL_DEBUG("Cluster was boosted with (bx,by,bz)=("
383 << aBoostVector.getX() << ", " << aBoostVector.getY() << ", " << aBoostVector.getZ() << "):"
384 << '\n' << print());
385
386 }
#define INCL_DEBUG(x)
std::string print() const
void boost(const ThreeVector &aBoostVector)

Referenced by G4INCL::ProjectileRemnant::ProjectileRemnant(), and G4INCL::StandardPropagationModel::shootCompositeAtrest().

◆ clearParticles()

void G4INCL::Cluster::clearParticles ( )
inline

Definition at line 159 of file G4INCLCluster.hh.

159{ particles.clear(); }

Referenced by deleteParticles().

◆ deleteParticles()

void G4INCL::Cluster::deleteParticles ( )
inline

Definition at line 152 of file G4INCLCluster.hh.

152 {
153 for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
154 delete (*p);
155 }
157 }

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

◆ freezeInternalMotion()

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 396 of file G4INCLCluster.hh.

396 {
397 const ThreeVector &normMomentum = theMomentum / getMass();
398 for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
399 const G4double pMass = (*p)->getMass();
400 const ThreeVector frozenMomentum = normMomentum * pMass;
401 const G4double frozenEnergy = std::sqrt(frozenMomentum.mag2()+pMass*pMass);
402 (*p)->setFrozenMomentum(frozenMomentum);
403 (*p)->setFrozenEnergy(frozenEnergy);
404 (*p)->freezePropagation();
405 }
406 }
double G4double
Definition G4Types.hh:83
G4double getMass() const
Get the cached particle mass.

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

◆ getAngularMomentum()

G4INCL::ThreeVector G4INCL::Cluster::getAngularMomentum ( ) const
inlinevirtual

Get the total angular momentum (orbital + spin).

Reimplemented from G4INCL::Particle.

Definition at line 457 of file G4INCLCluster.hh.

457 {
459 }
ThreeVector const & getSpin() const
Get the spin of the nucleus.
virtual G4INCL::ThreeVector getAngularMomentum() const

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

◆ getExcitationEnergy()

G4double G4INCL::Cluster::getExcitationEnergy ( ) const
inline

Get the excitation energy of the cluster.

Definition at line 171 of file G4INCLCluster.hh.

171{ return theExcitationEnergy; }

◆ getParticleList()

ParticleList G4INCL::Cluster::getParticleList ( ) const
inline

Returns the list of particles that make up the cluster.

Definition at line 235 of file G4INCLCluster.hh.

235{ return particles; }

◆ getParticles()

ParticleList const & G4INCL::Cluster::getParticles ( ) const
inline

◆ getSpecies()

ParticleSpecies G4INCL::Cluster::getSpecies ( ) const
inlinevirtual

Get the particle species.

Reimplemented from G4INCL::Particle.

Definition at line 148 of file G4INCLCluster.hh.

148 {
149 return ParticleSpecies(theA, theZ, theS);
150 }

◆ getSpin()

ThreeVector const & G4INCL::Cluster::getSpin ( ) const
inline

Get the spin of the nucleus.

Definition at line 451 of file G4INCLCluster.hh.

451{ return theSpin; }

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

◆ getTableMass()

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

Get the real particle mass.

Overloads the Particle method.

Reimplemented from G4INCL::Particle.

Definition at line 180 of file G4INCLCluster.hh.

180{ return getRealMass(); }
G4double getRealMass() const
Get the real particle mass.

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

◆ initializeParticles()

void G4INCL::Cluster::initializeParticles ( )
virtual

Initialise the NuclearDensity pointer and sample the particles.

Reimplemented in G4INCL::Nucleus.

Definition at line 43 of file G4INCLCluster.cc.

43 {
44// assert(theA>=2 || theA<=-2);
45 const ThreeVector oldPosition = thePosition;
46 theParticleSampler->sampleParticlesIntoList(thePosition, particles);
47#if !defined(NDEBUG) && !defined(INCLXX_IN_GEANT4_MODE)
48 const G4int theMassNumber = theA;
49 const G4int theChargeNumber = theZ;
50#endif
52 thePosition = oldPosition;
53// assert(theMassNumber==theA && theChargeNumber==theZ);
54 INCL_DEBUG("Cluster initialized:" << '\n' << print());
55 }
int G4int
Definition G4Types.hh:85

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

◆ internalBoostToCM()

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 270 of file G4INCLCluster.hh.

270 {
271
272 // First compute the current CM position and total momentum
273 ThreeVector theCMPosition, theTotalMomentum;
274 for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
275 theCMPosition += (*p)->getPosition();
276 theTotalMomentum += (*p)->getMomentum();
277 //theTotalEnergy += (*p)->getEnergy();
278 }
279 if(theA>=0){
280 theCMPosition /= theA;
281// assert((unsigned int)theA==particles.size());
282 } else if (theA < 0){
283 theCMPosition /= -theA;
284//assert(-theA==particles.size());
285 }
286
287 // Now determine the CM velocity of the particles
288 // commented out because currently unused, see below
289 // ThreeVector betaCM = theTotalMomentum / theTotalEnergy;
290
291 // The new particle positions and momenta are scaled by a factor of
292 // \f$\sqrt{A/(A-1)}\f$, so that the resulting density distributions in
293 // the CM have the same variance as the one we started with.
294 G4double rescaling;
295 if (theA>0)
296 rescaling = std::sqrt(((G4double)theA)/((G4double)(theA-1)));
297 else if (theA<0)
298 rescaling = std::sqrt(((G4double)(-theA))/((G4double)((-theA)-1)));
299 else
300 rescaling = 0 ;
301
302 // Loop again to boost and reposition
303 for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
304 // \bug{We should do the following, but the Fortran version actually
305 // does not!
306 // (*p)->boost(betaCM);
307 // Here is what the Fortran version does:}
308 if (theA>0)
309 (*p)->setMomentum(((*p)->getMomentum()-theTotalMomentum/theA)*rescaling);
310 else if (theA<0)
311 (*p)->setMomentum(((*p)->getMomentum()-theTotalMomentum/(-theA))*rescaling);
312
313 // Set the CM position of the particles
314 (*p)->setPosition(((*p)->getPosition()-theCMPosition)*rescaling);
315 }
316
317 // Set the global cluster kinematic variables
318 thePosition.setX(0.0);
319 thePosition.setY(0.0);
320 thePosition.setZ(0.0);
321 theMomentum.setX(0.0);
322 theMomentum.setY(0.0);
323 theMomentum.setZ(0.0);
324 theEnergy = getMass();
325
326 INCL_DEBUG("Cluster boosted to internal CM:" << '\n' << print());
327
328 }

Referenced by G4INCL::ProjectileRemnant::ProjectileRemnant(), and G4INCL::StandardPropagationModel::shootCompositeAtrest().

◆ makeParticipant()

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

Make all the components participants, too.

Reimplemented from G4INCL::Particle.

Definition at line 443 of file G4INCLCluster.hh.

443 {
445 for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
446 (*p)->makeParticipant();
447 }
448 }
virtual void makeParticipant()

◆ makeProjectileSpectator()

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

Make all the components projectile spectators, too.

Reimplemented from G4INCL::Particle.

Definition at line 427 of file G4INCLCluster.hh.

427 {
429 for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
430 (*p)->makeProjectileSpectator();
431 }
432 }
virtual void makeProjectileSpectator()

Referenced by G4INCL::ProjectileRemnant::ProjectileRemnant(), and G4INCL::StandardPropagationModel::shootCompositeAtrest().

◆ makeTargetSpectator()

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

Make all the components target spectators, too.

Reimplemented from G4INCL::Particle.

Definition at line 435 of file G4INCLCluster.hh.

435 {
437 for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
438 (*p)->makeTargetSpectator();
439 }
440 }
virtual void makeTargetSpectator()

◆ operator=()

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

Assignment operator.

Definition at line 130 of file G4INCLCluster.hh.

130 {
131 Cluster temporaryCluster(rhs);
132 Particle::operator=(temporaryCluster);
133 swap(temporaryCluster);
134 return *this;
135 }
Cluster(const G4int Z, const G4int A, const G4int S, const G4bool createParticleSampler=true)
Standard Cluster constructor.
void swap(Cluster &rhs)
Helper method for the assignment operator.
Particle & operator=(const Particle &rhs)
Assignment operator.

◆ print()

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

Definition at line 237 of file G4INCLCluster.hh.

237 {
238 std::stringstream ss;
239 ss << "Cluster (ID = " << ID << ") type = ";
241 ss << '\n'
242 << " A = " << theA << '\n'
243 << " Z = " << theZ << '\n'
244 << " S = " << theS << '\n'
245 << " mass = " << getMass() << '\n'
246 << " energy = " << theEnergy << '\n'
247 << " momentum = "
248 << theMomentum.print()
249 << '\n'
250 << " position = "
251 << thePosition.print()
252 << '\n'
253 << "Contains the following particles:"
254 << '\n';
255 for(ParticleIter i=particles.begin(), e=particles.end(); i!=e; ++i)
256 ss << (*i)->print();
257 ss << '\n';
258 return ss.str();
259 }
G4INCL::ParticleType theType
std::string getName(const ParticleType t)
Get the native INCL name of the particle.

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

◆ putParticlesOffShell()

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 335 of file G4INCLCluster.hh.

335 {
336 // Compute the dynamical potential
337 const G4double theDynamicalPotential = computeDynamicalPotential();
338 INCL_DEBUG("The dynamical potential is " << theDynamicalPotential << " MeV" << '\n');
339
340 for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
341 const G4double energy = (*p)->getEnergy() - theDynamicalPotential;
342 const ThreeVector &momentum = (*p)->getMomentum();
343 // Here particles are put off-shell so that we can satisfy the energy-
344 // and momentum-conservation laws
345 (*p)->setEnergy(energy);
346 (*p)->setMass(std::sqrt(energy*energy - momentum.mag2()));
347 }
348 INCL_DEBUG("Cluster components are now off shell:" << '\n'
349 << print());
350 }
G4double energy(const ThreeVector &p, const G4double m)

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

◆ removeParticle()

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

Remove a particle from the cluster components.

Definition at line 188 of file G4INCLCluster.hh.

188{ particles.remove(p); }

Referenced by G4INCL::ProjectileRemnant::removeParticle(), and G4INCL::StandardPropagationModel::shootCompositeAtrest().

◆ rotateMomentum()

void G4INCL::Cluster::rotateMomentum ( const G4double angle,
const ThreeVector & axis )
virtual

Rotate momentum of all the particles.

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

Parameters
anglethe rotation angle
axisa unit vector representing the rotation axis

Reimplemented from G4INCL::Particle.

Definition at line 64 of file G4INCLCluster.cc.

64 {
66 for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
67 (*p)->rotateMomentum(angle, axis);
68 }
69 }
virtual void rotateMomentum(const G4double angle, const ThreeVector &axis)
Rotate the particle momentum.
const axis_t axis_to_type< N >::axis
Definition pugixml.cc:9668

◆ rotatePosition()

void G4INCL::Cluster::rotatePosition ( const G4double angle,
const ThreeVector & axis )
virtual

Rotate position of all the particles.

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

Parameters
anglethe rotation angle
axisa unit vector representing the rotation axis

Reimplemented from G4INCL::Particle.

Definition at line 57 of file G4INCLCluster.cc.

57 {
59 for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
60 (*p)->rotatePosition(angle, axis);
61 }
62 }
virtual void rotatePosition(const G4double angle, const ThreeVector &axis)
Rotate the particle position.

◆ setA()

void G4INCL::Cluster::setA ( const G4int A)
inline

Set the mass number of the cluster.

Definition at line 165 of file G4INCLCluster.hh.

165{ theA = A; }

◆ setExcitationEnergy()

void G4INCL::Cluster::setExcitationEnergy ( const G4double e)
inline

Set the excitation energy of the cluster.

Definition at line 174 of file G4INCLCluster.hh.

◆ setPosition()

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 357 of file G4INCLCluster.hh.

357 {
358 ThreeVector shift(position-thePosition);
359 Particle::setPosition(position);
360 for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
361 (*p)->setPosition((*p)->getPosition()+shift);
362 }
363 }
virtual void setPosition(const G4INCL::ThreeVector &position)

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

◆ setS()

void G4INCL::Cluster::setS ( const G4int S)
inline

Set the strangess number of the cluster.

Definition at line 168 of file G4INCLCluster.hh.

168{ theS = S; }

◆ setSpin()

void G4INCL::Cluster::setSpin ( const ThreeVector & j)
inline

Set the spin of the nucleus.

Definition at line 454 of file G4INCLCluster.hh.

454{ theSpin = j; }

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

◆ setZ()

void G4INCL::Cluster::setZ ( const G4int Z)
inline

Set the charge number of the cluster.

Definition at line 162 of file G4INCLCluster.hh.

162{ theZ = Z; }

◆ swap()

void G4INCL::Cluster::swap ( Cluster & rhs)
inline

Helper method for the assignment operator.

Definition at line 138 of file G4INCLCluster.hh.

138 {
139 Particle::swap(rhs);
140 std::swap(theExcitationEnergy, rhs.theExcitationEnergy);
141 std::swap(theSpin, rhs.theSpin);
142 // std::swap is overloaded by std::list and guaranteed to operate in
143 // constant time
144 std::swap(particles, rhs.particles);
145 std::swap(theParticleSampler, rhs.theParticleSampler);
146 }
void swap(Particle &rhs)
Helper method for the assignment operator.

Referenced by operator=().

◆ updateClusterParameters()

void G4INCL::Cluster::updateClusterParameters ( )
inline

Set total cluster mass, energy, size, etc. from the particles.

Definition at line 207 of file G4INCLCluster.hh.

207 {
208 theEnergy = 0.;
210 theMomentum = ThreeVector();
211 thePosition = ThreeVector();
212 theA = 0;
213 theZ = 0;
214 theS = 0;
215 nCollisions = 0;
216 for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
217 theEnergy += (*p)->getEnergy();
218 thePotentialEnergy += (*p)->getPotentialEnergy();
219 theMomentum += (*p)->getMomentum();
220 thePosition += (*p)->getPosition();
221 theA += (*p)->getA();
222 theZ += (*p)->getZ();
223 theS += (*p)->getS();
224 nCollisions += (*p)->getNumberOfCollisions();
225 }
226 }

Referenced by addParticles(), and initializeParticles().

Member Data Documentation

◆ particles

◆ theExcitationEnergy

◆ theParticleSampler

ParticleSampler* G4INCL::Cluster::theParticleSampler
protected

◆ theSpin

ThreeVector G4INCL::Cluster::theSpin
protected

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