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 // CPA100 ionisation model class for electrons 27 // 28 // Based on the work of M. Terrissol and M. C. Bordage 29 // 30 // Users are requested to cite the following papers: 31 // - M. Terrissol, A. Baudre, Radiat. Prot. Dosim. 31 (1990) 175-177 32 // - M.C. Bordage, J. Bordes, S. Edel, M. Terrissol, X. Franceries, 33 // M. Bardies, N. Lampe, S. Incerti, Phys. Med. 32 (2016) 1833-1840 34 // 35 // Authors of this class: 36 // M.C. Bordage, M. Terrissol, S. Edel, J. Bordes, S. Incerti 37 // 38 // 15.01.2014: creation 39 // 40 // Based on the study by S. Zein et. al. Nucl. Inst. Meth. B 488 (2021) 70-82 41 // 1/2/2023 : Hoang added modification for DNA cross sections 42 43 #ifndef G4DNACPA100IonisationModel_h 44 #define G4DNACPA100IonisationModel_h 1 45 46 #include "G4DNACPA100IonisationStructure.hh" 47 #include "G4DNACrossSectionDataSet.hh" 48 #include "G4Electron.hh" 49 #include "G4LogLogInterpolation.hh" 50 #include "G4NistManager.hh" 51 #include "G4ParticleChangeForGamma.hh" 52 #include "G4ProductionCutsTable.hh" 53 #include "G4VAtomDeexcitation.hh" 54 #include "G4VDNAModel.hh" 55 56 class G4DNACPA100IonisationModel : public G4VDNAModel 57 { 58 using TriDimensionMap = 59 std::map<std::size_t, 60 std::map<const G4ParticleDefinition*, 61 std::map<G4int, std::map<G4double, std::map<G4double, G4double>>>>>; 62 using VecMap = 63 std::map<std::size_t, 64 std::map<const G4ParticleDefinition*, std::map<G4double, std::vector<G4double>>>>; 65 using VecMapWithShell = 66 std::map<std::size_t, 67 std::map<const G4ParticleDefinition*, 68 std::map<G4double, std::map<G4double, std::vector<G4double>>>>>; 69 using PartKineticInMat = 70 const std::tuple<std::size_t /*Mat*/, G4double /*Energy*/, G4int /*shell*/>&; 71 72 public: 73 explicit G4DNACPA100IonisationModel(const G4ParticleDefinition* p = nullptr, 74 const G4String& nam = "DNACPA100IonisationModel"); 75 76 ~G4DNACPA100IonisationModel() override = default; 77 78 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override; 79 80 G4double CrossSectionPerVolume(const G4Material* material, const G4ParticleDefinition* p, 81 G4double ekin, G4double emin, G4double emax) override; 82 83 void SampleSecondaries(std::vector<G4DynamicParticle*>*, const G4MaterialCutsCouple*, 84 const G4DynamicParticle*, G4double tmin, G4double maxEnergy) override; 85 86 G4double DifferentialCrossSection(PartKineticInMat info, const G4double& energyTransfer); 87 88 // G4double DifferentialCrossSection(const G4double& k, 89 // const G4double& energyTransfer, const G4int& 90 // ionizationLevelIndex, const std::size_t& materialID); 91 92 inline void SelectFasterComputation(G4bool input); 93 94 inline void SelectUseDcs(G4bool input); 95 96 inline void SelectStationary(G4bool input); 97 98 G4DNACPA100IonisationModel& operator=(const G4DNACPA100IonisationModel& right) = delete; 99 G4DNACPA100IonisationModel(const G4DNACPA100IonisationModel&) = delete; 100 void ReadDiffCSFile(const std::size_t& materialID, const G4ParticleDefinition* p, 101 const G4String& file, const G4double& scaleFactor) override; 102 103 protected: 104 G4ParticleChangeForGamma* fParticleChangeForGamma = nullptr; 105 106 private: 107 G4bool statCode = false; 108 G4bool fasterCode = true; 109 G4bool useDcs = false; 110 111 // const std::vector<G4double>* fpMolMaterialDensity; 112 113 // Deexcitation manager to produce fluo photons and e- 114 G4VAtomDeexcitation* fAtomDeexcitation = nullptr; 115 116 G4bool isInitialised = false; 117 G4int verboseLevel = 0; 118 119 G4DNACPA100IonisationStructure iStructure; 120 121 G4double RandomizeEjectedElectronEnergy(PartKineticInMat info); 122 123 G4double RandomizeEjectedElectronEnergyFromCumulatedDcs(PartKineticInMat info); 124 125 G4double RandomizeEjectedElectronEnergyFromanalytical(PartKineticInMat info); 126 127 G4double RandomTransferedEnergy(PartKineticInMat info); 128 129 void RandomizeEjectedElectronDirection(G4ParticleDefinition* aParticleDefinition, 130 G4double incomingParticleEnergy, 131 G4double outgoingParticleEnergy, G4double& cosTheta, 132 G4double& phi); 133 134 G4double Interpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2); 135 136 G4double QuadInterpolator(G4double e11, G4double e12, G4double e21, G4double e22, G4double x11, 137 G4double x12, G4double x21, G4double x22, G4double t1, G4double t2, 138 G4double t, G4double e); 139 140 TriDimensionMap diffCrossSectionData, fEnergySecondaryData; 141 std::map<std::size_t, std::map<const G4ParticleDefinition*, std::vector<G4double>>> 142 fTMapWithVec; 143 VecMap fEMapWithVector; 144 VecMapWithShell fProbaShellMap; 145 146 const G4Material* fpGuanine = nullptr; 147 const G4Material* fpG4_WATER = nullptr; 148 const G4Material* fpDeoxyribose = nullptr; 149 const G4Material* fpCytosine = nullptr; 150 const G4Material* fpThymine = nullptr; 151 const G4Material* fpAdenine = nullptr; 152 const G4Material* fpPhosphate = nullptr; 153 const G4ParticleDefinition* fpParticle = nullptr; 154 G4DNACPA100IonisationModel* fpModelData = nullptr; 155 }; 156 157 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 158 159 inline void G4DNACPA100IonisationModel::SelectFasterComputation(G4bool input) 160 { 161 fasterCode = input; 162 } 163 164 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 165 166 inline void G4DNACPA100IonisationModel::SelectUseDcs(G4bool input) 167 { 168 useDcs = input; 169 } 170 171 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 172 173 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 174 175 inline void G4DNACPA100IonisationModel::SelectStationary(G4bool input) 176 { 177 statCode = input; 178 } 179 180 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 181 182 #endif 183