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 G4EmBuilder 28 // 29 // Author V.Ivanchenko 22.05.2020 30 // 31 32 #include "G4EmBuilder.hh" 33 34 #include "G4SystemOfUnits.hh" 35 #include "G4ParticleDefinition.hh" 36 #include "G4PhysListUtil.hh" 37 #include "G4EmParameters.hh" 38 #include "G4VEnergyLossProcess.hh" 39 40 #include "G4eMultipleScattering.hh" 41 #include "G4MuMultipleScattering.hh" 42 #include "G4hMultipleScattering.hh" 43 #include "G4CoulombScattering.hh" 44 #include "G4WentzelVIModel.hh" 45 46 #include "G4ProcessManager.hh" 47 #include "G4TransportationWithMsc.hh" 48 #include "G4TransportationProcessType.hh" 49 50 #include "G4MuBremsstrahlungModel.hh" 51 #include "G4MuPairProductionModel.hh" 52 #include "G4hBremsstrahlungModel.hh" 53 #include "G4hPairProductionModel.hh" 54 55 #include "G4MuIonisation.hh" 56 #include "G4MuBremsstrahlung.hh" 57 #include "G4MuPairProduction.hh" 58 #include "G4hBremsstrahlung.hh" 59 #include "G4hPairProduction.hh" 60 61 #include "G4hIonisation.hh" 62 #include "G4ionIonisation.hh" 63 #include "G4NuclearStopping.hh" 64 65 #include "G4MuMultipleScattering.hh" 66 #include "G4hMultipleScattering.hh" 67 68 #include "G4ParticleTable.hh" 69 #include "G4Gamma.hh" 70 #include "G4Electron.hh" 71 #include "G4Positron.hh" 72 73 #include "G4ChargedGeantino.hh" 74 #include "G4Geantino.hh" 75 #include "G4NeutrinoMu.hh" 76 #include "G4AntiNeutrinoMu.hh" 77 #include "G4NeutrinoE.hh" 78 #include "G4AntiNeutrinoE.hh" 79 80 #include "G4MuonPlus.hh" 81 #include "G4MuonMinus.hh" 82 #include "G4PionPlus.hh" 83 #include "G4PionMinus.hh" 84 #include "G4PionZero.hh" 85 #include "G4KaonPlus.hh" 86 #include "G4KaonMinus.hh" 87 #include "G4Proton.hh" 88 #include "G4AntiProton.hh" 89 #include "G4Lambda.hh" 90 #include "G4AntiLambda.hh" 91 92 #include "G4Deuteron.hh" 93 #include "G4Triton.hh" 94 #include "G4He3.hh" 95 #include "G4Alpha.hh" 96 #include "G4GenericIon.hh" 97 98 #include "G4PhysicsListHelper.hh" 99 #include "G4HadParticles.hh" 100 #include "G4HadronicParameters.hh" 101 #include "G4LossTableManager.hh" 102 #include "G4UAtomicDeexcitation.hh" 103 104 void G4EmBuilder::ConstructBasicEmPhysics(G4hMultipleScattering* hmsc, 105 const std::vector<G4int>& partList) 106 { 107 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper(); 108 G4ParticleTable* table = G4ParticleTable::GetParticleTable(); 109 110 for( auto & pdg : partList ) { 111 auto part = table->FindParticle( pdg ); 112 if ( part == nullptr || part->GetPDGCharge() == 0.0 ) { continue; } 113 ph->RegisterProcess(hmsc, part); 114 ph->RegisterProcess(new G4hIonisation(), part); 115 } 116 } 117 118 void G4EmBuilder::ConstructIonEmPhysics(G4hMultipleScattering* hmsc, 119 G4NuclearStopping* nucStopping) 120 { 121 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper(); 122 123 G4ParticleDefinition* part = G4Deuteron::Deuteron(); 124 ph->RegisterProcess(hmsc, part); 125 ph->RegisterProcess(new G4hIonisation(), part); 126 127 part = G4Triton::Triton(); 128 ph->RegisterProcess(hmsc, part); 129 ph->RegisterProcess(new G4hIonisation(), part); 130 131 part = G4Alpha::Alpha(); 132 ph->RegisterProcess(new G4hMultipleScattering(), part); 133 ph->RegisterProcess(new G4ionIonisation(), part); 134 if( nucStopping != nullptr ) { 135 ph->RegisterProcess(nucStopping, part); 136 } 137 138 part = G4He3::He3(); 139 ph->RegisterProcess(new G4hMultipleScattering(), part); 140 ph->RegisterProcess(new G4ionIonisation(), part); 141 if( nucStopping != nullptr ) { 142 ph->RegisterProcess(nucStopping, part); 143 } 144 } 145 146 void G4EmBuilder::ConstructIonEmPhysicsSS() 147 { 148 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper(); 149 150 G4ParticleDefinition* part = G4Deuteron::Deuteron(); 151 ph->RegisterProcess(new G4hIonisation(), part); 152 ph->RegisterProcess(new G4CoulombScattering(false), part); 153 154 part = G4Triton::Triton(); 155 ph->RegisterProcess(new G4hIonisation(), part); 156 ph->RegisterProcess(new G4CoulombScattering(false), part); 157 158 part = G4Alpha::Alpha(); 159 ph->RegisterProcess(new G4ionIonisation(), part); 160 ph->RegisterProcess(new G4CoulombScattering(false), part); 161 162 part = G4He3::He3(); 163 ph->RegisterProcess(new G4ionIonisation(), part); 164 ph->RegisterProcess(new G4CoulombScattering(false), part); 165 } 166 167 void G4EmBuilder::ConstructLightHadrons(G4ParticleDefinition* part1, 168 G4ParticleDefinition* part2, 169 G4bool isHEP, G4bool isProton, 170 G4bool isWVI) 171 { 172 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper(); 173 174 G4hMultipleScattering* msc = new G4hMultipleScattering(); 175 if(isWVI) { msc->SetEmModel(new G4WentzelVIModel()); } 176 G4CoulombScattering* ss = ( isWVI ) ? new G4CoulombScattering() : nullptr; 177 178 ph->RegisterProcess(msc, part1); 179 ph->RegisterProcess(new G4hIonisation(), part1); 180 181 G4hBremsstrahlung* brem = ( isHEP ) ? new G4hBremsstrahlung() : nullptr; 182 G4hPairProduction* pair = ( isHEP ) ? new G4hPairProduction() : nullptr; 183 184 if( isHEP ) { 185 ph->RegisterProcess(brem, part1); 186 ph->RegisterProcess(pair, part1); 187 } 188 if( isWVI ) { ph->RegisterProcess(ss, part1); } 189 190 if( isProton ) { 191 msc = new G4hMultipleScattering(); 192 if(isWVI) { 193 msc->SetEmModel(new G4WentzelVIModel()); 194 ss = new G4CoulombScattering(); 195 } 196 } 197 ph->RegisterProcess(msc, part2); 198 ph->RegisterProcess(new G4hIonisation(), part2); 199 if( isHEP ) { 200 ph->RegisterProcess(brem, part2); 201 ph->RegisterProcess(pair, part2); 202 } 203 if( isWVI ) { ph->RegisterProcess(ss, part2); } 204 } 205 206 void G4EmBuilder::ConstructLightHadronsSS(G4ParticleDefinition* part1, 207 G4ParticleDefinition* part2, 208 G4bool isHEP) 209 { 210 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper(); 211 212 ph->RegisterProcess(new G4hIonisation(), part1); 213 214 G4hBremsstrahlung* brem = ( isHEP ) ? new G4hBremsstrahlung() : nullptr; 215 G4hPairProduction* pair = ( isHEP ) ? new G4hPairProduction() : nullptr; 216 217 if( isHEP ) { 218 ph->RegisterProcess(brem, part1); 219 ph->RegisterProcess(pair, part1); 220 } 221 auto ss = new G4CoulombScattering(false); 222 ph->RegisterProcess(ss, part1); 223 224 ph->RegisterProcess(new G4hIonisation(), part2); 225 if( isHEP ) { 226 ph->RegisterProcess(brem, part2); 227 ph->RegisterProcess(pair, part2); 228 } 229 ph->RegisterProcess(new G4CoulombScattering(false), part2); 230 } 231 232 void G4EmBuilder::ConstructCharged(G4hMultipleScattering* hmsc, 233 G4NuclearStopping* nucStopping, 234 G4bool isWVI) 235 { 236 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper(); 237 G4EmParameters* param = G4EmParameters::Instance(); 238 G4HadronicParameters* hpar = G4HadronicParameters::Instance(); 239 G4bool isHEP = ( param->MaxKinEnergy() > hpar->EnergyThresholdForHeavyHadrons() ); 240 241 // muon multiple and single scattering 242 G4MuMultipleScattering* mumsc = new G4MuMultipleScattering(); 243 if(isWVI) { mumsc->SetEmModel(new G4WentzelVIModel()); } 244 G4CoulombScattering* muss = ( isWVI ) ? new G4CoulombScattering() : nullptr; 245 246 // Add standard EM Processes 247 // mu+- 248 G4ParticleDefinition* part = G4MuonPlus::MuonPlus(); 249 ph->RegisterProcess(mumsc, part); 250 ph->RegisterProcess(new G4MuIonisation(), part); 251 252 // muon bremsstrahlung and pair production 253 G4MuBremsstrahlung* mub = ( isHEP ) ? new G4MuBremsstrahlung() : nullptr; 254 G4MuPairProduction* mup = ( isHEP ) ? new G4MuPairProduction() : nullptr; 255 256 if( isHEP ) { 257 ph->RegisterProcess(mub, part); 258 ph->RegisterProcess(mup, part); 259 } 260 if( isWVI ) { ph->RegisterProcess(muss, part); } 261 262 part = G4MuonMinus::MuonMinus(); 263 ph->RegisterProcess(mumsc, part); 264 ph->RegisterProcess(new G4MuIonisation(), part); 265 if( isHEP ) { 266 ph->RegisterProcess(mub, part); 267 ph->RegisterProcess(mup, part); 268 } 269 if( isWVI ) { ph->RegisterProcess(muss, part); } 270 271 // pi+- 272 ConstructLightHadrons(G4PionPlus::PionPlus(), G4PionMinus::PionMinus(), isHEP, false, isWVI); 273 274 // K+- 275 ConstructLightHadrons(G4KaonPlus::KaonPlus(), G4KaonMinus::KaonMinus(), isHEP, false, isWVI); 276 277 // p, pbar 278 ConstructLightHadrons(G4Proton::Proton(), G4AntiProton::AntiProton(), isHEP, true, isWVI); 279 if( nucStopping != nullptr ) { 280 ph->RegisterProcess(nucStopping, G4Proton::Proton()); 281 } 282 283 // ions 284 ConstructIonEmPhysics(hmsc, nucStopping); 285 286 // hyperons and anti particles 287 if( isHEP ) { 288 ConstructBasicEmPhysics(hmsc, G4HadParticles::GetHeavyChargedParticles()); 289 290 // b- and c- charged particles 291 if( hpar->EnableBCParticles() ) { 292 ConstructBasicEmPhysics(hmsc, G4HadParticles::GetBCChargedHadrons()); 293 } 294 // light hyper-nuclei 295 if( hpar->EnableHyperNuclei() ) { 296 ConstructBasicEmPhysics(hmsc, G4HadParticles::GetChargedHyperNuclei()); 297 } 298 } 299 } 300 301 void G4EmBuilder::ConstructChargedSS(G4hMultipleScattering* hmsc) 302 { 303 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper(); 304 G4EmParameters* param = G4EmParameters::Instance(); 305 G4HadronicParameters* hpar = G4HadronicParameters::Instance(); 306 G4bool isHEP = ( param->MaxKinEnergy() > hpar->EnergyThresholdForHeavyHadrons() ); 307 308 // muon multiple and single scattering 309 G4CoulombScattering* muss = new G4CoulombScattering(false); 310 311 // Add standard EM Processes 312 // mu+- 313 G4ParticleDefinition* part = G4MuonPlus::MuonPlus(); 314 ph->RegisterProcess(new G4MuIonisation(), part); 315 316 // muon bremsstrahlung and pair production 317 G4MuBremsstrahlung* mub = ( isHEP ) ? new G4MuBremsstrahlung() : nullptr; 318 G4MuPairProduction* mup = ( isHEP ) ? new G4MuPairProduction() : nullptr; 319 320 if( isHEP ) { 321 ph->RegisterProcess(mub, part); 322 ph->RegisterProcess(mup, part); 323 } 324 ph->RegisterProcess(muss, part); 325 326 part = G4MuonMinus::MuonMinus(); 327 ph->RegisterProcess(new G4MuIonisation(), part); 328 if( isHEP ) { 329 ph->RegisterProcess(mub, part); 330 ph->RegisterProcess(mup, part); 331 } 332 ph->RegisterProcess(muss, part); 333 334 // pi+- 335 ConstructLightHadronsSS(G4PionPlus::PionPlus(), G4PionMinus::PionMinus(), isHEP); 336 337 // K+- 338 ConstructLightHadronsSS(G4KaonPlus::KaonPlus(), G4KaonMinus::KaonMinus(), isHEP); 339 340 // p, pbar 341 ConstructLightHadronsSS(G4Proton::Proton(), G4AntiProton::AntiProton(), isHEP); 342 // ions 343 ConstructIonEmPhysicsSS(); 344 345 // hyperons and anti particles 346 if( isHEP ) { 347 ConstructBasicEmPhysics(hmsc, G4HadParticles::GetHeavyChargedParticles()); 348 349 // b- and c- charged particles 350 if( hpar->EnableBCParticles() ) { 351 ConstructBasicEmPhysics(hmsc, G4HadParticles::GetBCChargedHadrons()); 352 } 353 // light hyper-nuclei 354 if( hpar->EnableHyperNuclei() ) { 355 ConstructBasicEmPhysics(hmsc, G4HadParticles::GetChargedHyperNuclei()); 356 } 357 } 358 } 359 360 void G4EmBuilder::ConstructMinimalEmSet() 361 { 362 // instantiate singletones for physics 363 G4PhysListUtil::InitialiseParameters(); 364 // pseudo-particles 365 G4Geantino::GeantinoDefinition(); 366 G4ChargedGeantino::ChargedGeantinoDefinition(); 367 G4NeutrinoMu::NeutrinoMu(); 368 G4AntiNeutrinoMu::AntiNeutrinoMu(); 369 G4NeutrinoE::NeutrinoE(); 370 G4AntiNeutrinoE::AntiNeutrinoE(); 371 // gamma 372 G4Gamma::Gamma(); 373 // leptons 374 G4Electron::Electron(); 375 G4Positron::Positron(); 376 G4MuonPlus::MuonPlus(); 377 G4MuonMinus::MuonMinus(); 378 // mesons 379 G4PionPlus::PionPlus(); 380 G4PionMinus::PionMinus(); 381 G4PionZero::PionZero(); 382 G4KaonPlus::KaonPlus(); 383 G4KaonMinus::KaonMinus(); 384 // barions 385 G4Proton::Proton(); 386 G4AntiProton::AntiProton(); 387 G4Neutron::Neutron(); 388 G4AntiNeutron::AntiNeutron(); 389 G4Lambda::Lambda(); 390 G4AntiLambda::AntiLambda(); 391 // ions 392 G4Deuteron::Deuteron(); 393 G4Triton::Triton(); 394 G4He3::He3(); 395 G4Alpha::Alpha(); 396 G4GenericIon::GenericIon(); 397 } 398 399 void G4EmBuilder::PrepareEMPhysics() 400 { 401 G4LossTableManager* man = G4LossTableManager::Instance(); 402 G4VAtomDeexcitation* ad = man->AtomDeexcitation(); 403 if(nullptr == ad) { 404 ad = new G4UAtomicDeexcitation(); 405 man->SetAtomDeexcitation(ad); 406 } 407 } 408 409 void G4EmBuilder::ConstructElectronMscProcess(G4VMscModel* msc1, G4VMscModel* msc2, 410 G4ParticleDefinition* particle) 411 { 412 G4TransportationWithMscType type = 413 G4EmParameters::Instance()->TransportationWithMsc(); 414 G4ProcessManager* procManager = particle->GetProcessManager(); 415 auto plist = procManager->GetProcessList(); 416 G4int ptype = (0 < plist->size()) ? (*plist)[0]->GetProcessSubType() : 0; 417 if(type != G4TransportationWithMscType::fDisabled && ptype == TRANSPORTATION) { 418 // Remove default G4Transportation and replace with G4TransportationWithMsc. 419 procManager->RemoveProcess(0); 420 G4TransportationWithMsc* transportWithMsc = new G4TransportationWithMsc( 421 G4TransportationWithMsc::ScatteringType::MultipleScattering); 422 if(type == G4TransportationWithMscType::fMultipleSteps) { 423 transportWithMsc->SetMultipleSteps(true); 424 } 425 transportWithMsc->AddMscModel(msc1); 426 if(msc2 != nullptr) { 427 transportWithMsc->AddMscModel(msc2); 428 } 429 procManager->AddProcess(transportWithMsc, -1, 0, 0); 430 } else { 431 // Register as a separate process. 432 G4eMultipleScattering* msc = new G4eMultipleScattering; 433 msc->SetEmModel(msc1); 434 if(msc2 != nullptr) { 435 msc->SetEmModel(msc2); 436 } 437 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper(); 438 ph->RegisterProcess(msc, particle); 439 } 440 } 441 442 void G4EmBuilder::ConstructElectronSSProcess(G4VEmModel* ss, G4ParticleDefinition* particle) 443 { 444 G4TransportationWithMscType type = G4EmParameters::Instance()->TransportationWithMsc(); 445 G4ProcessManager* procManager = particle->GetProcessManager(); 446 auto plist = procManager->GetProcessList(); 447 G4int ptype = (0 < plist->size()) ? (*plist)[0]->GetProcessSubType() : 0; 448 if (type != G4TransportationWithMscType::fDisabled && ptype == TRANSPORTATION) { 449 // Remove default G4Transportation and replace with G4TransportationWithMsc. 450 procManager->RemoveProcess(0); 451 G4TransportationWithMsc* transportWithMsc = 452 new G4TransportationWithMsc(G4TransportationWithMsc::ScatteringType::SingleScattering); 453 if (type == G4TransportationWithMscType::fMultipleSteps) { 454 transportWithMsc->SetMultipleSteps(true); 455 } 456 transportWithMsc->AddSSModel(ss); 457 procManager->AddProcess(transportWithMsc, -1, 0, 0); 458 } 459 else { 460 // Register as a separate process. 461 G4CoulombScattering* ssProc = new G4CoulombScattering(false); 462 ssProc->SetEmModel(ss); 463 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper(); 464 ph->RegisterProcess(ssProc, particle); 465 } 466 } 467