Geant4 Cross Reference |
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 // Authors: S. Meylan and C. Villagrasa (IRSN, France) 27 // Models come from 28 // M. Bug et al, Rad. Phys and Chem. 130, 459-479 (2017) 29 // 30 31 #ifndef G4DNAPTBIONISATIONMODEL_h 32 #define G4DNAPTBIONISATIONMODEL_h 1 33 34 #include "G4DNACrossSectionDataSet.hh" 35 #include "G4DNAGenericIonsManager.hh" 36 #include "G4DNAPTBAugerModel.hh" 37 #include "G4DNAPTBIonisationStructure.hh" 38 #include "G4Electron.hh" 39 #include "G4LogLogInterpolation.hh" 40 #include "G4NistManager.hh" 41 #include "G4ParticleChangeForGamma.hh" 42 #include "G4ProductionCutsTable.hh" 43 #include "G4Proton.hh" 44 #include "G4VDNAModel.hh" 45 46 /*! 47 * \brief The G4DNAPTBIonisationModel class 48 * Implements the PTB ionisation model. 49 */ 50 class G4DNAPTBIonisationModel : public G4VDNAModel 51 { 52 public: 53 using TriDimensionMap = 54 std::map<std::size_t, std::map<const G4ParticleDefinition*, 55 std::map<G4double, std::map<G4double, std::map<G4double, G4double>>>>>; 56 using VecMap = std::map<std::size_t, 57 std::map<const G4ParticleDefinition*, std::map<G4double, std::vector<G4double>>>>; 58 using VecMapWithShell = 59 std::map<std::size_t, std::map<const G4ParticleDefinition*, 60 std::map<G4double, std::map<G4double, std::vector<G4double>>>>>; 61 /*! 62 * \brief G4DNAPTBIonisationModel 63 * Constructor 64 * \param applyToMaterial 65 * \param p 66 * \param nam 67 * \param isAuger 68 */ 69 explicit G4DNAPTBIonisationModel(const G4String& applyToMaterial = "all", 70 const G4ParticleDefinition* p = nullptr, const G4String& nam = "DNAPTBIonisationModel", 71 const G4bool isAuger = true); 72 73 /*! 74 * \brief ~G4DNAPTBIonisationModel 75 * Destructor 76 */ 77 ~G4DNAPTBIonisationModel() override = default; 78 79 // copy constructor and hide assignment operator 80 G4DNAPTBIonisationModel(const G4DNAPTBIonisationModel&) = delete; // prevent copy-construction 81 G4DNAPTBIonisationModel& operator=( 82 const G4DNAPTBIonisationModel& right) = delete; // prevent assignement 83 84 85 /*! 86 * \brief Initialise 87 * Method called once at the beginning of the simulation. It is used to setup the list of the 88 * materials managed by the model and the energy limits. All the materials are setup but only a 89 * part of them can be activated by the user through the constructor. 90 */ 91 void Initialise(const G4ParticleDefinition* particle, const G4DataVector& data) override; 92 93 /*! 94 * \brief CrossSectionPerVolume 95 * Mandatory for every model the CrossSectionPerVolume method is in charge of returning the 96 * cross section value corresponding to the material, particle and energy current values. 97 * \param material 98 * \param materialName 99 * \param p 100 * \param ekin 101 * \param emin 102 * \param emax 103 * \return the cross section value 104 */ 105 G4double CrossSectionPerVolume(const G4Material* material, const G4ParticleDefinition* p, 106 G4double ekin, G4double emin, G4double emax) override; 107 108 /*! 109 * \brief SampleSecondaries 110 * If the model is selected for the ModelInterface then SampleSecondaries will be called. 111 * The method sets the characteristics of the particles implied with the physical process after 112 * the ModelInterface (energy, momentum...). This method is mandatory for every model. \param 113 * materialName \param particleChangeForGamma \param tmin \param tmax 114 */ 115 void SampleSecondaries(std::vector<G4DynamicParticle*>*, const G4MaterialCutsCouple*, 116 const G4DynamicParticle*, G4double tmin, G4double tmax) override; 117 118 G4ParticleChangeForGamma* fParticleChangeForGamma = nullptr; 119 120 private: 121 std::unique_ptr<G4DNAPTBAugerModel> 122 fpDNAPTBAugerModel; ///< PTB Auger model instanciated in the constructor and deleted in the 123 ///< destructor of the class 124 G4int verboseLevel = 0; ///< verbose level 125 G4DNAPTBIonisationStructure 126 ptbStructure; /*!< ptbStructure class which contains the shell binding energies */ 127 TriDimensionMap diffCrossSectionData; 128 TriDimensionMap fEnergySecondaryData; 129 std::map<std::size_t, std::map<const G4ParticleDefinition*, std::vector<G4double>>> fTMapWithVec; 130 VecMap fEMapWithVector; 131 VecMapWithShell fProbaShellMap; 132 133 G4double RandomizeEjectedElectronEnergy(const G4ParticleDefinition* aP, 134 G4double incomingParticleEnergy, G4int shell, const std::size_t& materialName); 135 G4double DifferentialCrossSection(const G4ParticleDefinition* p, G4double k, 136 G4double energyTransfer, G4int shell, const std::size_t& materialName); 137 138 /*! 139 * \brief RandomizeEjectedElectronEnergyFromCumulated 140 * Uses the cumulated tables to find the energy of the ejected particle (electron) 141 * \param particleDefinition 142 * \param k 143 * \param shell 144 * \param materialName 145 * \return the ejected electron energy 146 */ 147 G4double RandomizeEjectedElectronEnergyFromCumulated( 148 const G4ParticleDefinition*, G4double k, G4int shell, const std::size_t& materialID); 149 150 /*! 151 * \brief RandomizeEjectedElectronDirection 152 * Method to calculate the ejected electron direction 153 * \param aParticleDefinition 154 * \param incomingParticleEnergy 155 * \param outgoingParticleEnergy 156 * \param cosTheta 157 * \param phi 158 */ 159 void RandomizeEjectedElectronDirection(const G4ParticleDefinition*, 160 G4double incomingParticleEnergy, G4double outgoingParticleEnergy, G4double& cosTheta, 161 G4double& phi); 162 /*! 163 * \brief ReadDiffCSFile 164 * Method to read the differential cross section files. 165 * \param materialName 166 * \param particleName 167 * \param file 168 * \param scaleFactor 169 */ 170 void ReadDiffCSFile(const std::size_t& materialName, const G4ParticleDefinition* p, 171 const G4String& file, const G4double& scaleFactor) override; 172 173 /*! 174 * \brief QuadInterpolator 175 * \param e11 176 * \param e12 177 * \param e21 178 * \param e22 179 * \param xs11 180 * \param xs12 181 * \param xs21 182 * \param xs22 183 * \param t1 184 * \param t2 185 * \param t 186 * \param e 187 * \return the interpolated value 188 */ 189 G4double QuadInterpolator(G4double e11, G4double e12, G4double e21, G4double e22, G4double xs11, 190 G4double xs12, G4double xs21, G4double xs22, G4double t1, G4double t2, G4double t, G4double e); 191 /*! 192 * \brief LogLogInterpolate 193 * \param e1 194 * \param e2 195 * \param e 196 * \param xs1 197 * \param xs2 198 * \return the interpolate value 199 */ 200 G4double LogLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2); 201 202 G4Material* fpGuanine_PU = nullptr; 203 G4Material* fpTHF = nullptr; 204 G4Material* fpPY = nullptr; 205 G4Material* fpPU = nullptr; 206 G4Material* fpTMP = nullptr; 207 G4Material* fpG4_WATER = nullptr; 208 G4Material* fpBackbone_THF = nullptr; 209 G4Material* fpCytosine_PY = nullptr; 210 G4Material* fpThymine_PY = nullptr; 211 G4Material* fpAdenine_PU = nullptr; 212 G4Material* fpBackbone_TMP = nullptr; 213 G4Material* fpN2 = nullptr; 214 G4DNAPTBIonisationModel* fpModelData = nullptr; 215 216 }; 217 218 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 219 220 #endif 221