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 // 27 // GEANT4 Class header file 28 // 29 // File name: G4AtimaEnergyLossModel 30 // 31 // Author: Jose Luis Rodriguez Sanchez on base of ATIMA code 32 // 33 // Creation date: 16.01.2018 34 // 35 // Modifications: 36 // 37 // 38 // Class Description: 39 // 40 // Implementation of ATIMA model of energy loss 41 // by heavy charged particles 42 43 // ------------------------------------------------------------------- 44 // 45 46 #ifndef G4AtimaEnergyLossModel_h 47 #define G4AtimaEnergyLossModel_h 1 48 49 #include <CLHEP/Units/SystemOfUnits.h> 50 51 #include "G4VEmModel.hh" 52 #include "G4NistManager.hh" 53 54 class G4EmCorrections; 55 class G4ParticleChangeForLoss; 56 57 class G4AtimaEnergyLossModel : public G4VEmModel 58 { 59 60 public: 61 62 explicit G4AtimaEnergyLossModel(const G4ParticleDefinition* p = nullptr, 63 const G4String& nam = "Atima"); 64 65 ~G4AtimaEnergyLossModel() override; 66 67 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override; 68 69 G4double MinEnergyCut(const G4ParticleDefinition*, 70 const G4MaterialCutsCouple* couple) override; 71 72 G4double ComputeCrossSectionPerElectron( 73 const G4ParticleDefinition*, 74 G4double kineticEnergy, 75 G4double cutEnergy, 76 G4double maxEnergy); 77 78 G4double ComputeCrossSectionPerAtom( 79 const G4ParticleDefinition*, 80 G4double kineticEnergy, 81 G4double Z, G4double A, 82 G4double cutEnergy, 83 G4double maxEnergy) override; 84 85 G4double CrossSectionPerVolume(const G4Material*, 86 const G4ParticleDefinition*, 87 G4double kineticEnergy, 88 G4double cutEnergy, 89 G4double maxEnergy) override; 90 91 G4double ComputeDEDXPerVolume(const G4Material*, 92 const G4ParticleDefinition*, 93 G4double kineticEnergy, 94 G4double) override; 95 96 G4double GetChargeSquareRatio(const G4ParticleDefinition* p, 97 const G4Material* mat, 98 G4double kineticEnergy) override; 99 100 G4double GetParticleCharge(const G4ParticleDefinition* p, 101 const G4Material* mat, 102 G4double kineticEnergy) override; 103 104 void CorrectionsAlongStep(const G4MaterialCutsCouple*, 105 const G4DynamicParticle*, 106 const G4double&, 107 G4double&) override; 108 109 void SampleSecondaries(std::vector<G4DynamicParticle*>*, 110 const G4MaterialCutsCouple*, 111 const G4DynamicParticle*, 112 G4double tmin, 113 G4double maxEnergy) override; 114 115 protected: 116 117 G4double MaxSecondaryEnergy(const G4ParticleDefinition*, 118 G4double kinEnergy) override; 119 120 inline G4double GetChargeSquareRatio() const; 121 122 inline void SetChargeSquareRatio(G4double val); 123 124 private: 125 126 void SetupParameters(); 127 128 inline void SetParticle(const G4ParticleDefinition* p); 129 130 inline void SetGenericIon(const G4ParticleDefinition* p); 131 132 G4double StoppingPower(G4double ap, G4double zp, G4double ep, G4double at, G4double zt); 133 G4double Bethek_dedx_e(G4double ap,G4double zp,G4double ep,G4double at,G4double zt); 134 G4double dedx_n(const G4double ap, const G4double zp, const G4double ep, const G4double at, const G4double zt); 135 G4double sezi_dedx_e(const G4double zp, const G4double ep, const G4double at, const G4double zt); 136 G4double sezi_p_se(const G4double energy, const G4double at, const G4double zt); 137 G4double EnergyTable_interpolate(G4double xval, const G4double* y); 138 139 // hide assignment operator 140 G4AtimaEnergyLossModel & operator=(const G4AtimaEnergyLossModel &right) = delete; 141 G4AtimaEnergyLossModel(const G4AtimaEnergyLossModel&) = delete; 142 143 const G4ParticleDefinition* particle; 144 const G4ParticleDefinition* theElectron; 145 G4EmCorrections* corr; 146 G4ParticleChangeForLoss* fParticleChange; 147 G4NistManager* nist; 148 G4Pow* g4calc; 149 150 G4double mass = 0.0; 151 G4double tlimit = DBL_MAX; 152 G4double spin = 0.0; 153 G4double magMoment2 = 0.0; 154 G4double chargeSquare = 1.0; 155 G4double ratio = 1.0; 156 G4double formfact = 0.0; 157 G4double corrFactor = 1.0; 158 G4double MLN10; 159 G4double atomic_mass_unit; 160 G4double dedx_constant; 161 G4double electron_mass; 162 G4double fine_structure; 163 G4double domega2dx_constant; 164 165 G4bool isIon = false; 166 167 static G4double stepE; 168 static G4double tableE[200]; 169 static const G4double element_atomic_weights[110]; 170 static const G4double ls_coefficients_a[110][200]; 171 static const G4double ls_coefficients_ahi[110][200]; 172 static const G4double proton_stopping_coef[92][8]; 173 static const G4double ionisation_potentials_z[121]; 174 175 static const G4double atima_vfermi[92]; 176 static const G4double atima_lambda_screening[92]; 177 static const G4double x0[92]; 178 static const G4double x1[92]; 179 static const G4double afermi[92]; 180 static const G4double c[92]; 181 static const G4double m0[92]; 182 static const G4double del_0[92]; 183 184 }; 185 186 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 187 188 inline void G4AtimaEnergyLossModel::SetParticle(const G4ParticleDefinition* p) 189 { 190 if(particle != p) { 191 particle = p; 192 if(p->GetBaryonNumber() > 3 || p->GetPDGCharge() > CLHEP::eplus) 193 { isIon = true; } 194 SetupParameters(); 195 } 196 } 197 198 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 199 200 inline void G4AtimaEnergyLossModel::SetGenericIon(const G4ParticleDefinition* p) 201 { 202 if(p && p->GetParticleName() == "GenericIon") { isIon = true; } 203 } 204 205 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 206 207 inline G4double G4AtimaEnergyLossModel::GetChargeSquareRatio() const 208 { 209 return chargeSquare; 210 } 211 212 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 213 214 inline void G4AtimaEnergyLossModel::SetChargeSquareRatio(G4double val) 215 { 216 chargeSquare = val; 217 } 218 219 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 220 221 #endif 222