Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4INCLCluster.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// INCL++ intra-nuclear cascade model
27// Alain Boudard, CEA-Saclay, France
28// Joseph Cugnon, University of Liege, Belgium
29// Jean-Christophe David, CEA-Saclay, France
30// Pekka Kaitaniemi, CEA-Saclay, France, and Helsinki Institute of Physics, Finland
31// Sylvie Leray, CEA-Saclay, France
32// Davide Mancusi, CEA-Saclay, France
33//
34#define INCLXX_IN_GEANT4_MODE 1
35
36#include "globals.hh"
37
38#ifndef G4INCLCluster_hh
39#define G4INCLCluster_hh 1
40
41#include "G4INCLParticle.hh"
45
46namespace G4INCL {
47
48 /**
49 * Cluster is a particle (inherits from the Particle class) that is
50 * actually a collection of elementary particles.
51 */
52 class Cluster : public Particle {
53 public:
54
55 /** \brief Standard Cluster constructor
56 *
57 * This constructor should mainly be used when constructing Nucleus or
58 * when constructing Clusters to be used as composite projectiles.
59 */
60 Cluster(const G4int Z, const G4int A, const G4int S, const G4bool createParticleSampler=true) :
61 Particle(),
63 theSpin(0.,0.,0.),
65 {
67 theZ = Z;
68 theA = A;
69 theS = S;
71 if(createParticleSampler)
73 }
74
75 /**
76 * A cluster can be directly built from a list of particles.
77 */
78 template<class Iterator>
79 Cluster(Iterator begin, Iterator end) :
80 Particle(),
82 theSpin(0.,0.,0.),
84 {
86 for(Iterator i = begin; i != end; ++i) {
87 addParticle(*i);
88 }
92 }
93
94 virtual ~Cluster() {
95 delete theParticleSampler;
96 }
97
98 /// \brief Copy constructor
99 Cluster(const Cluster &rhs) :
100 Particle(rhs),
102 theSpin(rhs.theSpin)
103 {
104 for(ParticleIter p=rhs.particles.begin(), e=rhs.particles.end(); p!=e; ++p) {
105 particles.push_back(new Particle(**p));
106 }
107 if(rhs.theParticleSampler)
109 else
110 theParticleSampler = NULL;
111 }
112
113 /// \brief Assignment operator
114 Cluster &operator=(const Cluster &rhs) {
115 Cluster temporaryCluster(rhs);
116 Particle::operator=(temporaryCluster);
117 swap(temporaryCluster);
118 return *this;
119 }
120
121 /// \brief Helper method for the assignment operator
122 void swap(Cluster &rhs) {
123 Particle::swap(rhs);
125 std::swap(theSpin, rhs.theSpin);
126 // std::swap is overloaded by std::list and guaranteed to operate in
127 // constant time
128 std::swap(particles, rhs.particles);
130 }
131
133 return ParticleSpecies(theA, theZ, theS);
134 }
135
137 for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
138 delete (*p);
139 }
141 }
142
143 void clearParticles() { particles.clear(); }
144
145 /// \brief Set the charge number of the cluster
146 void setZ(const G4int Z) { theZ = Z; }
147
148 /// \brief Set the mass number of the cluster
149 void setA(const G4int A) { theA = A; }
150
151 /// \brief Set the strangess number of the cluster
152 void setS(const G4int S) { theS = S; }
153
154 /// \brief Get the excitation energy of the cluster.
156
157 /// \brief Set the excitation energy of the cluster.
159
160 /** \brief Get the real particle mass.
161 *
162 * Overloads the Particle method.
163 */
164 inline virtual G4double getTableMass() const { return getRealMass(); }
165
166 /**
167 * Get the list of particles in the cluster.
168 */
169 ParticleList const &getParticles() const { return particles; }
170
171 /// \brief Remove a particle from the cluster components.
172 void removeParticle(Particle * const p) { particles.remove(p); }
173
174 /**
175 * Add one particle to the cluster. This updates the cluster mass,
176 * energy, size, etc.
177 */
178 void addParticle(Particle * const p) {
179 particles.push_back(p);
180 theEnergy += p->getEnergy();
182 theMomentum += p->getMomentum();
183 thePosition += p->getPosition();
184 theA += p->getA();
185 theZ += p->getZ();
186 theS += p->getS();
188 }
189
190 /// \brief Set total cluster mass, energy, size, etc. from the particles
192 theEnergy = 0.;
196 theA = 0;
197 theZ = 0;
198 theS = 0;
199 nCollisions = 0;
200 for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
201 theEnergy += (*p)->getEnergy();
202 thePotentialEnergy += (*p)->getPotentialEnergy();
203 theMomentum += (*p)->getMomentum();
204 thePosition += (*p)->getPosition();
205 theA += (*p)->getA();
206 theZ += (*p)->getZ();
207 theS += (*p)->getS();
208 nCollisions += (*p)->getNumberOfCollisions();
209 }
210 }
211
212 /// \brief Add a list of particles to the cluster
213 void addParticles(ParticleList const &pL) {
214 particles = pL;
216 }
217
218 /// \brief Returns the list of particles that make up the cluster
220
221 std::string print() const {
222 std::stringstream ss;
223 ss << "Cluster (ID = " << ID << ") type = ";
225 ss << '\n'
226 << " A = " << theA << '\n'
227 << " Z = " << theZ << '\n'
228 << " S = " << theS << '\n'
229 << " mass = " << getMass() << '\n'
230 << " energy = " << theEnergy << '\n'
231 << " momentum = "
232 << theMomentum.print()
233 << '\n'
234 << " position = "
235 << thePosition.print()
236 << '\n'
237 << "Contains the following particles:"
238 << '\n';
239 for(ParticleIter i=particles.begin(), e=particles.end(); i!=e; ++i)
240 ss << (*i)->print();
241 ss << '\n';
242 return ss.str();
243 }
244
245 /// \brief Initialise the NuclearDensity pointer and sample the particles
246 virtual void initializeParticles();
247
248 /** \brief Boost to the CM of the component particles
249 *
250 * The position of all particles in the particles list is shifted so that
251 * their centre of mass is in the origin and their total momentum is
252 * zero.
253 */
255
256 // First compute the current CM position and total momentum
257 ThreeVector theCMPosition, theTotalMomentum;
258 G4double theTotalEnergy = 0.0;
259 for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
260 theCMPosition += (*p)->getPosition();
261 theTotalMomentum += (*p)->getMomentum();
262 theTotalEnergy += (*p)->getEnergy();
263 }
264 theCMPosition /= theA;
265// assert((unsigned int)theA==particles.size());
266
267 // Now determine the CM velocity of the particles
268 // commented out because currently unused, see below
269 // ThreeVector betaCM = theTotalMomentum / theTotalEnergy;
270
271 // The new particle positions and momenta are scaled by a factor of
272 // \f$\sqrt{A/(A-1)}\f$, so that the resulting density distributions in
273 // the CM have the same variance as the one we started with.
274 const G4double rescaling = std::sqrt(((G4double)theA)/((G4double)(theA-1)));
275
276 // Loop again to boost and reposition
277 for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
278 // \bug{We should do the following, but the Fortran version actually
279 // does not!
280 // (*p)->boost(betaCM);
281 // Here is what the Fortran version does:}
282 (*p)->setMomentum(((*p)->getMomentum()-theTotalMomentum/theA)*rescaling);
283
284 // Set the CM position of the particles
285 (*p)->setPosition(((*p)->getPosition()-theCMPosition)*rescaling);
286 }
287
288 // Set the global cluster kinematic variables
289 thePosition.setX(0.0);
290 thePosition.setY(0.0);
291 thePosition.setZ(0.0);
292 theMomentum.setX(0.0);
293 theMomentum.setY(0.0);
294 theMomentum.setZ(0.0);
295 theEnergy = getMass();
296
297 INCL_DEBUG("Cluster boosted to internal CM:" << '\n' << print());
298
299 }
300
301 /** \brief Put the cluster components off shell
302 *
303 * The Cluster components are put off shell in such a way that their total
304 * energy equals the cluster mass.
305 */
307 // Compute the dynamical potential
308 const G4double theDynamicalPotential = computeDynamicalPotential();
309 INCL_DEBUG("The dynamical potential is " << theDynamicalPotential << " MeV" << '\n');
310
311 for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
312 const G4double energy = (*p)->getEnergy() - theDynamicalPotential;
313 const ThreeVector &momentum = (*p)->getMomentum();
314 // Here particles are put off-shell so that we can satisfy the energy-
315 // and momentum-conservation laws
316 (*p)->setEnergy(energy);
317 (*p)->setMass(std::sqrt(energy*energy - momentum.mag2()));
318 }
319 INCL_DEBUG("Cluster components are now off shell:" << '\n'
320 << print());
321 }
322
323 /** \brief Set the position of the cluster
324 *
325 * This overloads the Particle method to take into account that the
326 * positions of the cluster members must be updated as well.
327 */
331 for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
332 (*p)->setPosition((*p)->getPosition()+shift);
333 }
334 }
335
336 /** \brief Boost the cluster with the indicated velocity
337 *
338 * The Cluster is boosted as a whole, just like any Particle object;
339 * moreover, the internal components (particles list) are also boosted,
340 * according to Alain Boudard's off-shell recipe.
341 *
342 * \param aBoostVector the velocity to boost to [c]
343 */
344 void boost(const ThreeVector &aBoostVector) {
345 Particle::boost(aBoostVector);
346 for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
347 (*p)->boost(aBoostVector);
348 // Apply Lorentz contraction to the particle position
349 (*p)->lorentzContract(aBoostVector,thePosition);
350 (*p)->rpCorrelate();
351 }
352
353 INCL_DEBUG("Cluster was boosted with (bx,by,bz)=("
354 << aBoostVector.getX() << ", " << aBoostVector.getY() << ", " << aBoostVector.getZ() << "):"
355 << '\n' << print());
356
357 }
358
359 /** \brief Freeze the internal motion of the particles
360 *
361 * Each particle is assigned a frozen momentum four-vector determined by
362 * the collective cluster velocity. This is used for propagation, but not
363 * for dynamics. Normal propagation is restored by calling the
364 * Particle::thawPropagation() method, which should be done in
365 * InteractionAvatar::postInteraction.
366 */
368 const ThreeVector &normMomentum = theMomentum / getMass();
369 for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
370 const G4double pMass = (*p)->getMass();
371 const ThreeVector frozenMomentum = normMomentum * pMass;
372 const G4double frozenEnergy = std::sqrt(frozenMomentum.mag2()+pMass*pMass);
373 (*p)->setFrozenMomentum(frozenMomentum);
374 (*p)->setFrozenEnergy(frozenEnergy);
375 (*p)->freezePropagation();
376 }
377 }
378
379 /** \brief Rotate position of all the particles
380 *
381 * This includes the cluster components. Overloads Particle::rotateMomentum().
382 *
383 * \param angle the rotation angle
384 * \param axis a unit vector representing the rotation axis
385 */
386 virtual void rotatePosition(const G4double angle, const ThreeVector &axis);
387
388 /** \brief Rotate momentum of all the particles
389 *
390 * This includes the cluster components. Overloads Particle::rotateMomentum().
391 *
392 * \param angle the rotation angle
393 * \param axis a unit vector representing the rotation axis
394 */
395 virtual void rotateMomentum(const G4double angle, const ThreeVector &axis);
396
397 /// \brief Make all the components projectile spectators, too
398 virtual void makeProjectileSpectator() {
400 for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
401 (*p)->makeProjectileSpectator();
402 }
403 }
404
405 /// \brief Make all the components target spectators, too
406 virtual void makeTargetSpectator() {
408 for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
409 (*p)->makeTargetSpectator();
410 }
411 }
412
413 /// \brief Make all the components participants, too
414 virtual void makeParticipant() {
416 for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
417 (*p)->makeParticipant();
418 }
419 }
420
421 /// \brief Get the spin of the nucleus.
422 ThreeVector const &getSpin() const { return theSpin; }
423
424 /// \brief Set the spin of the nucleus.
425 void setSpin(const ThreeVector &j) { theSpin = j; }
426
427 /// \brief Get the total angular momentum (orbital + spin)
430 }
431
432 private:
433 /** \brief Compute the dynamical cluster potential
434 *
435 * Alain Boudard's boost prescription for low-energy beams requires to
436 * define a "dynamical potential" that allows us to conserve momentum and
437 * energy when boosting the projectile cluster.
438 */
439 G4double computeDynamicalPotential() {
440 G4double theDynamicalPotential = 0.0;
441 for(ParticleIter p=particles.begin(), e=particles.end(); p!=e; ++p) {
442 theDynamicalPotential += (*p)->getEnergy();
443 }
444 theDynamicalPotential -= getTableMass();
445 theDynamicalPotential /= theA;
446
447 return theDynamicalPotential;
448 }
449
450 protected:
455
457 };
458
459}
460
461#endif
double S(double temp)
double A(double temperature)
Singleton for recycling allocation of instances of a given class.
#define INCL_DECLARE_ALLOCATION_POOL(T)
#define INCL_DEBUG(x)
Class for sampling particles in a nucleus.
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
ThreeVector const & getSpin() const
Get the spin of the nucleus.
void boost(const ThreeVector &aBoostVector)
Boost the cluster with the indicated velocity.
void internalBoostToCM()
Boost to the CM of the component particles.
void removeParticle(Particle *const p)
Remove a particle from the cluster components.
ParticleSpecies getSpecies() const
Get the particle species.
virtual ~Cluster()
ParticleList particles
ParticleList const & getParticles() const
virtual void makeProjectileSpectator()
Make all the components projectile spectators, too.
void setExcitationEnergy(const G4double e)
Set the excitation energy of the cluster.
void setS(const G4int S)
Set the strangess number of the cluster.
virtual void rotateMomentum(const G4double angle, const ThreeVector &axis)
Rotate momentum of all the particles.
void setSpin(const ThreeVector &j)
Set the spin of the nucleus.
virtual void initializeParticles()
Initialise the NuclearDensity pointer and sample the particles.
ParticleSampler * theParticleSampler
G4INCL::ThreeVector getAngularMomentum() const
Get the total angular momentum (orbital + spin)
void addParticles(ParticleList const &pL)
Add a list of particles to the cluster.
Cluster(const G4int Z, const G4int A, const G4int S, const G4bool createParticleSampler=true)
Standard Cluster constructor.
virtual void makeTargetSpectator()
Make all the components target spectators, too.
ThreeVector theSpin
G4double getExcitationEnergy() const
Get the excitation energy of the cluster.
virtual G4double getTableMass() const
Get the real particle mass.
void swap(Cluster &rhs)
Helper method for the assignment operator.
virtual void rotatePosition(const G4double angle, const ThreeVector &axis)
Rotate position of all the particles.
Cluster & operator=(const Cluster &rhs)
Assignment operator.
void setPosition(const ThreeVector &position)
Set the position of the cluster.
void setA(const G4int A)
Set the mass number of the cluster.
Cluster(Iterator begin, Iterator end)
Cluster(const Cluster &rhs)
Copy constructor.
void addParticle(Particle *const p)
virtual void makeParticipant()
Make all the components participants, too.
G4double theExcitationEnergy
void updateClusterParameters()
Set total cluster mass, energy, size, etc. from the particles.
void freezeInternalMotion()
Freeze the internal motion of the particles.
void setZ(const G4int Z)
Set the charge number of the cluster.
ParticleList getParticleList() const
Returns the list of particles that make up the cluster.
std::string print() const
void putParticlesOffShell()
Put the cluster components off shell.
G4int getS() const
Returns the strangeness number.
G4INCL::ThreeVector theMomentum
virtual G4INCL::ThreeVector getAngularMomentum() const
G4double getEnergy() const
G4double getPotentialEnergy() const
Get the particle potential energy.
virtual void makeTargetSpectator()
G4double thePotentialEnergy
G4int getZ() const
Returns the charge number.
const G4INCL::ThreeVector & getPosition() const
G4INCL::ParticleType theType
G4int getNumberOfCollisions() const
Return the number of collisions undergone by the particle.
virtual void makeProjectileSpectator()
const ThreeVector & adjustMomentumFromEnergy()
Rescale the momentum to match the total energy.
void setINCLMass()
Set the mass of the Particle to its table mass.
G4double getRealMass() const
Get the real particle mass.
const G4INCL::ThreeVector & getMomentum() const
Particle & operator=(const Particle &rhs)
Assignment operator.
virtual void makeParticipant()
void setType(ParticleType t)
G4double getMass() const
Get the cached particle mass.
void boost(const ThreeVector &aBoostVector)
virtual void setPosition(const G4INCL::ThreeVector &position)
void swap(Particle &rhs)
Helper method for the assignment operator.
G4int getA() const
Returns the baryon number.
G4INCL::ThreeVector thePosition
void setY(G4double ay)
Set the y coordinate.
G4double getY() const
G4double getZ() const
void setZ(G4double az)
Set the z coordinate.
void setX(G4double ax)
Set the x coordinate.
std::string print() const
G4double mag2() const
G4double getX() const
std::string getName(const ParticleType t)
Get the native INCL name of the particle.
ParticleList::const_iterator ParticleIter