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 /// \file biasing/B03/src/B03PhysicsList.cc 27 /// \brief Implementation of the B03PhysicsList class 28 // 29 // 30 // 31 32 #include "B03PhysicsList.hh" 33 34 #include "G4BaryonConstructor.hh" 35 #include "G4BosonConstructor.hh" 36 #include "G4HadronicParameters.hh" 37 #include "G4IonConstructor.hh" 38 #include "G4LeptonConstructor.hh" 39 #include "G4Material.hh" 40 #include "G4MaterialTable.hh" 41 #include "G4MesonConstructor.hh" 42 #include "G4ParticleDefinition.hh" 43 #include "G4ParticleTable.hh" 44 #include "G4ParticleTypes.hh" 45 #include "G4ParticleWithCuts.hh" 46 #include "G4ProcessManager.hh" 47 #include "G4ProcessVector.hh" 48 #include "G4ShortLivedConstructor.hh" 49 #include "G4SystemOfUnits.hh" 50 #include "globals.hh" 51 52 #include <iomanip> 53 54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 55 56 B03PhysicsList::B03PhysicsList(G4String parallelname) 57 : G4VUserPhysicsList(), fBiasWorldName(parallelname) 58 { 59 fParaWorldName.clear(); 60 SetVerboseLevel(1); 61 } 62 63 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 64 65 B03PhysicsList::~B03PhysicsList() 66 { 67 fParaWorldName.clear(); 68 } 69 70 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 71 72 void B03PhysicsList::ConstructParticle() 73 { 74 // In this method, static member functions should be called 75 // for all particles which you want to use. 76 // This ensures that objects of these particle types will be 77 // created in the program. 78 79 ConstructAllBosons(); 80 ConstructAllLeptons(); 81 ConstructAllMesons(); 82 ConstructAllBaryons(); 83 ConstructAllIons(); 84 ConstructAllShortLiveds(); 85 } 86 87 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 88 89 void B03PhysicsList::ConstructAllBosons() 90 { 91 // Construct all bosons 92 G4BosonConstructor pConstructor; 93 pConstructor.ConstructParticle(); 94 } 95 96 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 97 98 void B03PhysicsList::ConstructAllLeptons() 99 { 100 // Construct all leptons 101 G4LeptonConstructor pConstructor; 102 pConstructor.ConstructParticle(); 103 } 104 105 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 106 107 void B03PhysicsList::ConstructAllMesons() 108 { 109 // Construct all mesons 110 G4MesonConstructor pConstructor; 111 pConstructor.ConstructParticle(); 112 } 113 114 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 115 116 void B03PhysicsList::ConstructAllBaryons() 117 { 118 // Construct all barions 119 G4BaryonConstructor pConstructor; 120 pConstructor.ConstructParticle(); 121 } 122 123 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 124 125 void B03PhysicsList::ConstructAllIons() 126 { 127 // Construct light ions 128 G4IonConstructor pConstructor; 129 pConstructor.ConstructParticle(); 130 } 131 132 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 133 134 void B03PhysicsList::ConstructAllShortLiveds() 135 { 136 // Construct resonaces and quarks 137 G4ShortLivedConstructor pConstructor; 138 pConstructor.ConstructParticle(); 139 } 140 141 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 142 143 void B03PhysicsList::ConstructProcess() 144 { 145 AddTransportation(); 146 AddScoringProcess(); 147 AddBiasingProcess(); 148 ConstructEM(); 149 ConstructLeptHad(); 150 ConstructHad(); 151 ConstructGeneral(); 152 } 153 154 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 155 156 #include "G4ComptonScattering.hh" 157 #include "G4GammaConversion.hh" 158 #include "G4MuBremsstrahlung.hh" 159 #include "G4MuIonisation.hh" 160 #include "G4MuMultipleScattering.hh" 161 #include "G4MuPairProduction.hh" 162 #include "G4PhotoElectricEffect.hh" 163 #include "G4eBremsstrahlung.hh" 164 #include "G4eIonisation.hh" 165 #include "G4eMultipleScattering.hh" 166 #include "G4eplusAnnihilation.hh" 167 #include "G4hIonisation.hh" 168 #include "G4hMultipleScattering.hh" 169 170 void B03PhysicsList::ConstructEM() 171 { 172 auto particleIterator = GetParticleIterator(); 173 particleIterator->reset(); 174 while ((*particleIterator)()) { 175 G4ParticleDefinition* particle = particleIterator->value(); 176 G4ProcessManager* pmanager = particle->GetProcessManager(); 177 G4String particleName = particle->GetParticleName(); 178 179 if (particleName == "gamma") { 180 // gamma 181 // Construct processes for gamma 182 pmanager->AddDiscreteProcess(new G4GammaConversion()); 183 pmanager->AddDiscreteProcess(new G4ComptonScattering()); 184 pmanager->AddDiscreteProcess(new G4PhotoElectricEffect()); 185 } 186 else if (particleName == "e-") { 187 // electron 188 // Construct processes for electron 189 pmanager->AddProcess(new G4eMultipleScattering(), -1, 1, 1); 190 pmanager->AddProcess(new G4eIonisation(), -1, 2, 2); 191 pmanager->AddProcess(new G4eBremsstrahlung(), -1, -1, 3); 192 } 193 else if (particleName == "e+") { 194 // positron 195 // Construct processes for positron 196 pmanager->AddProcess(new G4eMultipleScattering(), -1, 1, 1); 197 198 pmanager->AddProcess(new G4eIonisation(), -1, 2, 2); 199 pmanager->AddProcess(new G4eBremsstrahlung(), -1, -1, 3); 200 pmanager->AddProcess(new G4eplusAnnihilation(), 0, -1, 4); 201 } 202 else if (particleName == "mu+" || particleName == "mu-") { 203 // muon 204 // Construct processes for muon+ 205 pmanager->AddProcess(new G4MuMultipleScattering(), -1, 1, 1); 206 pmanager->AddProcess(new G4MuIonisation(), -1, 2, 2); 207 pmanager->AddProcess(new G4MuBremsstrahlung(), -1, -1, 3); 208 pmanager->AddProcess(new G4MuPairProduction(), -1, -1, 4); 209 } 210 else if (particleName == "GenericIon") { 211 pmanager->AddProcess(new G4hMultipleScattering(), -1, 1, 1); 212 pmanager->AddProcess(new G4hIonisation(), -1, 2, 2); 213 } 214 else { 215 if ((particle->GetPDGCharge() != 0.0) && (particle->GetParticleName() != "chargedgeantino") 216 && (!particle->IsShortLived())) 217 { 218 // all others charged particles except geantino 219 pmanager->AddProcess(new G4hMultipleScattering(), -1, 1, 1); 220 pmanager->AddProcess(new G4hIonisation(), -1, 2, 2); 221 } 222 } 223 } 224 } 225 226 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 227 228 // Hadron Processes 229 230 #include "G4HadronElasticProcess.hh" 231 #include "G4HadronInelasticProcess.hh" 232 #include "G4NeutronCaptureProcess.hh" 233 #include "G4NeutronFissionProcess.hh" 234 235 // Low-energy Models 236 237 #include "G4HadronElastic.hh" 238 #include "G4LFission.hh" 239 #include "G4NeutronRadCapture.hh" 240 241 // -- generator models 242 #include "G4BinaryLightIonReaction.hh" 243 #include "G4CascadeInterface.hh" 244 #include "G4CompetitiveFission.hh" 245 #include "G4ExcitationHandler.hh" 246 #include "G4ExcitedStringDecay.hh" 247 #include "G4FTFModel.hh" 248 #include "G4Fancy3DNucleus.hh" 249 #include "G4GeneratorPrecompoundInterface.hh" 250 #include "G4LundStringFragmentation.hh" 251 #include "G4PreCompoundModel.hh" 252 #include "G4QGSMFragmentation.hh" 253 #include "G4QMDReaction.hh" 254 #include "G4StringModel.hh" 255 #include "G4TheoFSGenerator.hh" 256 257 // Cross sections 258 #include "G4ComponentGGHadronNucleusXsc.hh" 259 #include "G4ComponentGGNuclNuclXsc.hh" 260 #include "G4CrossSectionElastic.hh" 261 #include "G4CrossSectionInelastic.hh" 262 #include "G4NeutronInelasticXS.hh" 263 264 // 265 // ConstructHad() 266 // 267 // Makes discrete physics processes for the hadrons 268 // The processes are: Elastic scattering, Inelastic scattering, 269 // Fission (for neutron only), and Capture (neutron). 270 // 271 272 void B03PhysicsList::ConstructHad() 273 { 274 // this will be the model class for high energies 275 G4TheoFSGenerator* theTheoModel = new G4TheoFSGenerator; 276 G4TheoFSGenerator* antiBHighEnergyModel = new G4TheoFSGenerator; 277 278 // Evaporation logic 279 G4ExcitationHandler* theHandler = new G4ExcitationHandler; 280 theHandler->SetMinEForMultiFrag(3 * MeV); 281 282 // Pre equilibrium stage 283 G4PreCompoundModel* thePreEquilib = new G4PreCompoundModel(theHandler); 284 285 // a no-cascade generator-precompound interaface 286 G4GeneratorPrecompoundInterface* theCascade = new G4GeneratorPrecompoundInterface; 287 theCascade->SetDeExcitation(thePreEquilib); 288 289 // Bertini cascade 290 G4CascadeInterface* bertini = new G4CascadeInterface; 291 bertini->SetMaxEnergy(22 * MeV); 292 293 // here come the high energy parts 294 G4VPartonStringModel* theStringModel; 295 theStringModel = new G4FTFModel; 296 theTheoModel->SetTransport(theCascade); 297 theTheoModel->SetHighEnergyGenerator(theStringModel); 298 theTheoModel->SetMinEnergy(19 * GeV); 299 theTheoModel->SetMaxEnergy(G4HadronicParameters::Instance()->GetMaxEnergy()); 300 301 G4VLongitudinalStringDecay* theFragmentation = new G4QGSMFragmentation; 302 G4ExcitedStringDecay* theStringDecay = new G4ExcitedStringDecay(theFragmentation); 303 theStringModel->SetFragmentationModel(theStringDecay); 304 305 // high energy model for anti-baryons 306 antiBHighEnergyModel = new G4TheoFSGenerator("ANTI-FTFP"); 307 G4FTFModel* antiBStringModel = new G4FTFModel; 308 G4ExcitedStringDecay* stringDecay = new G4ExcitedStringDecay(new G4LundStringFragmentation); 309 antiBStringModel->SetFragmentationModel(stringDecay); 310 311 G4GeneratorPrecompoundInterface* antiBCascade = new G4GeneratorPrecompoundInterface; 312 G4PreCompoundModel* preEquilib = new G4PreCompoundModel(new G4ExcitationHandler); 313 antiBCascade->SetDeExcitation(preEquilib); 314 315 antiBHighEnergyModel->SetTransport(antiBCascade); 316 antiBHighEnergyModel->SetHighEnergyGenerator(antiBStringModel); 317 antiBHighEnergyModel->SetMinEnergy(0.0); 318 antiBHighEnergyModel->SetMaxEnergy(20 * TeV); 319 320 // Light ion models 321 G4BinaryLightIonReaction* binaryCascade = new G4BinaryLightIonReaction; 322 binaryCascade->SetMinEnergy(0.0); 323 binaryCascade->SetMaxEnergy(110 * MeV); 324 325 G4QMDReaction* qmd = new G4QMDReaction; 326 qmd->SetMinEnergy(100 * MeV); 327 qmd->SetMaxEnergy(10 * GeV); 328 329 G4VCrossSectionDataSet* ionXS = new G4CrossSectionInelastic(new G4ComponentGGNuclNuclXsc); 330 331 G4ComponentGGHadronNucleusXsc* ggHNXsec = new G4ComponentGGHadronNucleusXsc(); 332 G4VCrossSectionDataSet* theGGHNEl = new G4CrossSectionElastic(ggHNXsec); 333 G4VCrossSectionDataSet* theGGHNInel = new G4CrossSectionInelastic(ggHNXsec); 334 335 // Elastic process 336 G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess; 337 theElasticProcess->AddDataSet(theGGHNEl); 338 G4HadronElastic* theElasticModel = new G4HadronElastic; 339 theElasticProcess->RegisterMe(theElasticModel); 340 341 auto particleIterator = GetParticleIterator(); 342 particleIterator->reset(); 343 while ((*particleIterator)()) { 344 G4ParticleDefinition* particle = particleIterator->value(); 345 G4ProcessManager* pmanager = particle->GetProcessManager(); 346 G4String particleName = particle->GetParticleName(); 347 348 if (particleName == "pi+") { 349 pmanager->AddDiscreteProcess(theElasticProcess); 350 G4HadronInelasticProcess* theInelasticProcess = 351 new G4HadronInelasticProcess("inelastic", G4PionPlus::Definition()); 352 theInelasticProcess->AddDataSet(theGGHNInel); 353 theInelasticProcess->RegisterMe(bertini); 354 theInelasticProcess->RegisterMe(theTheoModel); 355 pmanager->AddDiscreteProcess(theInelasticProcess); 356 } 357 else if (particleName == "pi-") { 358 pmanager->AddDiscreteProcess(theElasticProcess); 359 G4HadronInelasticProcess* theInelasticProcess = 360 new G4HadronInelasticProcess("inelastic", G4PionMinus::Definition()); 361 theInelasticProcess->AddDataSet(theGGHNInel); 362 theInelasticProcess->RegisterMe(bertini); 363 theInelasticProcess->RegisterMe(theTheoModel); 364 pmanager->AddDiscreteProcess(theInelasticProcess); 365 } 366 else if (particleName == "kaon+") { 367 pmanager->AddDiscreteProcess(theElasticProcess); 368 G4HadronInelasticProcess* theInelasticProcess = 369 new G4HadronInelasticProcess("inelastic", G4KaonPlus::Definition()); 370 theInelasticProcess->AddDataSet(theGGHNInel); 371 theInelasticProcess->RegisterMe(bertini); 372 theInelasticProcess->RegisterMe(theTheoModel); 373 pmanager->AddDiscreteProcess(theInelasticProcess); 374 } 375 else if (particleName == "kaon0S") { 376 pmanager->AddDiscreteProcess(theElasticProcess); 377 G4HadronInelasticProcess* theInelasticProcess = 378 new G4HadronInelasticProcess("inelastic", G4KaonZeroShort::Definition()); 379 theInelasticProcess->AddDataSet(theGGHNInel); 380 theInelasticProcess->RegisterMe(bertini); 381 theInelasticProcess->RegisterMe(theTheoModel); 382 pmanager->AddDiscreteProcess(theInelasticProcess); 383 } 384 else if (particleName == "kaon0L") { 385 pmanager->AddDiscreteProcess(theElasticProcess); 386 G4HadronInelasticProcess* theInelasticProcess = 387 new G4HadronInelasticProcess("inelastic", G4KaonZeroLong::Definition()); 388 theInelasticProcess->AddDataSet(theGGHNInel); 389 theInelasticProcess->RegisterMe(bertini); 390 theInelasticProcess->RegisterMe(theTheoModel); 391 pmanager->AddDiscreteProcess(theInelasticProcess); 392 } 393 else if (particleName == "kaon-") { 394 pmanager->AddDiscreteProcess(theElasticProcess); 395 G4HadronInelasticProcess* theInelasticProcess = 396 new G4HadronInelasticProcess("inelastic", G4KaonMinus::Definition()); 397 theInelasticProcess->AddDataSet(theGGHNInel); 398 theInelasticProcess->RegisterMe(bertini); 399 theInelasticProcess->RegisterMe(theTheoModel); 400 pmanager->AddDiscreteProcess(theInelasticProcess); 401 } 402 else if (particleName == "proton") { 403 pmanager->AddDiscreteProcess(theElasticProcess); 404 G4HadronInelasticProcess* theInelasticProcess = 405 new G4HadronInelasticProcess("inelastic", G4Proton::Definition()); 406 theInelasticProcess->AddDataSet(theGGHNInel); 407 theInelasticProcess->RegisterMe(bertini); 408 theInelasticProcess->RegisterMe(theTheoModel); 409 pmanager->AddDiscreteProcess(theInelasticProcess); 410 } 411 else if (particleName == "anti_proton") { 412 pmanager->AddDiscreteProcess(theElasticProcess); 413 G4HadronInelasticProcess* theInelasticProcess = 414 new G4HadronInelasticProcess("inelastic", G4AntiProton::Definition()); 415 theInelasticProcess->AddDataSet(theGGHNInel); 416 theInelasticProcess->RegisterMe(antiBHighEnergyModel); 417 pmanager->AddDiscreteProcess(theInelasticProcess); 418 } 419 else if (particleName == "neutron") { 420 // elastic scattering 421 pmanager->AddDiscreteProcess(theElasticProcess); 422 423 // inelastic scattering 424 G4HadronInelasticProcess* theInelasticProcess = 425 new G4HadronInelasticProcess("inelastic", G4Neutron::Definition()); 426 theInelasticProcess->AddDataSet(new G4NeutronInelasticXS()); 427 theInelasticProcess->RegisterMe(bertini); 428 theInelasticProcess->RegisterMe(theTheoModel); 429 pmanager->AddDiscreteProcess(theInelasticProcess); 430 431 // fission 432 G4NeutronFissionProcess* theFissionProcess = new G4NeutronFissionProcess; 433 G4LFission* theFissionModel = new G4LFission; 434 theFissionProcess->RegisterMe(theFissionModel); 435 pmanager->AddDiscreteProcess(theFissionProcess); 436 437 // capture 438 G4NeutronCaptureProcess* theCaptureProcess = new G4NeutronCaptureProcess; 439 G4NeutronRadCapture* theCaptureModel = new G4NeutronRadCapture; 440 theCaptureProcess->RegisterMe(theCaptureModel); 441 pmanager->AddDiscreteProcess(theCaptureProcess); 442 } 443 else if (particleName == "anti_neutron") { 444 pmanager->AddDiscreteProcess(theElasticProcess); 445 G4HadronInelasticProcess* theInelasticProcess = 446 new G4HadronInelasticProcess("inelastic", G4AntiNeutron::Definition()); 447 theInelasticProcess->AddDataSet(theGGHNInel); 448 theInelasticProcess->RegisterMe(antiBHighEnergyModel); 449 pmanager->AddDiscreteProcess(theInelasticProcess); 450 } 451 else if (particleName == "lambda") { 452 pmanager->AddDiscreteProcess(theElasticProcess); 453 G4HadronInelasticProcess* theInelasticProcess = 454 new G4HadronInelasticProcess("inelastic", G4Lambda::Definition()); 455 theInelasticProcess->AddDataSet(theGGHNInel); 456 theInelasticProcess->RegisterMe(bertini); 457 theInelasticProcess->RegisterMe(theTheoModel); 458 pmanager->AddDiscreteProcess(theInelasticProcess); 459 } 460 else if (particleName == "anti_lambda") { 461 pmanager->AddDiscreteProcess(theElasticProcess); 462 G4HadronInelasticProcess* theInelasticProcess = 463 new G4HadronInelasticProcess("inelastic", G4AntiLambda::Definition()); 464 theInelasticProcess->AddDataSet(theGGHNInel); 465 theInelasticProcess->RegisterMe(antiBHighEnergyModel); 466 pmanager->AddDiscreteProcess(theInelasticProcess); 467 } 468 else if (particleName == "sigma+") { 469 pmanager->AddDiscreteProcess(theElasticProcess); 470 G4HadronInelasticProcess* theInelasticProcess = 471 new G4HadronInelasticProcess("inelastic", G4SigmaPlus::Definition()); 472 theInelasticProcess->AddDataSet(theGGHNInel); 473 theInelasticProcess->RegisterMe(bertini); 474 theInelasticProcess->RegisterMe(theTheoModel); 475 pmanager->AddDiscreteProcess(theInelasticProcess); 476 } 477 else if (particleName == "sigma-") { 478 pmanager->AddDiscreteProcess(theElasticProcess); 479 G4HadronInelasticProcess* theInelasticProcess = 480 new G4HadronInelasticProcess("inelastic", G4SigmaMinus::Definition()); 481 theInelasticProcess->AddDataSet(theGGHNInel); 482 theInelasticProcess->RegisterMe(bertini); 483 theInelasticProcess->RegisterMe(theTheoModel); 484 pmanager->AddDiscreteProcess(theInelasticProcess); 485 } 486 else if (particleName == "anti_sigma+") { 487 pmanager->AddDiscreteProcess(theElasticProcess); 488 G4HadronInelasticProcess* theInelasticProcess = 489 new G4HadronInelasticProcess("inelastic", G4AntiSigmaPlus::Definition()); 490 theInelasticProcess->AddDataSet(theGGHNInel); 491 theInelasticProcess->RegisterMe(antiBHighEnergyModel); 492 pmanager->AddDiscreteProcess(theInelasticProcess); 493 } 494 else if (particleName == "anti_sigma-") { 495 pmanager->AddDiscreteProcess(theElasticProcess); 496 G4HadronInelasticProcess* theInelasticProcess = 497 new G4HadronInelasticProcess("inelastic", G4AntiSigmaMinus::Definition()); 498 theInelasticProcess->AddDataSet(theGGHNInel); 499 theInelasticProcess->RegisterMe(antiBHighEnergyModel); 500 pmanager->AddDiscreteProcess(theInelasticProcess); 501 } 502 else if (particleName == "xi0") { 503 pmanager->AddDiscreteProcess(theElasticProcess); 504 G4HadronInelasticProcess* theInelasticProcess = 505 new G4HadronInelasticProcess("inelastic", G4XiZero::Definition()); 506 theInelasticProcess->AddDataSet(theGGHNInel); 507 theInelasticProcess->RegisterMe(bertini); 508 theInelasticProcess->RegisterMe(theTheoModel); 509 pmanager->AddDiscreteProcess(theInelasticProcess); 510 } 511 else if (particleName == "xi-") { 512 pmanager->AddDiscreteProcess(theElasticProcess); 513 G4HadronInelasticProcess* theInelasticProcess = 514 new G4HadronInelasticProcess("inelastic", G4XiMinus::Definition()); 515 theInelasticProcess->AddDataSet(theGGHNInel); 516 theInelasticProcess->RegisterMe(bertini); 517 theInelasticProcess->RegisterMe(theTheoModel); 518 pmanager->AddDiscreteProcess(theInelasticProcess); 519 } 520 else if (particleName == "anti_xi0") { 521 pmanager->AddDiscreteProcess(theElasticProcess); 522 G4HadronInelasticProcess* theInelasticProcess = 523 new G4HadronInelasticProcess("inelastic", G4AntiXiZero::Definition()); 524 theInelasticProcess->AddDataSet(theGGHNInel); 525 theInelasticProcess->RegisterMe(antiBHighEnergyModel); 526 pmanager->AddDiscreteProcess(theInelasticProcess); 527 } 528 else if (particleName == "anti_xi-") { 529 pmanager->AddDiscreteProcess(theElasticProcess); 530 G4HadronInelasticProcess* theInelasticProcess = 531 new G4HadronInelasticProcess("inelastic", G4AntiXiMinus::Definition()); 532 theInelasticProcess->AddDataSet(theGGHNInel); 533 theInelasticProcess->RegisterMe(antiBHighEnergyModel); 534 pmanager->AddDiscreteProcess(theInelasticProcess); 535 } 536 else if (particleName == "deuteron") { 537 pmanager->AddDiscreteProcess(theElasticProcess); 538 G4HadronInelasticProcess* theInelasticProcess = 539 new G4HadronInelasticProcess("inelastic", G4Deuteron::Definition()); 540 theInelasticProcess->RegisterMe(binaryCascade); 541 theInelasticProcess->RegisterMe(qmd); 542 theInelasticProcess->AddDataSet(ionXS); 543 pmanager->AddDiscreteProcess(theInelasticProcess); 544 } 545 else if (particleName == "triton") { 546 pmanager->AddDiscreteProcess(theElasticProcess); 547 G4HadronInelasticProcess* theInelasticProcess = 548 new G4HadronInelasticProcess("inelastic", G4Triton::Definition()); 549 theInelasticProcess->RegisterMe(binaryCascade); 550 theInelasticProcess->RegisterMe(qmd); 551 theInelasticProcess->AddDataSet(ionXS); 552 pmanager->AddDiscreteProcess(theInelasticProcess); 553 } 554 else if (particleName == "alpha") { 555 pmanager->AddDiscreteProcess(theElasticProcess); 556 G4HadronInelasticProcess* theInelasticProcess = 557 new G4HadronInelasticProcess("inelastic", G4Alpha::Definition()); 558 theInelasticProcess->RegisterMe(binaryCascade); 559 theInelasticProcess->RegisterMe(qmd); 560 theInelasticProcess->AddDataSet(ionXS); 561 pmanager->AddDiscreteProcess(theInelasticProcess); 562 } 563 else if (particleName == "omega-") { 564 pmanager->AddDiscreteProcess(theElasticProcess); 565 G4HadronInelasticProcess* theInelasticProcess = 566 new G4HadronInelasticProcess("inelastic", G4OmegaMinus::Definition()); 567 theInelasticProcess->AddDataSet(theGGHNInel); 568 theInelasticProcess->RegisterMe(bertini); 569 theInelasticProcess->RegisterMe(theTheoModel); 570 pmanager->AddDiscreteProcess(theInelasticProcess); 571 } 572 else if (particleName == "anti_omega-") { 573 pmanager->AddDiscreteProcess(theElasticProcess); 574 G4HadronInelasticProcess* theInelasticProcess = 575 new G4HadronInelasticProcess("inelastic", G4AntiOmegaMinus::Definition()); 576 theInelasticProcess->AddDataSet(theGGHNInel); 577 theInelasticProcess->RegisterMe(antiBHighEnergyModel); 578 pmanager->AddDiscreteProcess(theInelasticProcess); 579 } 580 } 581 } 582 583 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 584 585 void B03PhysicsList::ConstructLeptHad() 586 { 587 ; 588 } 589 590 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 591 592 #include "G4Decay.hh" 593 void B03PhysicsList::ConstructGeneral() 594 { 595 G4Decay* theDecayProcess = new G4Decay(); 596 auto particleIterator = GetParticleIterator(); 597 particleIterator->reset(); 598 while ((*particleIterator)()) { 599 G4ParticleDefinition* particle = particleIterator->value(); 600 G4ProcessManager* pmanager = particle->GetProcessManager(); 601 if (theDecayProcess->IsApplicable(*particle)) { 602 pmanager->AddProcess(theDecayProcess); 603 pmanager->SetProcessOrdering(theDecayProcess, idxPostStep); 604 pmanager->SetProcessOrdering(theDecayProcess, idxAtRest); 605 } 606 } 607 } 608 609 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 610 611 void B03PhysicsList::SetCuts() 612 { 613 if (verboseLevel > 0) { 614 G4cout << "B03PhysicsList::SetCuts:"; 615 G4cout << "CutLength : " << defaultCutValue / mm << " (mm)" << G4endl; 616 } 617 // "G4VUserPhysicsList::SetCutsWithDefault" method sets 618 // the default cut value for all particle types 619 SetCutsWithDefault(); 620 } 621 622 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 623 624 #include "G4ParallelWorldProcess.hh" 625 void B03PhysicsList::AddScoringProcess() 626 { 627 G4int npw = fParaWorldName.size(); 628 for (G4int i = 0; i < npw; i++) { 629 G4String procName = "ParaWorldProc_" + fParaWorldName[i]; 630 G4ParallelWorldProcess* theParallelWorldProcess = new G4ParallelWorldProcess(procName); 631 theParallelWorldProcess->SetParallelWorld(fParaWorldName[i]); 632 633 auto particleIterator = GetParticleIterator(); 634 particleIterator->reset(); 635 while ((*particleIterator)()) { 636 G4ParticleDefinition* particle = particleIterator->value(); 637 G4ProcessManager* pmanager = particle->GetProcessManager(); 638 pmanager->AddProcess(theParallelWorldProcess); 639 if (theParallelWorldProcess->IsAtRestRequired(particle)) { 640 pmanager->SetProcessOrdering(theParallelWorldProcess, idxAtRest, 9900); 641 } 642 pmanager->SetProcessOrderingToSecond(theParallelWorldProcess, idxAlongStep); 643 pmanager->SetProcessOrdering(theParallelWorldProcess, idxPostStep, 9900); 644 } 645 } 646 } 647 648 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 649 650 #include "G4IStore.hh" 651 #include "G4ImportanceProcess.hh" 652 void B03PhysicsList::AddBiasingProcess() 653 { 654 G4cout << " Preparing Importance Sampling with GhostWorld " << fBiasWorldName << G4endl; 655 656 G4IStore* iStore = G4IStore::GetInstance(fBiasWorldName); 657 G4GeometrySampler fGeomSampler(fBiasWorldName, "neutron"); 658 fGeomSampler.SetParallel(true); // parallelworld 659 // fGeomSampler.SetWorld(iStore->GetParallelWorldVolumePointer()); 660 // fGeomSampler->PrepareImportanceSampling(G4IStore:: 661 // GetInstance(fBiasWorldName), 0); 662 static G4bool first = true; 663 if (first) { 664 fGeomSampler.PrepareImportanceSampling(iStore, 0); 665 666 fGeomSampler.Configure(); 667 G4cout << " GeomSampler Configured!!! " << G4endl; 668 first = false; 669 } 670 671 #ifdef G4MULTITHREADED 672 if (!G4Threading::IsMasterThread()) fGeomSampler.AddProcess(); 673 #else 674 G4cout << " Running in singlethreaded mode!!! " << G4endl; 675 #endif 676 } 677 678 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 679