|
||||
File indexing completed on 2025-01-18 09:58:29
0001 // 0002 // ******************************************************************** 0003 // * License and Disclaimer * 0004 // * * 0005 // * The Geant4 software is copyright of the Copyright Holders of * 0006 // * the Geant4 Collaboration. It is provided under the terms and * 0007 // * conditions of the Geant4 Software License, included in the file * 0008 // * LICENSE and available at http://cern.ch/geant4/license . These * 0009 // * include a list of copyright holders. * 0010 // * * 0011 // * Neither the authors of this software system, nor their employing * 0012 // * institutes,nor the agencies providing financial support for this * 0013 // * work make any representation or warranty, express or implied, * 0014 // * regarding this software system or assume any liability for its * 0015 // * use. Please see the license in the file LICENSE and URL above * 0016 // * for the full disclaimer and the limitation of liability. * 0017 // * * 0018 // * This code implementation is the result of the scientific and * 0019 // * technical work of the GEANT4 collaboration. * 0020 // * By using, copying, modifying or distributing the software (or * 0021 // * any work based on the software) you agree to acknowledge its * 0022 // * use in resulting scientific publications, and indicate your * 0023 // * acceptance of all terms of the Geant4 Software license. * 0024 // ******************************************************************** 0025 // 0026 // INCL++ intra-nuclear cascade model 0027 // Alain Boudard, CEA-Saclay, France 0028 // Joseph Cugnon, University of Liege, Belgium 0029 // Jean-Christophe David, CEA-Saclay, France 0030 // Pekka Kaitaniemi, CEA-Saclay, France, and Helsinki Institute of Physics, Finland 0031 // Sylvie Leray, CEA-Saclay, France 0032 // Davide Mancusi, CEA-Saclay, France 0033 // 0034 #define INCLXX_IN_GEANT4_MODE 1 0035 0036 #include "globals.hh" 0037 0038 /* 0039 * StandardPropagationModel.hh 0040 * 0041 * \date 4 June 2009 0042 * \author Pekka Kaitaniemi 0043 */ 0044 0045 #ifndef G4INCLStandardPropagationModel_hh 0046 #define G4INCLStandardPropagationModel_hh 1 0047 0048 #include "G4INCLNucleus.hh" 0049 #include "G4INCLIPropagationModel.hh" 0050 #include "G4INCLIAvatar.hh" 0051 #include "G4INCLConfigEnums.hh" 0052 0053 #include <iterator> 0054 0055 namespace G4INCL { 0056 0057 /** 0058 * Standard INCL4 particle propagation and avatar prediction 0059 * 0060 * This class implements the standard INCL4 avatar prediction and particle 0061 * propagation logic. The main idea is to predict all collisions between particles 0062 * and their reflections from the potential wall. After this we select the avatar 0063 * with the smallest time, propagate all particles to their positions at that time 0064 * and return the avatar to the INCL kernel @see G4INCL::Kernel. 0065 * 0066 * The particle trajectories in this propagation model are straight lines and all 0067 * particles are assumed to move with constant velocity. 0068 */ 0069 class StandardPropagationModel: public G4INCL::IPropagationModel { 0070 public: 0071 StandardPropagationModel(LocalEnergyType localEnergyType, LocalEnergyType localEnergyDeltaType, const G4double hTime = 0.0); 0072 virtual ~StandardPropagationModel(); 0073 0074 G4double getCurrentTime(); 0075 /** 0076 * Set the nucleus for this propagation model. 0077 */ 0078 void setNucleus(G4INCL::Nucleus *nucleus); 0079 0080 /** 0081 * Get the nucleus. 0082 */ 0083 G4INCL::Nucleus* getNucleus(); 0084 0085 G4double shoot(ParticleSpecies const &projectileSpecies, const G4double kineticEnergy, const G4double impactParameter, const G4double phi); 0086 G4double shootParticle(ParticleType const t, const G4double kineticEnergy, const G4double impactParameter, const G4double phi); 0087 G4double shootComposite(ParticleSpecies const &s, const G4double kineticEnergy, const G4double impactParameter, const G4double phi); 0088 G4double shootAtrest(ParticleType const t, const G4double kineticEnergy); 0089 0090 /** 0091 * Set the stopping time of the simulation. 0092 */ 0093 void setStoppingTime(G4double); 0094 0095 /** 0096 * Get the current stopping time. 0097 */ 0098 G4double getStoppingTime(); 0099 0100 /** 0101 * Add an avatar to the storage. 0102 */ 0103 void registerAvatar(G4INCL::IAvatar *anAvatar); 0104 0105 /** \brief Generate a two-particle avatar. 0106 * 0107 * Generate a two-particle avatar, if all the appropriate conditions are 0108 * met. 0109 */ 0110 IAvatar *generateBinaryCollisionAvatar(Particle * const p1, Particle * const p2); 0111 0112 /** \brief Get the reflection time. 0113 * 0114 * Returns the reflection time of a particle on the potential wall. 0115 * 0116 * \param aParticle pointer to the particle 0117 */ 0118 G4double getReflectionTime(G4INCL::Particle const * const aParticle); 0119 0120 /** 0121 * Get the predicted time of the collision between two particles. 0122 */ 0123 G4double getTime(G4INCL::Particle const * const particleA, 0124 G4INCL::Particle const * const particleB, G4double *minDistOfApproach) const; 0125 0126 /** \brief Generate and register collisions between a list of updated particles and all the other particles. 0127 * 0128 * This method does not generate collisions among the particles in 0129 * updatedParticles; in other words, it generates a collision between one 0130 * of the updatedParticles and one of the particles ONLY IF the latter 0131 * does not belong to updatedParticles. 0132 * 0133 * If you intend to generate all possible collisions among particles in a 0134 * list, use generateCollisions(). 0135 * 0136 * \param updatedParticles list of updated particles 0137 * \param particles list of particles 0138 */ 0139 void generateUpdatedCollisions(const ParticleList &updatedParticles, const ParticleList &particles); 0140 0141 /** \brief Generate and register collisions among particles in a list, except between those in another list. 0142 * 0143 * This method generates all possible collisions among the particles. 0144 * Each collision is generated only once. 0145 * 0146 * \param particles list of particles 0147 */ 0148 void generateCollisions(const ParticleList &particles); 0149 0150 /** \brief Generate and register collisions among particles in a list, except between those in another list. 0151 * 0152 * This method generates all possible collisions among the particles. 0153 * Each collision is generated only once. The collision is NOT generated 0154 * if BOTH collision partners belong to the except list. 0155 * 0156 * You should pass an empty list as the except parameter if you want to 0157 * generate all possible collisions among particles. 0158 * 0159 * \param particles list of particles 0160 * \param except list of excluded particles 0161 */ 0162 void generateCollisions(const ParticleList &particles, const ParticleList &except); 0163 0164 /** \brief Generate decays for particles that can decay. 0165 * 0166 * The list of particles given as an argument is allowed to contain also 0167 * stable particles. 0168 * 0169 * \param particles list of particles to (possibly) generate decays for 0170 */ 0171 void generateDecays(const ParticleList &particles); 0172 0173 /** 0174 * Update all avatars related to a particle. 0175 */ 0176 void updateAvatars(const ParticleList &particles); 0177 0178 /// \brief (Re)Generate all possible avatars. 0179 void generateAllAvatars(); 0180 0181 #ifdef INCL_REGENERATE_AVATARS 0182 /** \brief (Re)Generate all possible avatars. 0183 * 0184 * This method excludes collision avatars between updated particles. 0185 */ 0186 void generateAllAvatarsExceptUpdated(FinalState const * const fs); 0187 #endif 0188 0189 /** 0190 * Propagate all particles and return the first avatar. 0191 */ 0192 G4INCL::IAvatar* propagate(FinalState const * const fs); 0193 0194 private: 0195 G4INCL::Nucleus *theNucleus; 0196 G4double maximumTime; 0197 G4double currentTime; 0198 G4double hadronizationTime; 0199 G4bool firstAvatar; 0200 LocalEnergyType theLocalEnergyType, theLocalEnergyDeltaType; 0201 Particle backupParticle1, backupParticle2; 0202 }; 0203 0204 } 0205 0206 0207 #endif
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |