Geant4  10.01.p02
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 
55 namespace G4INCL {
56 
70  public:
71  StandardPropagationModel(LocalEnergyType localEnergyType, LocalEnergyType localEnergyDeltaType);
72  virtual ~StandardPropagationModel();
73 
78  void setNucleus(G4INCL::Nucleus *nucleus);
79 
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 
93 
98 
102  void registerAvatar(G4INCL::IAvatar *anAvatar);
103 
110 
117  G4double getReflectionTime(G4INCL::Particle const * const aParticle);
118 
122  G4double getTime(G4INCL::Particle const * const particleA,
123  G4INCL::Particle const * const particleB, G4double *minDistOfApproach) const;
124 
138  void generateUpdatedCollisions(const ParticleList &updatedParticles, const ParticleList &particles);
139 
147  void generateCollisions(const ParticleList &particles);
148 
161  void generateCollisions(const ParticleList &particles, const ParticleList &except);
162 
170  void generateDecays(const ParticleList &particles);
171 
175  void updateAvatars(const ParticleList &particles);
176 
178  void generateAllAvatars();
179 
180 #ifdef INCL_REGENERATE_AVATARS
181 
185  void generateAllAvatarsExceptUpdated(FinalState const * const fs);
186 #endif
187 
191  G4INCL::IAvatar* propagate(FinalState const * const fs);
192 
193  private:
200  };
201 
202 }
203 
204 
205 #endif
G4double shootComposite(ParticleSpecies const &s, const G4double kineticEnergy, const G4double impactParameter, const G4double phi)
G4double shoot(ParticleSpecies const &projectileSpecies, const G4double kineticEnergy, const G4double impactParameter, const G4double phi)
void registerAvatar(G4INCL::IAvatar *anAvatar)
Add an avatar to the storage.
void generateAllAvatars()
(Re)Generate all possible avatars.
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...
G4INCL::IAvatar * propagate(FinalState const *const fs)
Propagate all particles and return the first avatar.
G4double getReflectionTime(G4INCL::Particle const *const aParticle)
Get the reflection time.
void updateAvatars(const ParticleList &particles)
Update all avatars related to a particle.
void setNucleus(G4INCL::Nucleus *nucleus)
Set the nucleus for this propagation model.
void setStoppingTime(G4double)
Set the stopping time of the simulation.
Propagation model takes care of transporting the particles until something interesting (i...
G4double getTime(G4INCL::Particle const *const particleA, G4INCL::Particle const *const particleB, G4double *minDistOfApproach) const
Get the predicted time of the collision between two particles.
static const double s
Definition: G4SIunits.hh:150
Final state of an interaction.
bool G4bool
Definition: G4Types.hh:79
G4double getCurrentTime()
Returns the current global time of the system.
const G4double p2
const G4double p1
IAvatar * generateBinaryCollisionAvatar(Particle *const p1, Particle *const p2)
Generate a two-particle avatar.
Standard INCL4 particle propagation and avatar prediction.
void generateDecays(const ParticleList &particles)
Generate decays for particles that can decay.
G4INCL::Nucleus * getNucleus()
Get the nucleus.
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...
G4double getStoppingTime()
Get the current stopping time.