File indexing completed on 2025-12-16 09:30:08
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032 #include "B03PhysicsList.hh"
0033
0034 #include "G4BaryonConstructor.hh"
0035 #include "G4BosonConstructor.hh"
0036 #include "G4HadronicParameters.hh"
0037 #include "G4IonConstructor.hh"
0038 #include "G4LeptonConstructor.hh"
0039 #include "G4Material.hh"
0040 #include "G4MaterialTable.hh"
0041 #include "G4MesonConstructor.hh"
0042 #include "G4ParticleDefinition.hh"
0043 #include "G4ParticleTable.hh"
0044 #include "G4ParticleTypes.hh"
0045 #include "G4ParticleWithCuts.hh"
0046 #include "G4ProcessManager.hh"
0047 #include "G4ProcessVector.hh"
0048 #include "G4ShortLivedConstructor.hh"
0049 #include "G4SystemOfUnits.hh"
0050 #include "globals.hh"
0051
0052 #include <iomanip>
0053
0054
0055
0056 B03PhysicsList::B03PhysicsList(G4String parallelname)
0057 : G4VUserPhysicsList(), fBiasWorldName(parallelname)
0058 {
0059 fParaWorldName.clear();
0060 SetVerboseLevel(1);
0061 }
0062
0063
0064
0065 B03PhysicsList::~B03PhysicsList()
0066 {
0067 fParaWorldName.clear();
0068 }
0069
0070
0071
0072 void B03PhysicsList::ConstructParticle()
0073 {
0074
0075
0076
0077
0078
0079 ConstructAllBosons();
0080 ConstructAllLeptons();
0081 ConstructAllMesons();
0082 ConstructAllBaryons();
0083 ConstructAllIons();
0084 ConstructAllShortLiveds();
0085 }
0086
0087
0088
0089 void B03PhysicsList::ConstructAllBosons()
0090 {
0091
0092 G4BosonConstructor pConstructor;
0093 pConstructor.ConstructParticle();
0094 }
0095
0096
0097
0098 void B03PhysicsList::ConstructAllLeptons()
0099 {
0100
0101 G4LeptonConstructor pConstructor;
0102 pConstructor.ConstructParticle();
0103 }
0104
0105
0106
0107 void B03PhysicsList::ConstructAllMesons()
0108 {
0109
0110 G4MesonConstructor pConstructor;
0111 pConstructor.ConstructParticle();
0112 }
0113
0114
0115
0116 void B03PhysicsList::ConstructAllBaryons()
0117 {
0118
0119 G4BaryonConstructor pConstructor;
0120 pConstructor.ConstructParticle();
0121 }
0122
0123
0124
0125 void B03PhysicsList::ConstructAllIons()
0126 {
0127
0128 G4IonConstructor pConstructor;
0129 pConstructor.ConstructParticle();
0130 }
0131
0132
0133
0134 void B03PhysicsList::ConstructAllShortLiveds()
0135 {
0136
0137 G4ShortLivedConstructor pConstructor;
0138 pConstructor.ConstructParticle();
0139 }
0140
0141
0142
0143 void B03PhysicsList::ConstructProcess()
0144 {
0145 AddTransportation();
0146 AddScoringProcess();
0147 AddBiasingProcess();
0148 ConstructEM();
0149 ConstructLeptHad();
0150 ConstructHad();
0151 ConstructGeneral();
0152 }
0153
0154
0155
0156 #include "G4ComptonScattering.hh"
0157 #include "G4GammaConversion.hh"
0158 #include "G4MuBremsstrahlung.hh"
0159 #include "G4MuIonisation.hh"
0160 #include "G4MuMultipleScattering.hh"
0161 #include "G4MuPairProduction.hh"
0162 #include "G4PhotoElectricEffect.hh"
0163 #include "G4eBremsstrahlung.hh"
0164 #include "G4eIonisation.hh"
0165 #include "G4eMultipleScattering.hh"
0166 #include "G4eplusAnnihilation.hh"
0167 #include "G4hIonisation.hh"
0168 #include "G4hMultipleScattering.hh"
0169
0170 void B03PhysicsList::ConstructEM()
0171 {
0172 auto particleIterator = GetParticleIterator();
0173 particleIterator->reset();
0174 while ((*particleIterator)()) {
0175 G4ParticleDefinition* particle = particleIterator->value();
0176 G4ProcessManager* pmanager = particle->GetProcessManager();
0177 G4String particleName = particle->GetParticleName();
0178
0179 if (particleName == "gamma") {
0180
0181
0182 pmanager->AddDiscreteProcess(new G4GammaConversion());
0183 pmanager->AddDiscreteProcess(new G4ComptonScattering());
0184 pmanager->AddDiscreteProcess(new G4PhotoElectricEffect());
0185 }
0186 else if (particleName == "e-") {
0187
0188
0189 pmanager->AddProcess(new G4eMultipleScattering(), -1, 1, 1);
0190 pmanager->AddProcess(new G4eIonisation(), -1, 2, 2);
0191 pmanager->AddProcess(new G4eBremsstrahlung(), -1, -1, 3);
0192 }
0193 else if (particleName == "e+") {
0194
0195
0196 pmanager->AddProcess(new G4eMultipleScattering(), -1, 1, 1);
0197
0198 pmanager->AddProcess(new G4eIonisation(), -1, 2, 2);
0199 pmanager->AddProcess(new G4eBremsstrahlung(), -1, -1, 3);
0200 pmanager->AddProcess(new G4eplusAnnihilation(), 0, -1, 4);
0201 }
0202 else if (particleName == "mu+" || particleName == "mu-") {
0203
0204
0205 pmanager->AddProcess(new G4MuMultipleScattering(), -1, 1, 1);
0206 pmanager->AddProcess(new G4MuIonisation(), -1, 2, 2);
0207 pmanager->AddProcess(new G4MuBremsstrahlung(), -1, -1, 3);
0208 pmanager->AddProcess(new G4MuPairProduction(), -1, -1, 4);
0209 }
0210 else if (particleName == "GenericIon") {
0211 pmanager->AddProcess(new G4hMultipleScattering(), -1, 1, 1);
0212 pmanager->AddProcess(new G4hIonisation(), -1, 2, 2);
0213 }
0214 else {
0215 if ((particle->GetPDGCharge() != 0.0) && (particle->GetParticleName() != "chargedgeantino")
0216 && (!particle->IsShortLived()))
0217 {
0218
0219 pmanager->AddProcess(new G4hMultipleScattering(), -1, 1, 1);
0220 pmanager->AddProcess(new G4hIonisation(), -1, 2, 2);
0221 }
0222 }
0223 }
0224 }
0225
0226
0227
0228
0229
0230 #include "G4HadronElasticProcess.hh"
0231 #include "G4HadronInelasticProcess.hh"
0232 #include "G4NeutronCaptureProcess.hh"
0233 #include "G4NeutronFissionProcess.hh"
0234
0235
0236
0237 #include "G4HadronElastic.hh"
0238 #include "G4LFission.hh"
0239 #include "G4NeutronRadCapture.hh"
0240
0241
0242 #include "G4BinaryLightIonReaction.hh"
0243 #include "G4CascadeInterface.hh"
0244 #include "G4CompetitiveFission.hh"
0245 #include "G4ExcitationHandler.hh"
0246 #include "G4ExcitedStringDecay.hh"
0247 #include "G4FTFModel.hh"
0248 #include "G4Fancy3DNucleus.hh"
0249 #include "G4GeneratorPrecompoundInterface.hh"
0250 #include "G4LundStringFragmentation.hh"
0251 #include "G4PreCompoundModel.hh"
0252 #include "G4QGSMFragmentation.hh"
0253 #include "G4QMDReaction.hh"
0254 #include "G4StringModel.hh"
0255 #include "G4TheoFSGenerator.hh"
0256
0257
0258 #include "G4ComponentGGHadronNucleusXsc.hh"
0259 #include "G4ComponentGGNuclNuclXsc.hh"
0260 #include "G4CrossSectionElastic.hh"
0261 #include "G4CrossSectionInelastic.hh"
0262 #include "G4NeutronInelasticXS.hh"
0263
0264
0265
0266
0267
0268
0269
0270
0271
0272 void B03PhysicsList::ConstructHad()
0273 {
0274
0275 G4TheoFSGenerator* theTheoModel = new G4TheoFSGenerator;
0276 G4TheoFSGenerator* antiBHighEnergyModel = new G4TheoFSGenerator;
0277
0278
0279 G4ExcitationHandler* theHandler = new G4ExcitationHandler;
0280 theHandler->SetMinEForMultiFrag(3 * MeV);
0281
0282
0283 G4PreCompoundModel* thePreEquilib = new G4PreCompoundModel(theHandler);
0284
0285
0286 G4GeneratorPrecompoundInterface* theCascade = new G4GeneratorPrecompoundInterface;
0287 theCascade->SetDeExcitation(thePreEquilib);
0288
0289
0290 G4CascadeInterface* bertini = new G4CascadeInterface;
0291 bertini->SetMaxEnergy(22 * MeV);
0292
0293
0294 G4VPartonStringModel* theStringModel;
0295 theStringModel = new G4FTFModel;
0296 theTheoModel->SetTransport(theCascade);
0297 theTheoModel->SetHighEnergyGenerator(theStringModel);
0298 theTheoModel->SetMinEnergy(19 * GeV);
0299 theTheoModel->SetMaxEnergy(G4HadronicParameters::Instance()->GetMaxEnergy());
0300
0301 G4VLongitudinalStringDecay* theFragmentation = new G4QGSMFragmentation;
0302 G4ExcitedStringDecay* theStringDecay = new G4ExcitedStringDecay(theFragmentation);
0303 theStringModel->SetFragmentationModel(theStringDecay);
0304
0305
0306 antiBHighEnergyModel = new G4TheoFSGenerator("ANTI-FTFP");
0307 G4FTFModel* antiBStringModel = new G4FTFModel;
0308 G4ExcitedStringDecay* stringDecay = new G4ExcitedStringDecay(new G4LundStringFragmentation);
0309 antiBStringModel->SetFragmentationModel(stringDecay);
0310
0311 G4GeneratorPrecompoundInterface* antiBCascade = new G4GeneratorPrecompoundInterface;
0312 G4PreCompoundModel* preEquilib = new G4PreCompoundModel(new G4ExcitationHandler);
0313 antiBCascade->SetDeExcitation(preEquilib);
0314
0315 antiBHighEnergyModel->SetTransport(antiBCascade);
0316 antiBHighEnergyModel->SetHighEnergyGenerator(antiBStringModel);
0317 antiBHighEnergyModel->SetMinEnergy(0.0);
0318 antiBHighEnergyModel->SetMaxEnergy(20 * TeV);
0319
0320
0321 G4BinaryLightIonReaction* binaryCascade = new G4BinaryLightIonReaction;
0322 binaryCascade->SetMinEnergy(0.0);
0323 binaryCascade->SetMaxEnergy(110 * MeV);
0324
0325 G4QMDReaction* qmd = new G4QMDReaction;
0326 qmd->SetMinEnergy(100 * MeV);
0327 qmd->SetMaxEnergy(10 * GeV);
0328
0329 G4VCrossSectionDataSet* ionXS = new G4CrossSectionInelastic(new G4ComponentGGNuclNuclXsc);
0330
0331 G4ComponentGGHadronNucleusXsc* ggHNXsec = new G4ComponentGGHadronNucleusXsc();
0332 G4VCrossSectionDataSet* theGGHNEl = new G4CrossSectionElastic(ggHNXsec);
0333 G4VCrossSectionDataSet* theGGHNInel = new G4CrossSectionInelastic(ggHNXsec);
0334
0335
0336 G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
0337 theElasticProcess->AddDataSet(theGGHNEl);
0338 G4HadronElastic* theElasticModel = new G4HadronElastic;
0339 theElasticProcess->RegisterMe(theElasticModel);
0340
0341 auto particleIterator = GetParticleIterator();
0342 particleIterator->reset();
0343 while ((*particleIterator)()) {
0344 G4ParticleDefinition* particle = particleIterator->value();
0345 G4ProcessManager* pmanager = particle->GetProcessManager();
0346 G4String particleName = particle->GetParticleName();
0347
0348 if (particleName == "pi+") {
0349 pmanager->AddDiscreteProcess(theElasticProcess);
0350 G4HadronInelasticProcess* theInelasticProcess =
0351 new G4HadronInelasticProcess("inelastic", G4PionPlus::Definition());
0352 theInelasticProcess->AddDataSet(theGGHNInel);
0353 theInelasticProcess->RegisterMe(bertini);
0354 theInelasticProcess->RegisterMe(theTheoModel);
0355 pmanager->AddDiscreteProcess(theInelasticProcess);
0356 }
0357 else if (particleName == "pi-") {
0358 pmanager->AddDiscreteProcess(theElasticProcess);
0359 G4HadronInelasticProcess* theInelasticProcess =
0360 new G4HadronInelasticProcess("inelastic", G4PionMinus::Definition());
0361 theInelasticProcess->AddDataSet(theGGHNInel);
0362 theInelasticProcess->RegisterMe(bertini);
0363 theInelasticProcess->RegisterMe(theTheoModel);
0364 pmanager->AddDiscreteProcess(theInelasticProcess);
0365 }
0366 else if (particleName == "kaon+") {
0367 pmanager->AddDiscreteProcess(theElasticProcess);
0368 G4HadronInelasticProcess* theInelasticProcess =
0369 new G4HadronInelasticProcess("inelastic", G4KaonPlus::Definition());
0370 theInelasticProcess->AddDataSet(theGGHNInel);
0371 theInelasticProcess->RegisterMe(bertini);
0372 theInelasticProcess->RegisterMe(theTheoModel);
0373 pmanager->AddDiscreteProcess(theInelasticProcess);
0374 }
0375 else if (particleName == "kaon0S") {
0376 pmanager->AddDiscreteProcess(theElasticProcess);
0377 G4HadronInelasticProcess* theInelasticProcess =
0378 new G4HadronInelasticProcess("inelastic", G4KaonZeroShort::Definition());
0379 theInelasticProcess->AddDataSet(theGGHNInel);
0380 theInelasticProcess->RegisterMe(bertini);
0381 theInelasticProcess->RegisterMe(theTheoModel);
0382 pmanager->AddDiscreteProcess(theInelasticProcess);
0383 }
0384 else if (particleName == "kaon0L") {
0385 pmanager->AddDiscreteProcess(theElasticProcess);
0386 G4HadronInelasticProcess* theInelasticProcess =
0387 new G4HadronInelasticProcess("inelastic", G4KaonZeroLong::Definition());
0388 theInelasticProcess->AddDataSet(theGGHNInel);
0389 theInelasticProcess->RegisterMe(bertini);
0390 theInelasticProcess->RegisterMe(theTheoModel);
0391 pmanager->AddDiscreteProcess(theInelasticProcess);
0392 }
0393 else if (particleName == "kaon-") {
0394 pmanager->AddDiscreteProcess(theElasticProcess);
0395 G4HadronInelasticProcess* theInelasticProcess =
0396 new G4HadronInelasticProcess("inelastic", G4KaonMinus::Definition());
0397 theInelasticProcess->AddDataSet(theGGHNInel);
0398 theInelasticProcess->RegisterMe(bertini);
0399 theInelasticProcess->RegisterMe(theTheoModel);
0400 pmanager->AddDiscreteProcess(theInelasticProcess);
0401 }
0402 else if (particleName == "proton") {
0403 pmanager->AddDiscreteProcess(theElasticProcess);
0404 G4HadronInelasticProcess* theInelasticProcess =
0405 new G4HadronInelasticProcess("inelastic", G4Proton::Definition());
0406 theInelasticProcess->AddDataSet(theGGHNInel);
0407 theInelasticProcess->RegisterMe(bertini);
0408 theInelasticProcess->RegisterMe(theTheoModel);
0409 pmanager->AddDiscreteProcess(theInelasticProcess);
0410 }
0411 else if (particleName == "anti_proton") {
0412 pmanager->AddDiscreteProcess(theElasticProcess);
0413 G4HadronInelasticProcess* theInelasticProcess =
0414 new G4HadronInelasticProcess("inelastic", G4AntiProton::Definition());
0415 theInelasticProcess->AddDataSet(theGGHNInel);
0416 theInelasticProcess->RegisterMe(antiBHighEnergyModel);
0417 pmanager->AddDiscreteProcess(theInelasticProcess);
0418 }
0419 else if (particleName == "neutron") {
0420
0421 pmanager->AddDiscreteProcess(theElasticProcess);
0422
0423
0424 G4HadronInelasticProcess* theInelasticProcess =
0425 new G4HadronInelasticProcess("inelastic", G4Neutron::Definition());
0426 theInelasticProcess->AddDataSet(new G4NeutronInelasticXS());
0427 theInelasticProcess->RegisterMe(bertini);
0428 theInelasticProcess->RegisterMe(theTheoModel);
0429 pmanager->AddDiscreteProcess(theInelasticProcess);
0430
0431
0432 G4NeutronFissionProcess* theFissionProcess = new G4NeutronFissionProcess;
0433 G4LFission* theFissionModel = new G4LFission;
0434 theFissionProcess->RegisterMe(theFissionModel);
0435 pmanager->AddDiscreteProcess(theFissionProcess);
0436
0437
0438 G4NeutronCaptureProcess* theCaptureProcess = new G4NeutronCaptureProcess;
0439 G4NeutronRadCapture* theCaptureModel = new G4NeutronRadCapture;
0440 theCaptureProcess->RegisterMe(theCaptureModel);
0441 pmanager->AddDiscreteProcess(theCaptureProcess);
0442 }
0443 else if (particleName == "anti_neutron") {
0444 pmanager->AddDiscreteProcess(theElasticProcess);
0445 G4HadronInelasticProcess* theInelasticProcess =
0446 new G4HadronInelasticProcess("inelastic", G4AntiNeutron::Definition());
0447 theInelasticProcess->AddDataSet(theGGHNInel);
0448 theInelasticProcess->RegisterMe(antiBHighEnergyModel);
0449 pmanager->AddDiscreteProcess(theInelasticProcess);
0450 }
0451 else if (particleName == "lambda") {
0452 pmanager->AddDiscreteProcess(theElasticProcess);
0453 G4HadronInelasticProcess* theInelasticProcess =
0454 new G4HadronInelasticProcess("inelastic", G4Lambda::Definition());
0455 theInelasticProcess->AddDataSet(theGGHNInel);
0456 theInelasticProcess->RegisterMe(bertini);
0457 theInelasticProcess->RegisterMe(theTheoModel);
0458 pmanager->AddDiscreteProcess(theInelasticProcess);
0459 }
0460 else if (particleName == "anti_lambda") {
0461 pmanager->AddDiscreteProcess(theElasticProcess);
0462 G4HadronInelasticProcess* theInelasticProcess =
0463 new G4HadronInelasticProcess("inelastic", G4AntiLambda::Definition());
0464 theInelasticProcess->AddDataSet(theGGHNInel);
0465 theInelasticProcess->RegisterMe(antiBHighEnergyModel);
0466 pmanager->AddDiscreteProcess(theInelasticProcess);
0467 }
0468 else if (particleName == "sigma+") {
0469 pmanager->AddDiscreteProcess(theElasticProcess);
0470 G4HadronInelasticProcess* theInelasticProcess =
0471 new G4HadronInelasticProcess("inelastic", G4SigmaPlus::Definition());
0472 theInelasticProcess->AddDataSet(theGGHNInel);
0473 theInelasticProcess->RegisterMe(bertini);
0474 theInelasticProcess->RegisterMe(theTheoModel);
0475 pmanager->AddDiscreteProcess(theInelasticProcess);
0476 }
0477 else if (particleName == "sigma-") {
0478 pmanager->AddDiscreteProcess(theElasticProcess);
0479 G4HadronInelasticProcess* theInelasticProcess =
0480 new G4HadronInelasticProcess("inelastic", G4SigmaMinus::Definition());
0481 theInelasticProcess->AddDataSet(theGGHNInel);
0482 theInelasticProcess->RegisterMe(bertini);
0483 theInelasticProcess->RegisterMe(theTheoModel);
0484 pmanager->AddDiscreteProcess(theInelasticProcess);
0485 }
0486 else if (particleName == "anti_sigma+") {
0487 pmanager->AddDiscreteProcess(theElasticProcess);
0488 G4HadronInelasticProcess* theInelasticProcess =
0489 new G4HadronInelasticProcess("inelastic", G4AntiSigmaPlus::Definition());
0490 theInelasticProcess->AddDataSet(theGGHNInel);
0491 theInelasticProcess->RegisterMe(antiBHighEnergyModel);
0492 pmanager->AddDiscreteProcess(theInelasticProcess);
0493 }
0494 else if (particleName == "anti_sigma-") {
0495 pmanager->AddDiscreteProcess(theElasticProcess);
0496 G4HadronInelasticProcess* theInelasticProcess =
0497 new G4HadronInelasticProcess("inelastic", G4AntiSigmaMinus::Definition());
0498 theInelasticProcess->AddDataSet(theGGHNInel);
0499 theInelasticProcess->RegisterMe(antiBHighEnergyModel);
0500 pmanager->AddDiscreteProcess(theInelasticProcess);
0501 }
0502 else if (particleName == "xi0") {
0503 pmanager->AddDiscreteProcess(theElasticProcess);
0504 G4HadronInelasticProcess* theInelasticProcess =
0505 new G4HadronInelasticProcess("inelastic", G4XiZero::Definition());
0506 theInelasticProcess->AddDataSet(theGGHNInel);
0507 theInelasticProcess->RegisterMe(bertini);
0508 theInelasticProcess->RegisterMe(theTheoModel);
0509 pmanager->AddDiscreteProcess(theInelasticProcess);
0510 }
0511 else if (particleName == "xi-") {
0512 pmanager->AddDiscreteProcess(theElasticProcess);
0513 G4HadronInelasticProcess* theInelasticProcess =
0514 new G4HadronInelasticProcess("inelastic", G4XiMinus::Definition());
0515 theInelasticProcess->AddDataSet(theGGHNInel);
0516 theInelasticProcess->RegisterMe(bertini);
0517 theInelasticProcess->RegisterMe(theTheoModel);
0518 pmanager->AddDiscreteProcess(theInelasticProcess);
0519 }
0520 else if (particleName == "anti_xi0") {
0521 pmanager->AddDiscreteProcess(theElasticProcess);
0522 G4HadronInelasticProcess* theInelasticProcess =
0523 new G4HadronInelasticProcess("inelastic", G4AntiXiZero::Definition());
0524 theInelasticProcess->AddDataSet(theGGHNInel);
0525 theInelasticProcess->RegisterMe(antiBHighEnergyModel);
0526 pmanager->AddDiscreteProcess(theInelasticProcess);
0527 }
0528 else if (particleName == "anti_xi-") {
0529 pmanager->AddDiscreteProcess(theElasticProcess);
0530 G4HadronInelasticProcess* theInelasticProcess =
0531 new G4HadronInelasticProcess("inelastic", G4AntiXiMinus::Definition());
0532 theInelasticProcess->AddDataSet(theGGHNInel);
0533 theInelasticProcess->RegisterMe(antiBHighEnergyModel);
0534 pmanager->AddDiscreteProcess(theInelasticProcess);
0535 }
0536 else if (particleName == "deuteron") {
0537 pmanager->AddDiscreteProcess(theElasticProcess);
0538 G4HadronInelasticProcess* theInelasticProcess =
0539 new G4HadronInelasticProcess("inelastic", G4Deuteron::Definition());
0540 theInelasticProcess->RegisterMe(binaryCascade);
0541 theInelasticProcess->RegisterMe(qmd);
0542 theInelasticProcess->AddDataSet(ionXS);
0543 pmanager->AddDiscreteProcess(theInelasticProcess);
0544 }
0545 else if (particleName == "triton") {
0546 pmanager->AddDiscreteProcess(theElasticProcess);
0547 G4HadronInelasticProcess* theInelasticProcess =
0548 new G4HadronInelasticProcess("inelastic", G4Triton::Definition());
0549 theInelasticProcess->RegisterMe(binaryCascade);
0550 theInelasticProcess->RegisterMe(qmd);
0551 theInelasticProcess->AddDataSet(ionXS);
0552 pmanager->AddDiscreteProcess(theInelasticProcess);
0553 }
0554 else if (particleName == "alpha") {
0555 pmanager->AddDiscreteProcess(theElasticProcess);
0556 G4HadronInelasticProcess* theInelasticProcess =
0557 new G4HadronInelasticProcess("inelastic", G4Alpha::Definition());
0558 theInelasticProcess->RegisterMe(binaryCascade);
0559 theInelasticProcess->RegisterMe(qmd);
0560 theInelasticProcess->AddDataSet(ionXS);
0561 pmanager->AddDiscreteProcess(theInelasticProcess);
0562 }
0563 else if (particleName == "omega-") {
0564 pmanager->AddDiscreteProcess(theElasticProcess);
0565 G4HadronInelasticProcess* theInelasticProcess =
0566 new G4HadronInelasticProcess("inelastic", G4OmegaMinus::Definition());
0567 theInelasticProcess->AddDataSet(theGGHNInel);
0568 theInelasticProcess->RegisterMe(bertini);
0569 theInelasticProcess->RegisterMe(theTheoModel);
0570 pmanager->AddDiscreteProcess(theInelasticProcess);
0571 }
0572 else if (particleName == "anti_omega-") {
0573 pmanager->AddDiscreteProcess(theElasticProcess);
0574 G4HadronInelasticProcess* theInelasticProcess =
0575 new G4HadronInelasticProcess("inelastic", G4AntiOmegaMinus::Definition());
0576 theInelasticProcess->AddDataSet(theGGHNInel);
0577 theInelasticProcess->RegisterMe(antiBHighEnergyModel);
0578 pmanager->AddDiscreteProcess(theInelasticProcess);
0579 }
0580 }
0581 }
0582
0583
0584
0585 void B03PhysicsList::ConstructLeptHad()
0586 {
0587 ;
0588 }
0589
0590
0591
0592 #include "G4Decay.hh"
0593 void B03PhysicsList::ConstructGeneral()
0594 {
0595 G4Decay* theDecayProcess = new G4Decay();
0596 auto particleIterator = GetParticleIterator();
0597 particleIterator->reset();
0598 while ((*particleIterator)()) {
0599 G4ParticleDefinition* particle = particleIterator->value();
0600 G4ProcessManager* pmanager = particle->GetProcessManager();
0601 if (theDecayProcess->IsApplicable(*particle)) {
0602 pmanager->AddProcess(theDecayProcess);
0603 pmanager->SetProcessOrdering(theDecayProcess, idxPostStep);
0604 pmanager->SetProcessOrdering(theDecayProcess, idxAtRest);
0605 }
0606 }
0607 }
0608
0609
0610
0611 void B03PhysicsList::SetCuts()
0612 {
0613 if (verboseLevel > 0) {
0614 G4cout << "B03PhysicsList::SetCuts:";
0615 G4cout << "CutLength : " << defaultCutValue / mm << " (mm)" << G4endl;
0616 }
0617
0618
0619 SetCutsWithDefault();
0620 }
0621
0622
0623
0624 #include "G4ParallelWorldProcess.hh"
0625 void B03PhysicsList::AddScoringProcess()
0626 {
0627 G4int npw = fParaWorldName.size();
0628 for (G4int i = 0; i < npw; i++) {
0629 G4String procName = "ParaWorldProc_" + fParaWorldName[i];
0630 G4ParallelWorldProcess* theParallelWorldProcess = new G4ParallelWorldProcess(procName);
0631 theParallelWorldProcess->SetParallelWorld(fParaWorldName[i]);
0632
0633 auto particleIterator = GetParticleIterator();
0634 particleIterator->reset();
0635 while ((*particleIterator)()) {
0636 G4ParticleDefinition* particle = particleIterator->value();
0637 G4ProcessManager* pmanager = particle->GetProcessManager();
0638 pmanager->AddProcess(theParallelWorldProcess);
0639 if (theParallelWorldProcess->IsAtRestRequired(particle)) {
0640 pmanager->SetProcessOrdering(theParallelWorldProcess, idxAtRest, 9900);
0641 }
0642 pmanager->SetProcessOrderingToSecond(theParallelWorldProcess, idxAlongStep);
0643 pmanager->SetProcessOrdering(theParallelWorldProcess, idxPostStep, 9900);
0644 }
0645 }
0646 }
0647
0648
0649
0650 #include "G4IStore.hh"
0651 #include "G4ImportanceProcess.hh"
0652 void B03PhysicsList::AddBiasingProcess()
0653 {
0654 G4cout << " Preparing Importance Sampling with GhostWorld " << fBiasWorldName << G4endl;
0655
0656 G4IStore* iStore = G4IStore::GetInstance(fBiasWorldName);
0657 G4GeometrySampler fGeomSampler(fBiasWorldName, "neutron");
0658 fGeomSampler.SetParallel(true);
0659
0660
0661
0662 static G4bool first = true;
0663 if (first) {
0664 fGeomSampler.PrepareImportanceSampling(iStore, 0);
0665
0666 fGeomSampler.Configure();
0667 G4cout << " GeomSampler Configured!!! " << G4endl;
0668 first = false;
0669 }
0670
0671 #ifdef G4MULTITHREADED
0672 if (!G4Threading::IsMasterThread()) fGeomSampler.AddProcess();
0673 #else
0674 G4cout << " Running in singlethreaded mode!!! " << G4endl;
0675 #endif
0676 }
0677
0678