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