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 // ------------------------------------------------------------------- 28 // 29 // GEANT4 Class header file 30 // 31 // 32 // File name: G4EmCorrections 33 // 34 // Author: Vladimir Ivanchenko 35 // 36 // Creation date: 13.01.2005 37 // 38 // Modifications: 39 // 28.04.2006 General cleanup, add finite size corrections (V.Ivanchenko) 40 // 13.05.2006 Add corrections for ion stopping (V.Ivanhcenko) 41 // 20.05.2008 Removed Finite Size correction (V.Ivanchenko) 42 // 12.09.2008 Added inlined interfaces to effective charge (V.Ivanchenko) 43 // 19.04.2012 Fix reproducibility problem (A.Ribon) 44 // 45 // Class Description: 46 // 47 // This class provides calculation of EM corrections to ionisation 48 // 49 50 // ------------------------------------------------------------------- 51 // 52 53 #ifndef G4EmCorrections_h 54 #define G4EmCorrections_h 1 55 56 #include "globals.hh" 57 #include "G4ionEffectiveCharge.hh" 58 #include "G4Material.hh" 59 #include "G4ParticleDefinition.hh" 60 61 class G4VEmModel; 62 class G4PhysicsVector; 63 class G4IonTable; 64 class G4MaterialCutsCouple; 65 class G4PhysicsFreeVector; 66 class G4Pow; 67 68 class G4EmCorrections 69 { 70 71 public: 72 73 explicit G4EmCorrections(G4int verb); 74 75 ~G4EmCorrections(); 76 77 G4double HighOrderCorrections(const G4ParticleDefinition*, 78 const G4Material*, 79 const G4double kineticEnergy, 80 const G4double cutEnergy); 81 82 G4double IonHighOrderCorrections(const G4ParticleDefinition*, 83 const G4MaterialCutsCouple*, 84 const G4double kineticEnergy); 85 86 G4double ComputeIonCorrections(const G4ParticleDefinition*, 87 const G4Material*, 88 const G4double kineticEnergy); 89 90 G4double IonBarkasCorrection(const G4ParticleDefinition*, 91 const G4Material*, 92 const G4double kineticEnergy); 93 94 G4double Bethe(const G4ParticleDefinition*, 95 const G4Material*, 96 const G4double kineticEnergy); 97 98 G4double SpinCorrection(const G4ParticleDefinition*, 99 const G4Material*, 100 const G4double kineticEnergy); 101 102 G4double KShellCorrection(const G4ParticleDefinition*, 103 const G4Material*, 104 const G4double kineticEnergy); 105 106 G4double LShellCorrection(const G4ParticleDefinition*, 107 const G4Material*, 108 const G4double kineticEnergy); 109 110 G4double ShellCorrection(const G4ParticleDefinition*, 111 const G4Material*, 112 const G4double kineticEnergy); 113 114 G4double ShellCorrectionSTD(const G4ParticleDefinition*, 115 const G4Material*, 116 const G4double kineticEnergy); 117 118 G4double DensityCorrection(const G4ParticleDefinition*, 119 const G4Material*, 120 const G4double kineticEnergy); 121 122 G4double BarkasCorrection(const G4ParticleDefinition*, 123 const G4Material*, 124 const G4double kineticEnergy, 125 const G4bool isInitialized = false); 126 127 G4double BlochCorrection(const G4ParticleDefinition*, 128 const G4Material*, 129 const G4double kineticEnergy, 130 const G4bool isInitialized = false); 131 132 G4double MottCorrection(const G4ParticleDefinition*, 133 const G4Material*, 134 const G4double kineticEnergy, 135 const G4bool isInitialized = false); 136 137 void AddStoppingData(const G4int Z, const G4int A, 138 const G4String& materialName, 139 G4PhysicsVector* dVector); 140 141 void InitialiseForNewRun(); 142 143 // effective charge correction using stopping power data 144 G4double EffectiveChargeCorrection(const G4ParticleDefinition*, 145 const G4Material*, 146 const G4double kineticEnergy); 147 148 // effective charge of an ion 149 inline G4double GetParticleCharge(const G4ParticleDefinition*, 150 const G4Material*, 151 const G4double kineticEnergy); 152 153 inline 154 G4double EffectiveChargeSquareRatio(const G4ParticleDefinition*, 155 const G4Material*, 156 const G4double kineticEnergy); 157 158 // ionisation models for ions 159 inline void SetIonisationModels(G4VEmModel* mod1 = nullptr, 160 G4VEmModel* mod2 = nullptr); 161 162 inline G4int GetNumberOfStoppingVectors() const; 163 164 inline void SetVerbose(G4int verb); 165 166 // hide assignment operator 167 G4EmCorrections & operator=(const G4EmCorrections &right) = delete; 168 G4EmCorrections(const G4EmCorrections&) = delete; 169 170 private: 171 172 void Initialise(); 173 174 void BuildCorrectionVector(); 175 176 void SetupKinematics(const G4ParticleDefinition*, 177 const G4Material*, 178 const G4double kineticEnergy); 179 180 G4double KShell(const G4double theta, const G4double eta); 181 182 G4double LShell(const G4double theta, const G4double eta); 183 184 G4int Index(const G4double x, const G4double* y, const G4int n) const; 185 186 G4double Value(const G4double xv, const G4double x1, const G4double x2, 187 const G4double y1, const G4double y2) const; 188 189 G4double Value2(const G4double xv, const G4double yv, 190 const G4double x1, const G4double x2, 191 const G4double y1, const G4double y2, 192 const G4double z11, const G4double z21, 193 const G4double z12, const G4double z22) const; 194 195 G4Pow* g4calc; 196 G4IonTable* ionTable; 197 198 const G4ParticleDefinition* particle = nullptr; 199 const G4ParticleDefinition* curParticle = nullptr; 200 const G4Material* material = nullptr; 201 const G4Material* curMaterial = nullptr; 202 const G4ElementVector* theElementVector = nullptr; 203 const G4double* atomDensity = nullptr; 204 205 G4PhysicsVector* curVector = nullptr; 206 207 G4VEmModel* ionLEModel = nullptr; 208 G4VEmModel* ionHEModel = nullptr; 209 210 G4double kinEnergy = 0.0; 211 G4double mass = 0.0; 212 G4double massFactor = 1.0; 213 G4double tau = 0.0; 214 G4double gamma = 1.0; 215 G4double bg2 = 0.0; 216 G4double beta2 = 0.0; 217 G4double beta = 0.0; 218 G4double ba2 = 0.0; 219 G4double tmax = 0.0; 220 G4double charge = 0.0; 221 G4double q2 = 0.0; 222 G4double eth; 223 G4double eCorrMin; 224 G4double eCorrMax; 225 226 std::size_t ncouples = 0; 227 std::size_t idxBarkas = 0; 228 G4int nK = 20; 229 G4int nL = 26; 230 G4int nEtaK = 29; 231 G4int nEtaL = 28; 232 G4int nbinCorr = 52; 233 G4int numberOfElements = 0; 234 235 // Ion stopping data 236 G4int nIons = 0; 237 G4int idx = 0; 238 G4int currentZ = 0; 239 240 G4int verbose; 241 G4bool isInitializer = false; 242 243 std::vector<G4int> Zion; 244 std::vector<G4int> Aion; 245 std::vector<G4String> materialName; 246 std::vector<const G4ParticleDefinition*> ionList; 247 248 std::map< G4int, std::vector<G4double> > thcorr; 249 250 std::vector<const G4Material*> currmat; 251 std::vector<const G4Material*> materialList; 252 std::vector<G4PhysicsVector*> stopData; 253 254 G4ionEffectiveCharge effCharge; 255 256 static const G4double ZD[11]; 257 static const G4double UK[20]; 258 static const G4double VK[20]; 259 static G4double ZK[20]; 260 static const G4double Eta[29]; 261 static G4double CK[20][29]; 262 static G4double CL[26][28]; 263 static const G4double UL[26]; 264 static G4double VL[26]; 265 266 static G4double sWmaxBarkas; 267 static G4PhysicsFreeVector* sBarkasCorr; 268 static G4PhysicsFreeVector* sThetaK; 269 static G4PhysicsFreeVector* sThetaL; 270 }; 271 272 inline G4int 273 G4EmCorrections::Index(const G4double x, const G4double* y, const G4int n) const 274 { 275 G4int iddd = n-1; 276 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko 277 do {--iddd;} while (iddd>0 && x<y[iddd]); 278 return iddd; 279 } 280 281 inline G4double G4EmCorrections::Value(const G4double xv, const G4double x1, 282 const G4double x2, 283 const G4double y1, const G4double y2) const 284 { 285 return y1 + (y2 - y1)*(xv - x1)/(x2 - x1); 286 } 287 288 inline G4double G4EmCorrections::Value2(const G4double xv, const G4double yv, 289 const G4double x1, const G4double x2, 290 const G4double y1, const G4double y2, 291 const G4double z11, const G4double z21, 292 const G4double z12, const G4double z22) const 293 { 294 return ( z11*(x2-xv)*(y2-yv) + z22*(xv-x1)*(yv-y1) + 295 z12*(x2-xv)*(yv-y1) + z21*(xv-x1)*(y2-yv) ) 296 / ((x2-x1)*(y2-y1)); 297 } 298 299 inline void 300 G4EmCorrections::SetIonisationModels(G4VEmModel* mod1, G4VEmModel* mod2) 301 { 302 if(nullptr != mod1) { ionLEModel = mod1; } 303 if(nullptr != mod2) { ionHEModel = mod2; } 304 } 305 306 inline G4int G4EmCorrections::GetNumberOfStoppingVectors() const 307 { 308 return nIons; 309 } 310 311 inline G4double 312 G4EmCorrections::GetParticleCharge(const G4ParticleDefinition* p, 313 const G4Material* mat, 314 const G4double kineticEnergy) 315 { 316 return effCharge.EffectiveCharge(p,mat,kineticEnergy); 317 } 318 319 inline G4double 320 G4EmCorrections::EffectiveChargeSquareRatio(const G4ParticleDefinition* p, 321 const G4Material* mat, 322 const G4double kineticEnergy) 323 { 324 return effCharge.EffectiveChargeSquareRatio(p,mat,kineticEnergy); 325 } 326 327 inline void G4EmCorrections::SetVerbose(G4int verb) 328 { 329 verbose = verb; 330 } 331 332 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 333 334 #endif 335