|
||||
File indexing completed on 2025-01-18 09:58:08
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 // 0027 // Contact authors: S. Meylan, C. Villagrasa 0028 // 0029 // email: sylvain.meylan@symalgo-tech.com, carmen.villagrasa@irsn.fr 0030 // updated : Hoang Tran : 6/1/2023 clean code 0031 0032 #ifndef G4DNAMODELINTERFACE_HH 0033 #define G4DNAMODELINTERFACE_HH 0034 0035 #include "G4VEmModel.hh" 0036 0037 #include <map> 0038 class G4ParticleChangeForGamma; 0039 class G4VDNAModel; 0040 class G4DNAModelInterface : public G4VEmModel 0041 { 0042 using MaterialParticleModelTable = 0043 std::map<std::size_t /*MatID*/, std::map<const G4ParticleDefinition*, G4VEmModel*>>; 0044 //should have only one model 0045 public: 0046 /*! 0047 * \brief G4DNAModelManager 0048 * Constructor 0049 * \param nam 0050 */ 0051 explicit G4DNAModelInterface(const G4String& nam); 0052 0053 /*! 0054 * \brief ~G4DNAModelManager 0055 * Destructor 0056 */ 0057 ~G4DNAModelInterface() override = default; 0058 0059 G4DNAModelInterface(const G4DNAModelInterface&) = delete; // prevent copy-construction 0060 G4DNAModelInterface& operator=(const G4DNAModelInterface& right) = delete; // prevent assignement 0061 0062 /*! 0063 * \brief Initialise 0064 * Initialise method to call all the initialise methods of the registered models 0065 * \param particle 0066 * \param cuts 0067 */ 0068 void Initialise(const G4ParticleDefinition* particle, const G4DataVector& cuts) override; 0069 0070 /*! 0071 * \brief CrossSectionPerVolume 0072 * Method called by the process and used to call the CrossSectionPerVolume method of the 0073 * registered models. The method also calculates through G4DNAMolecularMaterial the number of 0074 * molecule per volume unit for the current material or (component of a composite material). 0075 * \param material 0076 * \param p 0077 * \param ekin 0078 * \param emin 0079 * \param emax 0080 * \return the final cross section value times with the number of molecule per volume unit 0081 */ 0082 G4double CrossSectionPerVolume(const G4Material* material, const G4ParticleDefinition* p, 0083 G4double ekin, G4double emin, G4double emax) override; 0084 0085 /*! 0086 * \brief SampleSecondaries 0087 * Used to call the SampleSecondaries method of the registered models. A sampling is done to 0088 * select a component if the material is a composite one. \param fVect \param couple \param 0089 * aDynamicElectron \param tmin \param tmax 0090 */ 0091 void SampleSecondaries(std::vector<G4DynamicParticle*>* fVect, const G4MaterialCutsCouple* couple, 0092 const G4DynamicParticle* aDynamicElectron, G4double tmin, G4double tmax) override; 0093 0094 /*! 0095 * \brief RegisterModel 0096 * Method used to associate a model with the interaction 0097 * \param model 0098 */ 0099 void RegisterModel(G4VEmModel* model); 0100 /*! 0101 * \brief GetSelectedMaterial 0102 * To allow the user to retrieve the selected material in case of a composite material. 0103 * \return the last selected material by SampleSecondaries. 0104 */ 0105 inline std::size_t GetSelectedMaterial() 0106 { 0107 return fSampledMat; 0108 } 0109 0110 void StreamInfo(std::ostream& os) const; 0111 0112 private: 0113 /*! 0114 * \brief BuildMaterialParticleModelTable 0115 * Method used to build a map allowing the code to quickly retrieve the good model for a 0116 * particle/material couple \param p 0117 */ 0118 void BuildMaterialParticleModelTable(const G4ParticleDefinition* p); 0119 0120 void BuildMaterialMolPerVolTable(); 0121 0122 /*! 0123 * \brief InsertModelInTable 0124 * Used to put a model in the table after performing some checks. 0125 * \param matName 0126 * \param pName 0127 */ 0128 void InsertModelInTable(const std::size_t& matID, const G4ParticleDefinition* p); 0129 0130 /*! 0131 * \brief GetDNAModel 0132 * \param material 0133 * \param particle 0134 * \param ekin 0135 * \return G4VDNAModel* 0136 * Return the model corresponding to the material, particle and energy specified. 0137 * This method will check the energy range of the models to find to good one for the current ekin. 0138 */ 0139 G4VEmModel* SelectModel( 0140 const std::size_t& material, const G4ParticleDefinition* particle, const G4double& ekin); 0141 0142 G4double GetNumMoleculePerVolumeUnitForMaterial(const G4Material* mat); 0143 G4double GetNumMolPerVolUnitForComponentInComposite( 0144 const G4Material* component, const G4Material* composite); 0145 0146 const G4String fName; ///< name of the interaction 0147 G4ParticleChangeForGamma* fpParticleChangeForGamma = nullptr; 0148 std::vector<G4VEmModel*> fRegisteredModels; ///< vector containing all the registered models 0149 0150 std::map<std::size_t, G4double> fMaterialCS; ///< map used to share information between 0151 ///< CrossSectionPerVolume and SampleSecondaries 0152 G4double fCSsumTot = 0; ///< value which contains the sum of all the component cross sections in 0153 ///< case of a composite material 0154 std::size_t fSampledMat = 0; ///< for the user to retrieve selected material/component 0155 MaterialParticleModelTable fMaterialParticleModelTable; 0156 std::map<std::size_t, const std::vector<G4double>*> fMaterialMolPerVol; 0157 G4Material* fpG4_WATER = nullptr; 0158 }; 0159 0160 #endif // G4DNAMODELINTERFACE_HH
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |