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 // ClassName: G4EnergyLossForExtrapolator 30 // 31 // Description: This class provide calculation of energy loss, fluctuation, 32 // and msc angle 33 // 34 // Author: 09.12.04 V.Ivanchenko 35 // 36 // Modification: 37 // 08-04-05 Rename Propogator -> Extrapolator 38 // 16-03-06 Add muon tables 39 // 21-03-06 Add verbosity defined in the constructor and Initialisation 40 // start only when first public method is called (V.Ivanchenko) 41 // 03-05-06 Remove unused pointer G4Material* from number of methods (VI) 42 // 28-07-07 Add maxEnergyTransfer for computation of energy loss (VI) 43 // 44 //---------------------------------------------------------------------------- 45 // 46 47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 48 49 #ifndef G4EnergyLossForExtrapolator_h 50 #define G4EnergyLossForExtrapolator_h 1 51 52 #include <vector> 53 #include <CLHEP/Units/PhysicalConstants.h> 54 55 #include "globals.hh" 56 #include "G4PhysicsTable.hh" 57 #include "G4TablesForExtrapolator.hh" 58 #include "G4Log.hh" 59 #include "G4Threading.hh" 60 61 class G4ParticleDefinition; 62 class G4Material; 63 class G4MaterialCutsCouple; 64 65 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 66 67 class G4EnergyLossForExtrapolator 68 { 69 public: 70 71 explicit G4EnergyLossForExtrapolator(G4int verb = 1); 72 73 ~G4EnergyLossForExtrapolator(); 74 75 void Initialisation(); 76 77 G4double ComputeDEDX(G4double kinEnergy, const G4ParticleDefinition*, 78 const G4Material*); 79 80 G4double ComputeRange(G4double kinEnergy, const G4ParticleDefinition*, 81 const G4Material*); 82 83 G4double ComputeEnergy(G4double range, const G4ParticleDefinition*, 84 const G4Material*); 85 86 G4double EnergyAfterStep(G4double kinEnergy, G4double step, 87 const G4Material*, const G4ParticleDefinition*); 88 89 G4double EnergyBeforeStep(G4double kinEnergy, G4double step, 90 const G4Material*, const G4ParticleDefinition*); 91 92 G4double TrueStepLength(G4double kinEnergy, G4double step, 93 const G4Material*, const G4ParticleDefinition* part); 94 95 inline G4double EnergyAfterStep(G4double kinEnergy, G4double step, 96 const G4Material*, 97 const G4String& particleName); 98 99 inline G4double EnergyBeforeStep(G4double kinEnergy, G4double step, 100 const G4Material*, 101 const G4String& particleName); 102 103 G4double AverageScatteringAngle(G4double kinEnergy, G4double step, 104 const G4Material*, 105 const G4ParticleDefinition* part); 106 107 inline G4double AverageScatteringAngle(G4double kinEnergy, G4double step, 108 const G4Material*, 109 const G4String& particleName); 110 111 inline G4double ComputeTrueStep(const G4Material*, 112 const G4ParticleDefinition* part, 113 G4double kinEnergy, G4double stepLength); 114 115 G4double EnergyDispersion(G4double kinEnergy, G4double step, 116 const G4Material*, 117 const G4ParticleDefinition*); 118 119 inline G4double EnergyDispersion(G4double kinEnergy, G4double step, 120 const G4Material*, 121 const G4String& particleName); 122 123 inline void SetVerbose(G4int val); 124 125 inline void SetMinKinEnergy(G4double); 126 127 inline void SetMaxKinEnergy(G4double); 128 129 inline void SetMaxEnergyTransfer(G4double); 130 131 // hide assignment operator 132 G4EnergyLossForExtrapolator & operator= 133 (const G4EnergyLossForExtrapolator &right) = delete; 134 G4EnergyLossForExtrapolator(const G4EnergyLossForExtrapolator&) = delete; 135 136 private: 137 138 G4bool SetupKinematics(const G4ParticleDefinition*, const G4Material*, 139 G4double kinEnergy); 140 141 const G4ParticleDefinition* FindParticle(const G4String& name); 142 143 inline G4double ComputeValue(G4double x, const G4PhysicsTable* table, 144 size_t idxMat); 145 146 inline const G4PhysicsTable* GetPhysicsTable(ExtTableType type) const; 147 148 #ifdef G4MULTITHREADED 149 static G4Mutex extrMutex; 150 #endif 151 static G4TablesForExtrapolator* tables; 152 153 const G4ParticleDefinition* currentParticle = nullptr; 154 const G4ParticleDefinition* electron = nullptr; 155 const G4ParticleDefinition* positron = nullptr; 156 const G4ParticleDefinition* muonPlus = nullptr; 157 const G4ParticleDefinition* muonMinus= nullptr; 158 const G4ParticleDefinition* proton = nullptr; 159 const G4Material* currentMaterial = nullptr; 160 161 G4double electronDensity = 0.0; 162 G4double radLength = 0.0; 163 G4double charge2 = 0.0; 164 G4double kineticEnergy = 0.0; 165 G4double gam = 1.0; 166 G4double bg2 = 0.0; 167 G4double beta2 = 0.0; 168 G4double tmax = 0.0; 169 170 G4double linLossLimit = 0.01; 171 G4double emin = 0.0; 172 G4double emax = 0.0; 173 G4double maxEnergyTransfer = 0.0; 174 175 size_t index = 0; 176 size_t nmat = 0; 177 G4int nbins = 80; 178 G4int verbose = 0; 179 180 G4bool isMaster = false; 181 }; 182 183 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 184 185 inline const G4PhysicsTable* 186 G4EnergyLossForExtrapolator::GetPhysicsTable(ExtTableType type) const 187 { 188 return tables->GetPhysicsTable(type); 189 } 190 191 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 192 193 inline G4double 194 G4EnergyLossForExtrapolator::EnergyAfterStep(G4double kinEnergy, 195 G4double step, 196 const G4Material* mat, 197 const G4String& name) 198 { 199 return EnergyAfterStep(kinEnergy,step,mat,FindParticle(name)); 200 } 201 202 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 203 204 inline G4double 205 G4EnergyLossForExtrapolator::EnergyBeforeStep(G4double kinEnergy, 206 G4double step, 207 const G4Material* mat, 208 const G4String& name) 209 { 210 return EnergyBeforeStep(kinEnergy,step,mat,FindParticle(name)); 211 } 212 213 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 214 215 inline G4double 216 G4EnergyLossForExtrapolator::AverageScatteringAngle(G4double kinEnergy, 217 G4double step, 218 const G4Material* mat, 219 const G4String& name) 220 { 221 return AverageScatteringAngle(kinEnergy,step,mat,FindParticle(name)); 222 } 223 224 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 225 226 inline G4double 227 G4EnergyLossForExtrapolator::EnergyDispersion(G4double kinEnergy, 228 G4double step, 229 const G4Material* mat, 230 const G4String& name) 231 { 232 return EnergyDispersion(kinEnergy,step,mat,FindParticle(name)); 233 } 234 235 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 236 237 inline G4double 238 G4EnergyLossForExtrapolator::ComputeTrueStep(const G4Material* mat, 239 const G4ParticleDefinition* part, 240 G4double kinEnergy, 241 G4double stepLength) 242 { 243 G4double theta = AverageScatteringAngle(kinEnergy,stepLength,mat,part); 244 return stepLength*std::sqrt(1.0 + 0.625*theta*theta); 245 } 246 247 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 248 249 inline G4double 250 G4EnergyLossForExtrapolator::ComputeValue(G4double x, 251 const G4PhysicsTable* table, 252 size_t idxMat) 253 { 254 return (nullptr != table) ? ((*table)[idxMat])->Value(x, index) : 0.0; 255 } 256 257 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 258 259 inline void G4EnergyLossForExtrapolator::SetVerbose(G4int val) 260 { 261 verbose = val; 262 } 263 264 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 265 266 inline void G4EnergyLossForExtrapolator::SetMinKinEnergy(G4double val) 267 { 268 emin = val; 269 } 270 271 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 272 273 inline void G4EnergyLossForExtrapolator::SetMaxKinEnergy(G4double val) 274 { 275 emax = val; 276 } 277 278 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 279 280 inline void G4EnergyLossForExtrapolator::SetMaxEnergyTransfer(G4double val) 281 { 282 maxEnergyTransfer = val; 283 } 284 285 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 286 287 #endif 288 289