Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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// 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/*
39 * StandardPropagationModel.hh
40 *
41 * \date 4 June 2009
42 * \author Pekka Kaitaniemi
43 */
44
45#ifndef G4INCLStandardPropagationModel_hh
46#define G4INCLStandardPropagationModel_hh 1
47
48#include "G4INCLNucleus.hh"
50#include "G4INCLIAvatar.hh"
51#include "G4INCLConfigEnums.hh"
52
53#include <iterator>
54
55namespace G4INCL {
56
57 /**
58 * Standard INCL4 particle propagation and avatar prediction
59 *
60 * This class implements the standard INCL4 avatar prediction and particle
61 * propagation logic. The main idea is to predict all collisions between particles
62 * and their reflections from the potential wall. After this we select the avatar
63 * with the smallest time, propagate all particles to their positions at that time
64 * and return the avatar to the INCL kernel @see G4INCL::Kernel.
65 *
66 * The particle trajectories in this propagation model are straight lines and all
67 * particles are assumed to move with constant velocity.
68 */
70 public:
71 StandardPropagationModel(LocalEnergyType localEnergyType, LocalEnergyType localEnergyDeltaType, const G4double hTime = 0.0);
73
75 /**
76 * Set the nucleus for this propagation model.
77 */
78 void setNucleus(G4INCL::Nucleus *nucleus);
79
80 /**
81 * Get the nucleus.
82 */
84
85 G4double shoot(ParticleSpecies const &projectileSpecies, const G4double kineticEnergy, const G4double impactParameter, const G4double phi);
86 G4double shootParticle(ParticleType const t, const G4double kineticEnergy, const G4double impactParameter, const G4double phi);
87 G4double shootComposite(ParticleSpecies const &s, const G4double kineticEnergy, const G4double impactParameter, const G4double phi);
88 G4double shootAtrest(ParticleType const t, const G4double kineticEnergy);
89 G4double shootCompositeAtrest(ParticleSpecies const &s, const G4double kineticEnergy);
90
91 /**
92 * Set the stopping time of the simulation.
93 */
95
96 /**
97 * Get the current stopping time.
98 */
100
101 /**
102 * Add an avatar to the storage.
103 */
104 void registerAvatar(G4INCL::IAvatar *anAvatar);
105
106 /** \brief Generate a two-particle avatar.
107 *
108 * Generate a two-particle avatar, if all the appropriate conditions are
109 * met.
110 */
112
113 /** \brief Get the reflection time.
114 *
115 * Returns the reflection time of a particle on the potential wall.
116 *
117 * \param aParticle pointer to the particle
118 */
119 G4double getReflectionTime(G4INCL::Particle const * const aParticle);
120
121 /**
122 * Get the predicted time of the collision between two particles.
123 */
124 G4double getTime(G4INCL::Particle const * const particleA,
125 G4INCL::Particle const * const particleB, G4double *minDistOfApproach) const;
126
127 /** \brief Generate and register collisions between a list of updated particles and all the other particles.
128 *
129 * This method does not generate collisions among the particles in
130 * updatedParticles; in other words, it generates a collision between one
131 * of the updatedParticles and one of the particles ONLY IF the latter
132 * does not belong to updatedParticles.
133 *
134 * If you intend to generate all possible collisions among particles in a
135 * list, use generateCollisions().
136 *
137 * \param updatedParticles list of updated particles
138 * \param particles list of particles
139 */
140 void generateUpdatedCollisions(const ParticleList &updatedParticles, const ParticleList &particles);
141
142 /** \brief Generate and register collisions among particles in a list, except between those in another list.
143 *
144 * This method generates all possible collisions among the particles.
145 * Each collision is generated only once.
146 *
147 * \param particles list of particles
148 */
149 void generateCollisions(const ParticleList &particles);
150
151 /** \brief Generate and register collisions among particles in a list, except between those in another list.
152 *
153 * This method generates all possible collisions among the particles.
154 * Each collision is generated only once. The collision is NOT generated
155 * if BOTH collision partners belong to the except list.
156 *
157 * You should pass an empty list as the except parameter if you want to
158 * generate all possible collisions among particles.
159 *
160 * \param particles list of particles
161 * \param except list of excluded particles
162 */
163 void generateCollisions(const ParticleList &particles, const ParticleList &except);
164
165 /** \brief Generate decays for particles that can decay.
166 *
167 * The list of particles given as an argument is allowed to contain also
168 * stable particles.
169 *
170 * \param particles list of particles to (possibly) generate decays for
171 */
172 void generateDecays(const ParticleList &particles);
173
174 /**
175 * Update all avatars related to a particle.
176 */
177 void updateAvatars(const ParticleList &particles);
178
179 /// \brief (Re)Generate all possible avatars.
180 void generateAllAvatars();
181
182#ifdef INCL_REGENERATE_AVATARS
183 /** \brief (Re)Generate all possible avatars.
184 *
185 * This method excludes collision avatars between updated particles.
186 */
187 void generateAllAvatarsExceptUpdated(FinalState const * const fs);
188#endif
189
190 /**
191 * Propagate all particles and return the first avatar.
192 */
193 G4INCL::IAvatar* propagate(FinalState const * const fs);
194
195 private:
196 G4INCL::Nucleus *theNucleus;
197 G4double maximumTime;
198 G4double currentTime;
199 G4double hadronizationTime;
200 G4bool firstAvatar;
201 LocalEnergyType theLocalEnergyType, theLocalEnergyDeltaType;
202 Particle backupParticle1, backupParticle2;
203 };
204
205}
206
207
208#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
G4double getReflectionTime(G4INCL::Particle const *const aParticle)
Get the reflection time.
G4INCL::IAvatar * propagate(FinalState const *const fs)
G4double getTime(G4INCL::Particle const *const particleA, G4INCL::Particle const *const particleB, G4double *minDistOfApproach) const
void generateCollisions(const ParticleList &particles)
Generate and register collisions among particles in a list, except between those in another list.
G4double shootAtrest(ParticleType const t, const G4double kineticEnergy)
StandardPropagationModel(LocalEnergyType localEnergyType, LocalEnergyType localEnergyDeltaType, const G4double hTime=0.0)
void registerAvatar(G4INCL::IAvatar *anAvatar)
void updateAvatars(const ParticleList &particles)
void generateAllAvatars()
(Re)Generate all possible avatars.
G4double shoot(ParticleSpecies const &projectileSpecies, const G4double kineticEnergy, const G4double impactParameter, const G4double phi)
IAvatar * generateBinaryCollisionAvatar(Particle *const p1, Particle *const p2)
Generate a two-particle avatar.
G4double shootCompositeAtrest(ParticleSpecies const &s, const G4double kineticEnergy)
void generateUpdatedCollisions(const ParticleList &updatedParticles, const ParticleList &particles)
Generate and register collisions between a list of updated particles and all the other particles.
G4double shootComposite(ParticleSpecies const &s, const G4double kineticEnergy, const G4double impactParameter, const G4double phi)
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)