Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4INCLStandardPropagationModel.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// Pekka Kaitaniemi, CEA and Helsinki Institute of Physics
28// Davide Mancusi, CEA
29// Alain Boudard, CEA
30// Sylvie Leray, CEA
31// Joseph Cugnon, University of Liege
32//
33// INCL++ revision: v5.1.8
34//
35#define INCLXX_IN_GEANT4_MODE 1
36
37#include "globals.hh"
38
39/*
40 * StandardPropagationModel.hh
41 *
42 * \date 4 June 2009
43 * \author Pekka Kaitaniemi
44 */
45
46#ifndef G4INCLStandardPropagationModel_hh
47#define G4INCLStandardPropagationModel_hh 1
48
49#include "G4INCLNucleus.hh"
51#include "G4INCLIAvatar.hh"
52#include "G4INCLConfigEnums.hh"
53
54#include <iterator>
55
56namespace G4INCL {
57
58 /**
59 * Standard INCL4 particle propagation and avatar prediction
60 *
61 * This class implements the standard INCL4 avatar prediction and particle
62 * propagation logic. The main idea is to predict all collisions between particles
63 * and their reflections from the potential wall. After this we select the avatar
64 * with the smallest time, propagate all particles to their positions at that time
65 * and return the avatar to the INCL kernel @see G4INCL::Kernel.
66 *
67 * The particle trajectories in this propagation model are straight lines and all
68 * particles are assumed to move with constant velocity.
69 */
71 public:
72 StandardPropagationModel(LocalEnergyType localEnergyType, LocalEnergyType localEnergyDeltaType);
74
76 /**
77 * Set the nucleus for this propagation model.
78 */
79 void setNucleus(G4INCL::Nucleus *nucleus);
80
81 /**
82 * Get the nucleus.
83 */
85
86 G4double shoot(ParticleSpecies const projectileSpecies, const G4double kineticEnergy, const G4double impactParameter, const G4double phi);
87 G4double shootParticle(ParticleType const t, const G4double kineticEnergy, const G4double impactParameter, const G4double phi);
88 G4double shootComposite(ParticleSpecies const s, const G4double kineticEnergy, const G4double impactParameter, const G4double phi);
89
90 /**
91 * Set the stopping time of the simulation.
92 */
94
95 /**
96 * Get the current stopping time.
97 */
99
100 /**
101 * Add an avatar to the storage.
102 */
103 void registerAvatar(G4INCL::IAvatar *anAvatar);
104
105 /** \brief Generate a two-particle avatar.
106 *
107 * Generate a two-particle avatar, if all the appropriate conditions are
108 * met.
109 */
110 IAvatar *generateBinaryCollisionAvatar(Particle * const p1, Particle * const p2) const;
111
112 /** \brief Get the reflection time.
113 *
114 * Returns the reflection time of a particle on the potential wall.
115 *
116 * \param aParticle pointer to the particle
117 */
118 G4double getReflectionTime(G4INCL::Particle const * const aParticle);
119
120 /**
121 * Get the predicted time of the collision between two particles.
122 */
123 G4double getTime(G4INCL::Particle const * const particleA,
124 G4INCL::Particle const * const particleB, G4double *minDistOfApproach) const;
125
126 /** \brief Generate and register collisions between a list of updated particles and all the other particles.
127 *
128 * This method does not generate collisions among the particles in
129 * updatedParticles; in other words, it generates a collision between one
130 * of the updatedParticles and one of the particles ONLY IF the latter
131 * does not belong to updatedParticles.
132 *
133 * If you intend to generate all possible collisions among particles in a
134 * list, use generateCollisions().
135 *
136 * \param updatedParticles list of updated particles
137 * \param particles list of particles
138 */
139 void generateUpdatedCollisions(const ParticleList &updatedParticles, const ParticleList &particles);
140
141 /** \brief Generate and register collisions among particles in a list, except between those in another list.
142 *
143 * This method generates all possible collisions among the particles.
144 * Each collision is generated only once. The collision is NOT generated
145 * if BOTH collision partners belong to the except list.
146 *
147 * You should pass an empty list as the except parameter if you want to
148 * generate all possible collisions among particles.
149 *
150 * \param particles list of particles
151 * \param except list of excluded particles
152 */
153 void generateCollisions(const ParticleList &particles, const ParticleList &except);
154
155 /** \brief Generate decays for particles that can decay.
156 *
157 * The list of particles given as an argument is allowed to contain also
158 * stable particles.
159 *
160 * \param particles list of particles to (possibly) generate decays for
161 */
162 void generateDecays(const ParticleList &particles);
163
164 /**
165 * Update all avatars related to a particle.
166 */
167 void updateAvatars(const ParticleList &particles);
168
169 /** \brief (Re)Generate all possible avatars.
170 *
171 * \param excludeUpdated exclude collisions between updated particles.
172 */
173 void generateAllAvatars(G4bool excludeUpdated=false);
174
175 /**
176 * Propagate all particles and return the first avatar.
177 */
179
180 private:
181 G4INCL::Nucleus *theNucleus;
182 G4double maximumTime;
183 G4double currentTime;
184 G4bool firstAvatar;
185 LocalEnergyType theLocalEnergyType, theLocalEnergyDeltaType;
186
187 /// \brief Put spectators on shell by extracting energy from the participants.
188 void putSpectatorsOnShell(IAvatarList const &entryAvatars, ParticleList const &spectators);
189 };
190
191}
192
193
194#endif
double G4double
Definition: G4Types.hh:64
bool G4bool
Definition: G4Types.hh:67
G4double shootComposite(ParticleSpecies const s, const G4double kineticEnergy, const G4double impactParameter, const G4double phi)
G4double getReflectionTime(G4INCL::Particle const *const aParticle)
Get the reflection time.
G4double getTime(G4INCL::Particle const *const particleA, G4INCL::Particle const *const particleB, G4double *minDistOfApproach) const
void generateAllAvatars(G4bool excludeUpdated=false)
(Re)Generate all possible avatars.
void registerAvatar(G4INCL::IAvatar *anAvatar)
void updateAvatars(const ParticleList &particles)
IAvatar * generateBinaryCollisionAvatar(Particle *const p1, Particle *const p2) const
Generate a two-particle avatar.
G4double shoot(ParticleSpecies const projectileSpecies, const G4double kineticEnergy, const G4double impactParameter, const G4double phi)
void generateUpdatedCollisions(const ParticleList &updatedParticles, const ParticleList &particles)
Generate and register collisions between a list of updated particles and all the other particles.
void generateCollisions(const ParticleList &particles, const ParticleList &except)
Generate and register collisions among particles in a list, except between those in another list.
void generateDecays(const ParticleList &particles)
Generate decays for particles that can decay.
G4double shootParticle(ParticleType const t, const G4double kineticEnergy, const G4double impactParameter, const G4double phi)
std::list< IAvatar * > IAvatarList
std::list< G4INCL::Particle * > ParticleList