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