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 #include "G4EmLivermorePhysics.hh" 28 #include "G4ParticleDefinition.hh" 29 #include "G4SystemOfUnits.hh" 30 31 // *** Processes and models 32 // gamma 33 #include "G4PhotoElectricEffect.hh" 34 #include "G4LivermorePhotoElectricModel.hh" 35 #include "G4PhotoElectricAngularGeneratorPolarized.hh" 36 37 #include "G4ComptonScattering.hh" 38 #include "G4LivermoreComptonModel.hh" 39 #include "G4LivermorePolarizedComptonModel.hh" 40 #include "G4KleinNishinaModel.hh" 41 #include "G4LowEPPolarizedComptonModel.hh" 42 43 #include "G4GammaConversion.hh" 44 #include "G4LivermoreGammaConversionModel.hh" 45 #include "G4LivermoreGammaConversion5DModel.hh" 46 47 #include "G4RayleighScattering.hh" 48 #include "G4LivermoreRayleighModel.hh" 49 #include "G4LivermorePolarizedRayleighModel.hh" 50 51 #include "G4PEEffectFluoModel.hh" 52 #include "G4KleinNishinaModel.hh" 53 54 // e+- 55 #include "G4eMultipleScattering.hh" 56 #include "G4eIonisation.hh" 57 #include "G4LivermoreIonisationModel.hh" 58 59 #include "G4eBremsstrahlung.hh" 60 #include "G4LivermoreBremsstrahlungModel.hh" 61 #include "G4Generator2BS.hh" 62 #include "G4SeltzerBergerModel.hh" 63 #include "G4BetheHeitler5DModel.hh" 64 65 // e+ 66 #include "G4eplusAnnihilation.hh" 67 #include "G4eplusTo2or3GammaModel.hh" 68 69 // hadrons 70 #include "G4hMultipleScattering.hh" 71 #include "G4MscStepLimitType.hh" 72 73 #include "G4ePairProduction.hh" 74 75 #include "G4hIonisation.hh" 76 #include "G4ionIonisation.hh" 77 #include "G4IonParametrisedLossModel.hh" 78 #include "G4NuclearStopping.hh" 79 #include "G4LindhardSorensenIonModel.hh" 80 81 // msc models 82 #include "G4UrbanMscModel.hh" 83 #include "G4WentzelVIModel.hh" 84 #include "G4GoudsmitSaundersonMscModel.hh" 85 #include "G4CoulombScattering.hh" 86 #include "G4eCoulombScatteringModel.hh" 87 88 // interfaces 89 #include "G4EmParameters.hh" 90 #include "G4EmBuilder.hh" 91 #include "G4EmStandUtil.hh" 92 93 // particles 94 95 #include "G4Gamma.hh" 96 #include "G4Electron.hh" 97 #include "G4Positron.hh" 98 #include "G4GenericIon.hh" 99 100 // 101 #include "G4PhysicsListHelper.hh" 102 #include "G4BuilderType.hh" 103 #include "G4EmModelActivator.hh" 104 105 // factory 106 #include "G4PhysicsConstructorFactory.hh" 107 // 108 G4_DECLARE_PHYSCONSTR_FACTORY(G4EmLivermorePhysics); 109 110 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 111 112 G4EmLivermorePhysics::G4EmLivermorePhysics(G4int ver, const G4String& pname) 113 : G4VPhysicsConstructor(pname) 114 { 115 SetVerboseLevel(ver); 116 G4EmParameters* param = G4EmParameters::Instance(); 117 param->SetDefaults(); 118 param->SetVerbose(ver); 119 param->SetMinEnergy(100*CLHEP::eV); 120 param->SetLowestElectronEnergy(100*CLHEP::eV); 121 param->SetNumberOfBinsPerDecade(20); 122 param->ActivateAngularGeneratorForIonisation(true); 123 param->SetStepFunction(0.2, 10*CLHEP::um); 124 param->SetStepFunctionMuHad(0.1, 50*CLHEP::um); 125 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um); 126 param->SetStepFunctionIons(0.1, 1*CLHEP::um); 127 param->SetUseMottCorrection(true); 128 param->SetMscStepLimitType(fUseSafetyPlus); 129 param->SetMscSkin(3); 130 param->SetMscRangeFactor(0.08); 131 param->SetMuHadLateralDisplacement(true); 132 param->SetFluo(true); 133 param->SetUseICRU90Data(true); 134 param->SetFluctuationType(fUrbanFluctuation); 135 param->SetMaxNIELEnergy(1*CLHEP::MeV); 136 param->SetPositronAtRestModelType(fAllisonPositronium); 137 SetPhysicsType(bElectromagnetic); 138 } 139 140 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 141 142 G4EmLivermorePhysics::~G4EmLivermorePhysics() = default; 143 144 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 145 146 void G4EmLivermorePhysics::ConstructParticle() 147 { 148 // minimal set of particles for EM physics 149 G4EmBuilder::ConstructMinimalEmSet(); 150 } 151 152 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 153 154 void G4EmLivermorePhysics::ConstructProcess() 155 { 156 if(verboseLevel > 1) { 157 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl; 158 } 159 G4EmBuilder::PrepareEMPhysics(); 160 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper(); 161 G4EmParameters* param = G4EmParameters::Instance(); 162 163 // processes used by several particles 164 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc"); 165 166 // high energy limit for e+- scattering models 167 G4double highEnergyLimit= param->MscEnergyLimit(); 168 G4double livEnergyLimit = 1*CLHEP::GeV; 169 170 // nuclear stopping 171 G4double nielEnergyLimit = param->MaxNIELEnergy(); 172 G4NuclearStopping* pnuc = nullptr; 173 if(nielEnergyLimit > 0.0) { 174 pnuc = new G4NuclearStopping(); 175 pnuc->SetMaxKinEnergy(nielEnergyLimit); 176 } 177 178 // Add Livermore EM Processes 179 180 // Add gamma EM Processes 181 G4ParticleDefinition* particle = G4Gamma::Gamma(); 182 G4bool polar = param->EnablePolarisation(); 183 184 // photoelectric effect - Livermore model only 185 G4PhotoElectricEffect* pe = new G4PhotoElectricEffect(); 186 G4VEmModel* peModel = new G4LivermorePhotoElectricModel(); 187 pe->SetEmModel(peModel); 188 if(polar) { 189 peModel->SetAngularDistribution(new G4PhotoElectricAngularGeneratorPolarized()); 190 } 191 192 // Compton scattering - Livermore model only 193 G4ComptonScattering* cs = new G4ComptonScattering; 194 cs->SetEmModel(new G4KleinNishinaModel()); 195 G4VEmModel* cModel = nullptr; 196 if(polar) { 197 cModel = new G4LivermorePolarizedComptonModel(); 198 cModel->SetHighEnergyLimit(20*CLHEP::MeV); 199 } else { 200 cModel = new G4LivermoreComptonModel(); 201 cModel->SetHighEnergyLimit(livEnergyLimit); 202 } 203 cs->AddEmModel(0, cModel); 204 205 // gamma conversion 206 G4GammaConversion* gc = new G4GammaConversion(); 207 G4VEmModel* convLiv = new G4LivermoreGammaConversion5DModel(); 208 gc->SetEmModel(convLiv); 209 210 // Rayleigh scattering 211 G4RayleighScattering* rl = new G4RayleighScattering(); 212 if(polar) { 213 rl->SetEmModel(new G4LivermorePolarizedRayleighModel()); 214 } 215 216 ph->RegisterProcess(pe, particle); 217 ph->RegisterProcess(cs, particle); 218 ph->RegisterProcess(gc, particle); 219 ph->RegisterProcess(rl, particle); 220 221 // e- 222 particle = G4Electron::Electron(); 223 224 // multiple and single scattering 225 G4GoudsmitSaundersonMscModel* msc1 = new G4GoudsmitSaundersonMscModel(); 226 G4WentzelVIModel* msc2 = new G4WentzelVIModel(); 227 msc1->SetHighEnergyLimit(highEnergyLimit); 228 msc2->SetLowEnergyLimit(highEnergyLimit); 229 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle); 230 231 G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel(); 232 G4CoulombScattering* ss = new G4CoulombScattering(); 233 ss->SetEmModel(ssm); 234 ss->SetMinKinEnergy(highEnergyLimit); 235 ssm->SetLowEnergyLimit(highEnergyLimit); 236 ssm->SetActivationLowEnergyLimit(highEnergyLimit); 237 238 // Ionisation - Livermore should be used only for low energies 239 G4eIonisation* eioni = new G4eIonisation(); 240 eioni->SetFluctModel(G4EmStandUtil::ModelOfFluctuations()); 241 G4VEmModel* theIoniLiv = new G4LivermoreIonisationModel(); 242 theIoniLiv->SetHighEnergyLimit(0.1*CLHEP::MeV); 243 eioni->AddEmModel(0, theIoniLiv); 244 245 // Bremsstrahlung 246 G4eBremsstrahlung* brem = new G4eBremsstrahlung(); 247 G4SeltzerBergerModel* br1 = new G4SeltzerBergerModel(); 248 G4eBremsstrahlungRelModel* br2 = new G4eBremsstrahlungRelModel(); 249 br1->SetAngularDistribution(new G4Generator2BS()); 250 br2->SetAngularDistribution(new G4Generator2BS()); 251 brem->SetEmModel(br1); 252 brem->SetEmModel(br2); 253 br1->SetHighEnergyLimit(GeV); 254 255 G4ePairProduction* ee = new G4ePairProduction(); 256 257 // register processes 258 ph->RegisterProcess(eioni, particle); 259 ph->RegisterProcess(brem, particle); 260 ph->RegisterProcess(ee, particle); 261 ph->RegisterProcess(ss, particle); 262 263 // e+ 264 particle = G4Positron::Positron(); 265 266 // multiple and single scattering 267 msc1 = new G4GoudsmitSaundersonMscModel(); 268 msc2 = new G4WentzelVIModel(); 269 msc1->SetHighEnergyLimit(highEnergyLimit); 270 msc2->SetLowEnergyLimit(highEnergyLimit); 271 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle); 272 273 ssm = new G4eCoulombScatteringModel(); 274 ss = new G4CoulombScattering(); 275 ss->SetEmModel(ssm); 276 ss->SetMinKinEnergy(highEnergyLimit); 277 ssm->SetLowEnergyLimit(highEnergyLimit); 278 ssm->SetActivationLowEnergyLimit(highEnergyLimit); 279 280 // ionisation 281 eioni = new G4eIonisation(); 282 283 // Bremsstrahlung from standard 284 brem = new G4eBremsstrahlung(); 285 br1 = new G4SeltzerBergerModel(); 286 br2 = new G4eBremsstrahlungRelModel(); 287 br1->SetAngularDistribution(new G4Generator2BS()); 288 br2->SetAngularDistribution(new G4Generator2BS()); 289 brem->SetEmModel(br1); 290 brem->SetEmModel(br2); 291 br1->SetHighEnergyLimit(GeV); 292 293 // annihilation 294 auto anni = new G4eplusAnnihilation(); 295 if (param->Use3GammaAnnihilationOnFly()) { 296 anni->SetEmModel(new G4eplusTo2or3GammaModel()); 297 } 298 299 // register processes 300 ph->RegisterProcess(eioni, particle); 301 ph->RegisterProcess(brem, particle); 302 ph->RegisterProcess(ee, particle); 303 ph->RegisterProcess(anni, particle); 304 ph->RegisterProcess(ss, particle); 305 306 // generic ion 307 particle = G4GenericIon::GenericIon(); 308 G4ionIonisation* ionIoni = new G4ionIonisation(); 309 ionIoni->SetEmModel(new G4LindhardSorensenIonModel()); 310 ph->RegisterProcess(hmsc, particle); 311 ph->RegisterProcess(ionIoni, particle); 312 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); } 313 314 // muons, hadrons, ions 315 G4EmBuilder::ConstructCharged(hmsc, pnuc); 316 317 // extra configuration 318 G4EmModelActivator mact(GetPhysicsName()); 319 320 } 321 322 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 323