Warning, file /geant4/examples/extended/hadronic/FlukaCern/ProcessLevel/FinalState/src/HadronicGenerator.cc was not indexed
or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).
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
0033
0034
0035
0036
0037
0038 #include "HadronicGenerator.hh"
0039
0040 #include "G4AblaInterface.hh"
0041 #include "G4Alpha.hh"
0042 #include "G4AntiAlpha.hh"
0043 #include "G4AntiBMesonZero.hh"
0044 #include "G4AntiBsMesonZero.hh"
0045 #include "G4AntiDMesonZero.hh"
0046 #include "G4AntiDeuteron.hh"
0047 #include "G4AntiHe3.hh"
0048 #include "G4AntiLambda.hh"
0049 #include "G4AntiLambdab.hh"
0050 #include "G4AntiLambdacPlus.hh"
0051 #include "G4AntiNeutron.hh"
0052 #include "G4AntiOmegaMinus.hh"
0053 #include "G4AntiOmegabMinus.hh"
0054 #include "G4AntiOmegacZero.hh"
0055 #include "G4AntiProton.hh"
0056 #include "G4AntiSigmaMinus.hh"
0057 #include "G4AntiSigmaPlus.hh"
0058 #include "G4AntiSigmaZero.hh"
0059 #include "G4AntiTriton.hh"
0060 #include "G4AntiXiMinus.hh"
0061 #include "G4AntiXiZero.hh"
0062 #include "G4AntiXibMinus.hh"
0063 #include "G4AntiXibZero.hh"
0064 #include "G4AntiXicPlus.hh"
0065 #include "G4AntiXicZero.hh"
0066 #include "G4BGGNucleonInelasticXS.hh"
0067 #include "G4BGGPionInelasticXS.hh"
0068 #include "G4BMesonMinus.hh"
0069 #include "G4BMesonPlus.hh"
0070 #include "G4BMesonZero.hh"
0071 #include "G4BcMesonMinus.hh"
0072 #include "G4BcMesonPlus.hh"
0073 #include "G4BinaryCascade.hh"
0074 #include "G4BinaryLightIonReaction.hh"
0075 #include "G4Box.hh"
0076 #include "G4BsMesonZero.hh"
0077 #include "G4CascadeInterface.hh"
0078 #include "G4ChipsHyperonInelasticXS.hh"
0079 #include "G4ComponentAntiNuclNuclearXS.hh"
0080 #include "G4ComponentGGHadronNucleusXsc.hh"
0081 #include "G4ComponentGGNuclNuclXsc.hh"
0082 #include "G4CrossSectionInelastic.hh"
0083 #include "G4DMesonMinus.hh"
0084 #include "G4DMesonPlus.hh"
0085 #include "G4DMesonZero.hh"
0086 #include "G4DecayPhysics.hh"
0087 #include "G4Deuteron.hh"
0088 #include "G4DsMesonMinus.hh"
0089 #include "G4DsMesonPlus.hh"
0090 #include "G4DynamicParticle.hh"
0091 #include "G4ExcitationHandler.hh"
0092 #include "G4ExcitedStringDecay.hh"
0093 #include "G4FTFModel.hh"
0094 #include "G4GeneratorPrecompoundInterface.hh"
0095 #include "G4GenericIon.hh"
0096 #include "G4HadronInelasticProcess.hh"
0097 #include "G4HadronicParameters.hh"
0098 #include "G4He3.hh"
0099 #include "G4INCLXXInterface.hh"
0100 #include "G4IonTable.hh"
0101 #include "G4KaonMinus.hh"
0102 #include "G4KaonPlus.hh"
0103 #include "G4KaonZeroLong.hh"
0104 #include "G4KaonZeroShort.hh"
0105 #include "G4Lambda.hh"
0106 #include "G4Lambdab.hh"
0107 #include "G4LambdacPlus.hh"
0108 #include "G4LundStringFragmentation.hh"
0109 #include "G4Material.hh"
0110 #include "G4Neutron.hh"
0111 #include "G4NeutronInelasticXS.hh"
0112 #include "G4OmegaMinus.hh"
0113 #include "G4OmegabMinus.hh"
0114 #include "G4OmegacZero.hh"
0115 #include "G4PVPlacement.hh"
0116 #include "G4ParticleTable.hh"
0117 #include "G4PhysicalConstants.hh"
0118 #include "G4PionMinus.hh"
0119 #include "G4PionPlus.hh"
0120 #include "G4PreCompoundModel.hh"
0121 #include "G4ProcessManager.hh"
0122 #include "G4Proton.hh"
0123 #include "G4QGSMFragmentation.hh"
0124 #include "G4QGSModel.hh"
0125 #include "G4QGSParticipants.hh"
0126 #include "G4QuasiElasticChannel.hh"
0127 #include "G4SigmaMinus.hh"
0128 #include "G4SigmaPlus.hh"
0129 #include "G4SigmaZero.hh"
0130 #include "G4StateManager.hh"
0131 #include "G4Step.hh"
0132 #include "G4SystemOfUnits.hh"
0133 #include "G4TheoFSGenerator.hh"
0134 #include "G4TouchableHistory.hh"
0135 #include "G4TransportationManager.hh"
0136 #include "G4Triton.hh"
0137 #include "G4UnitsTable.hh"
0138 #include "G4VCrossSectionDataSet.hh"
0139 #include "G4VParticleChange.hh"
0140 #include "G4XiMinus.hh"
0141 #include "G4XiZero.hh"
0142 #include "G4XibMinus.hh"
0143 #include "G4XibZero.hh"
0144 #include "G4XicPlus.hh"
0145 #include "G4XicZero.hh"
0146 #include "G4ios.hh"
0147 #include "globals.hh"
0148
0149 #include <iomanip>
0150
0151
0152 #ifdef G4_USE_FLUKA
0153 # include "FLUKAInelasticScatteringXS.hh"
0154 # include "FLUKANuclearInelasticModel.hh"
0155 # include "FLUKAParticleTable.hh"
0156 # include "fluka_interface.hh"
0157 #endif
0158
0159
0160
0161 HadronicGenerator::HadronicGenerator(const G4String physicsCase)
0162 : fPhysicsCase(physicsCase),
0163 fPhysicsCaseIsSupported(false),
0164 fLastHadronicProcess(nullptr),
0165 fPartTable(nullptr)
0166 {
0167
0168
0169
0170
0171
0172
0173
0174
0175
0176
0177
0178
0179
0180
0181
0182
0183
0184 G4GenericIon* gion = G4GenericIon::Definition();
0185 gion->SetProcessManager(new G4ProcessManager(gion));
0186 G4DecayPhysics* decays = new G4DecayPhysics;
0187 decays->ConstructParticle();
0188 fPartTable = G4ParticleTable::GetParticleTable();
0189 fPartTable->SetReadiness();
0190 G4IonTable* ions = fPartTable->GetIonTable();
0191 ions->CreateAllIon();
0192 ions->CreateAllIsomer();
0193
0194
0195 if (fPhysicsCase == "CFLUKAHI") {
0196 #ifdef G4_USE_FLUKA
0197
0198
0199
0200
0201
0202
0203 const G4bool activateCoalescence = true;
0204 const G4bool activateHeavyFragmentsEvaporation = true;
0205 fluka_interface::initialize(activateCoalescence, activateHeavyFragmentsEvaporation);
0206
0207
0208 fluka_particle_table::initialize();
0209 #else
0210 G4Exception("HadronicGenerator.cc", "Wrong compilation mode.", FatalException,
0211 "Requested G4_HP_CernFLUKAHadronInelastic physics list.\n"
0212 "This requires COMPILATION in FLUKA mode.\n"
0213 "Please fully recompile the example with G4_USE_FLUKA=yes.\n"
0214 "For example:\n"
0215 "source geant4/examples/extended/hadronic/FlukaCern/FlukaInterface/"
0216 "env_FLUKA_G4_interface.sh\n"
0217 "cd geant4/examples/extended/hadronic/FlukaCern/ProcessLevel/CrossSection/\n"
0218 "mkdir build\n"
0219 "cd build\n"
0220 "cmake -DGeant4_DIR=your_path_to_geant4 -DG4_USE_FLUKA=1 .. \n"
0221 "make -j8 G4_USE_FLUKA=1\n"
0222 "NB: First time use of FLUKA interface:\n"
0223 "Do not forget to first compile the FLUKA interface itself.\n"
0224 "For example: cd geant4/examples/extended/hadronic/FlukaCern/FlukaInterface/ "
0225 "&& make interface && make env\n"
0226 "FlukaInterface/env_FLUKA_G4_interface.sh then needs to be sourced\n"
0227 "in whichever terminal you want to use the FLUKA interface.\n");
0228 #endif
0229 }
0230
0231
0232 G4CascadeInterface* theBERTmodel = new G4CascadeInterface;
0233
0234
0235 G4BinaryCascade* theBICmodel = new G4BinaryCascade;
0236 G4PreCompoundModel* thePreEquilib = new G4PreCompoundModel(new G4ExcitationHandler);
0237 theBICmodel->SetDeExcitation(thePreEquilib);
0238
0239
0240 G4PreCompoundModel* thePreEquilibBis = new G4PreCompoundModel(new G4ExcitationHandler);
0241 G4BinaryLightIonReaction* theIonBICmodel = new G4BinaryLightIonReaction(thePreEquilibBis);
0242
0243
0244 G4INCLXXInterface* theINCLmodel = new G4INCLXXInterface;
0245 const G4bool useAblaDeExcitation = false;
0246
0247 if (theINCLmodel && useAblaDeExcitation) {
0248 G4AblaInterface* theAblaInterface = new G4AblaInterface;
0249 theINCLmodel->SetDeExcitation(theAblaInterface);
0250 }
0251
0252
0253
0254
0255
0256
0257 G4TheoFSGenerator* theFTFPmodel = new G4TheoFSGenerator("FTFP");
0258 theFTFPmodel->SetMaxEnergy(G4HadronicParameters::Instance()->GetMaxEnergy());
0259 G4GeneratorPrecompoundInterface* theCascade = new G4GeneratorPrecompoundInterface;
0260 theCascade->SetDeExcitation(thePreEquilib);
0261 theFTFPmodel->SetTransport(theCascade);
0262 G4LundStringFragmentation* theLundFragmentation = new G4LundStringFragmentation;
0263 G4ExcitedStringDecay* theStringDecay = new G4ExcitedStringDecay(theLundFragmentation);
0264 G4FTFModel* theStringModel = new G4FTFModel;
0265 theStringModel->SetFragmentationModel(theStringDecay);
0266
0267
0268
0269
0270
0271
0272 theFTFPmodel->SetHighEnergyGenerator(theStringModel);
0273
0274
0275
0276 G4TheoFSGenerator* theFTFPmodel_aboveThreshold = new G4TheoFSGenerator("FTFP");
0277 theFTFPmodel_aboveThreshold->SetMaxEnergy(G4HadronicParameters::Instance()->GetMaxEnergy());
0278 theFTFPmodel_aboveThreshold->SetTransport(theCascade);
0279 theFTFPmodel_aboveThreshold->SetHighEnergyGenerator(theStringModel);
0280
0281
0282 G4TheoFSGenerator* theFTFPmodel_constrained = new G4TheoFSGenerator("FTFP");
0283 theFTFPmodel_constrained->SetMaxEnergy(G4HadronicParameters::Instance()->GetMaxEnergy());
0284 theFTFPmodel_constrained->SetTransport(theCascade);
0285 theFTFPmodel_constrained->SetHighEnergyGenerator(theStringModel);
0286
0287
0288
0289 G4TheoFSGenerator* theFTFPmodel_belowThreshold = new G4TheoFSGenerator("FTFP");
0290 theFTFPmodel_belowThreshold->SetMaxEnergy(G4HadronicParameters::Instance()->GetMaxEnergy());
0291 theFTFPmodel_belowThreshold->SetTransport(theCascade);
0292 theFTFPmodel_belowThreshold->SetHighEnergyGenerator(theStringModel);
0293
0294
0295 G4TheoFSGenerator* theQGSPmodel = new G4TheoFSGenerator("QGSP");
0296 theQGSPmodel->SetMaxEnergy(G4HadronicParameters::Instance()->GetMaxEnergy());
0297 theQGSPmodel->SetTransport(theCascade);
0298 G4QGSMFragmentation* theQgsmFragmentation = new G4QGSMFragmentation;
0299 G4ExcitedStringDecay* theQgsmStringDecay = new G4ExcitedStringDecay(theQgsmFragmentation);
0300 G4VPartonStringModel* theQgsmStringModel = new G4QGSModel<G4QGSParticipants>;
0301 theQgsmStringModel->SetFragmentationModel(theQgsmStringDecay);
0302 theQGSPmodel->SetHighEnergyGenerator(theQgsmStringModel);
0303 G4QuasiElasticChannel* theQuasiElastic = new G4QuasiElasticChannel;
0304 theQGSPmodel->SetQuasiElasticChannel(theQuasiElastic);
0305
0306
0307
0308
0309
0310
0311
0312 if (fPhysicsCase == "FTFP_BERT_ATL" || fPhysicsCase == "FTFP_BERT"
0313 || fPhysicsCase == "FTFP_INCLXX" || fPhysicsCase == "QGSP_BERT" || fPhysicsCase == "QGSP_BIC")
0314 {
0315 const G4double ftfpMinE = G4HadronicParameters::Instance()->GetMinEnergyTransitionFTF_Cascade();
0316 const G4double bertMaxE = G4HadronicParameters::Instance()->GetMaxEnergyTransitionFTF_Cascade();
0317 const G4double ftfpMinE_ATL = 9.0 * CLHEP::GeV;
0318 const G4double bertMaxE_ATL = 12.0 * CLHEP::GeV;
0319 const G4double ftfpMaxE = G4HadronicParameters::Instance()->GetMaxEnergyTransitionQGS_FTF();
0320 const G4double qgspMinE = G4HadronicParameters::Instance()->GetMinEnergyTransitionQGS_FTF();
0321 theFTFPmodel->SetMinEnergy(0.0);
0322 theFTFPmodel_belowThreshold->SetMinEnergy(0.0);
0323 if (fPhysicsCase == "FTFP_BERT_ATL") {
0324 theBERTmodel->SetMaxEnergy(bertMaxE_ATL);
0325 theIonBICmodel->SetMaxEnergy(bertMaxE_ATL);
0326 theFTFPmodel_aboveThreshold->SetMinEnergy(ftfpMinE_ATL);
0327 theFTFPmodel_constrained->SetMinEnergy(ftfpMinE_ATL);
0328 }
0329 else {
0330 theBERTmodel->SetMaxEnergy(bertMaxE);
0331 theIonBICmodel->SetMaxEnergy(bertMaxE);
0332 theFTFPmodel_aboveThreshold->SetMinEnergy(ftfpMinE);
0333 theFTFPmodel_constrained->SetMinEnergy(ftfpMinE);
0334 }
0335 if (fPhysicsCase == "FTFP_INCLXX") {
0336 theINCLmodel->SetMaxEnergy(bertMaxE);
0337 }
0338 if (fPhysicsCase == "QGSP_BERT" || fPhysicsCase == "QGSP_BIC") {
0339 theFTFPmodel_constrained->SetMaxEnergy(ftfpMaxE);
0340 theFTFPmodel_belowThreshold->SetMaxEnergy(ftfpMaxE);
0341 theQGSPmodel->SetMinEnergy(qgspMinE);
0342 theBICmodel->SetMaxEnergy(bertMaxE);
0343 }
0344 }
0345
0346
0347 G4VCrossSectionDataSet* thePionMinusXSdata = new G4BGGPionInelasticXS(G4PionMinus::Definition());
0348 thePionMinusXSdata->BuildPhysicsTable(*(G4PionMinus::Definition()));
0349 G4VCrossSectionDataSet* thePionPlusXSdata = new G4BGGPionInelasticXS(G4PionPlus::Definition());
0350 thePionPlusXSdata->BuildPhysicsTable(*(G4PionPlus::Definition()));
0351 G4VCrossSectionDataSet* theKaonXSdata =
0352 new G4CrossSectionInelastic(new G4ComponentGGHadronNucleusXsc);
0353 theKaonXSdata->BuildPhysicsTable(*(G4KaonMinus::Definition()));
0354 theKaonXSdata->BuildPhysicsTable(*(G4KaonPlus::Definition()));
0355 theKaonXSdata->BuildPhysicsTable(*(G4KaonZeroLong::Definition()));
0356 theKaonXSdata->BuildPhysicsTable(*(G4KaonZeroShort::Definition()));
0357 G4VCrossSectionDataSet* theProtonXSdata = new G4BGGNucleonInelasticXS(G4Proton::Proton());
0358 theProtonXSdata->BuildPhysicsTable(*(G4Proton::Definition()));
0359 G4VCrossSectionDataSet* theNeutronXSdata = new G4NeutronInelasticXS;
0360 theNeutronXSdata->BuildPhysicsTable(*(G4Neutron::Definition()));
0361
0362
0363
0364 G4VCrossSectionDataSet* theHyperonsXSdata =
0365 new G4CrossSectionInelastic(new G4ComponentGGHadronNucleusXsc);
0366 G4VCrossSectionDataSet* theAntibaryonsXSdata =
0367 new G4CrossSectionInelastic(new G4ComponentAntiNuclNuclearXS);
0368 G4VCrossSectionDataSet* theNuclNuclXSdata =
0369 new G4CrossSectionInelastic(new G4ComponentGGNuclNuclXsc);
0370
0371
0372
0373 typedef std::pair<G4ParticleDefinition*, G4HadronicProcess*> ProcessPair;
0374 G4HadronicProcess* thePionMinusInelasticProcess =
0375 new G4HadronInelasticProcess("pi-Inelastic", G4PionMinus::Definition());
0376 fProcessMap.insert(ProcessPair(G4PionMinus::Definition(), thePionMinusInelasticProcess));
0377 G4HadronicProcess* thePionPlusInelasticProcess =
0378 new G4HadronInelasticProcess("pi+Inelastic", G4PionPlus::Definition());
0379 fProcessMap.insert(ProcessPair(G4PionPlus::Definition(), thePionPlusInelasticProcess));
0380 G4HadronicProcess* theKaonMinusInelasticProcess =
0381 new G4HadronInelasticProcess("kaon-Inelastic", G4KaonMinus::Definition());
0382 fProcessMap.insert(ProcessPair(G4KaonMinus::Definition(), theKaonMinusInelasticProcess));
0383 G4HadronicProcess* theKaonPlusInelasticProcess =
0384 new G4HadronInelasticProcess("kaon+Inelastic", G4KaonPlus::Definition());
0385 fProcessMap.insert(ProcessPair(G4KaonPlus::Definition(), theKaonPlusInelasticProcess));
0386 G4HadronicProcess* theKaonZeroLInelasticProcess =
0387 new G4HadronInelasticProcess("kaon0LInelastic", G4KaonZeroLong::Definition());
0388 fProcessMap.insert(ProcessPair(G4KaonZeroLong::Definition(), theKaonZeroLInelasticProcess));
0389 G4HadronicProcess* theKaonZeroSInelasticProcess =
0390 new G4HadronInelasticProcess("kaon0SInelastic", G4KaonZeroShort::Definition());
0391 fProcessMap.insert(ProcessPair(G4KaonZeroShort::Definition(), theKaonZeroSInelasticProcess));
0392 G4HadronicProcess* theProtonInelasticProcess =
0393 new G4HadronInelasticProcess("protonInelastic", G4Proton::Definition());
0394 fProcessMap.insert(ProcessPair(G4Proton::Definition(), theProtonInelasticProcess));
0395 G4HadronicProcess* theNeutronInelasticProcess =
0396 new G4HadronInelasticProcess("neutronInelastic", G4Neutron::Definition());
0397 fProcessMap.insert(ProcessPair(G4Neutron::Definition(), theNeutronInelasticProcess));
0398 G4HadronicProcess* theDeuteronInelasticProcess =
0399 new G4HadronInelasticProcess("dInelastic", G4Deuteron::Definition());
0400 fProcessMap.insert(ProcessPair(G4Deuteron::Definition(), theDeuteronInelasticProcess));
0401 G4HadronicProcess* theTritonInelasticProcess =
0402 new G4HadronInelasticProcess("tInelastic", G4Triton::Definition());
0403 fProcessMap.insert(ProcessPair(G4Triton::Definition(), theTritonInelasticProcess));
0404 G4HadronicProcess* theHe3InelasticProcess =
0405 new G4HadronInelasticProcess("he3Inelastic", G4He3::Definition());
0406 fProcessMap.insert(ProcessPair(G4He3::Definition(), theHe3InelasticProcess));
0407 G4HadronicProcess* theAlphaInelasticProcess =
0408 new G4HadronInelasticProcess("alphaInelastic", G4Alpha::Definition());
0409 fProcessMap.insert(ProcessPair(G4Alpha::Definition(), theAlphaInelasticProcess));
0410 G4HadronicProcess* theIonInelasticProcess =
0411 new G4HadronInelasticProcess("ionInelastic", G4GenericIon::Definition());
0412 fProcessMap.insert(ProcessPair(G4GenericIon::Definition(), theIonInelasticProcess));
0413 G4HadronicProcess* theLambdaInelasticProcess =
0414 new G4HadronInelasticProcess("lambdaInelastic", G4Lambda::Definition());
0415 fProcessMap.insert(ProcessPair(G4Lambda::Definition(), theLambdaInelasticProcess));
0416 G4HadronicProcess* theSigmaMinusInelasticProcess =
0417 new G4HadronInelasticProcess("sigma-Inelastic", G4SigmaMinus::Definition());
0418 fProcessMap.insert(ProcessPair(G4SigmaMinus::Definition(), theSigmaMinusInelasticProcess));
0419 G4HadronicProcess* theSigmaPlusInelasticProcess =
0420 new G4HadronInelasticProcess("sigma+Inelastic", G4SigmaPlus::Definition());
0421 fProcessMap.insert(ProcessPair(G4SigmaPlus::Definition(), theSigmaPlusInelasticProcess));
0422 G4HadronicProcess* theXiMinusInelasticProcess =
0423 new G4HadronInelasticProcess("xi-Inelastic", G4XiMinus::Definition());
0424 fProcessMap.insert(ProcessPair(G4XiMinus::Definition(), theXiMinusInelasticProcess));
0425 G4HadronicProcess* theXiZeroInelasticProcess =
0426 new G4HadronInelasticProcess("xi0Inelastic", G4XiZero::Definition());
0427 fProcessMap.insert(ProcessPair(G4XiZero::Definition(), theXiZeroInelasticProcess));
0428 G4HadronicProcess* theOmegaMinusInelasticProcess =
0429 new G4HadronInelasticProcess("omega-Inelastic", G4OmegaMinus::Definition());
0430 fProcessMap.insert(ProcessPair(G4OmegaMinus::Definition(), theOmegaMinusInelasticProcess));
0431 G4HadronicProcess* theAntiProtonInelasticProcess =
0432 new G4HadronInelasticProcess("anti_protonInelastic", G4AntiProton::Definition());
0433 fProcessMap.insert(ProcessPair(G4AntiProton::Definition(), theAntiProtonInelasticProcess));
0434 G4HadronicProcess* theAntiNeutronInelasticProcess =
0435 new G4HadronInelasticProcess("anti_neutronInelastic", G4AntiNeutron::Definition());
0436 fProcessMap.insert(ProcessPair(G4AntiNeutron::Definition(), theAntiNeutronInelasticProcess));
0437 G4HadronicProcess* theAntiDeuteronInelasticProcess =
0438 new G4HadronInelasticProcess("anti_deuteronInelastic", G4AntiDeuteron::Definition());
0439 fProcessMap.insert(ProcessPair(G4AntiDeuteron::Definition(), theAntiDeuteronInelasticProcess));
0440 G4HadronicProcess* theAntiTritonInelasticProcess =
0441 new G4HadronInelasticProcess("anti_tritonInelastic", G4AntiTriton::Definition());
0442 fProcessMap.insert(ProcessPair(G4AntiTriton::Definition(), theAntiTritonInelasticProcess));
0443 G4HadronicProcess* theAntiHe3InelasticProcess =
0444 new G4HadronInelasticProcess("anti_He3Inelastic", G4AntiHe3::Definition());
0445 fProcessMap.insert(ProcessPair(G4AntiHe3::Definition(), theAntiHe3InelasticProcess));
0446 G4HadronicProcess* theAntiAlphaInelasticProcess =
0447 new G4HadronInelasticProcess("anti_alphaInelastic", G4AntiAlpha::Definition());
0448 fProcessMap.insert(ProcessPair(G4AntiAlpha::Definition(), theAntiAlphaInelasticProcess));
0449 G4HadronicProcess* theAntiLambdaInelasticProcess =
0450 new G4HadronInelasticProcess("anti-lambdaInelastic", G4AntiLambda::Definition());
0451 fProcessMap.insert(ProcessPair(G4AntiLambda::Definition(), theAntiLambdaInelasticProcess));
0452 G4HadronicProcess* theAntiSigmaMinusInelasticProcess =
0453 new G4HadronInelasticProcess("anti_sigma-Inelastic", G4AntiSigmaMinus::Definition());
0454 fProcessMap.insert(
0455 ProcessPair(G4AntiSigmaMinus::Definition(), theAntiSigmaMinusInelasticProcess));
0456 G4HadronicProcess* theAntiSigmaPlusInelasticProcess =
0457 new G4HadronInelasticProcess("anti_sigma+Inelastic", G4AntiSigmaPlus::Definition());
0458 fProcessMap.insert(ProcessPair(G4AntiSigmaPlus::Definition(), theAntiSigmaPlusInelasticProcess));
0459 G4HadronicProcess* theAntiXiMinusInelasticProcess =
0460 new G4HadronInelasticProcess("anti_xi-Inelastic", G4AntiXiMinus::Definition());
0461 fProcessMap.insert(ProcessPair(G4AntiXiMinus::Definition(), theAntiXiMinusInelasticProcess));
0462 G4HadronicProcess* theAntiXiZeroInelasticProcess =
0463 new G4HadronInelasticProcess("anti_xi0Inelastic", G4AntiXiZero::Definition());
0464 fProcessMap.insert(ProcessPair(G4AntiXiZero::Definition(), theAntiXiZeroInelasticProcess));
0465 G4HadronicProcess* theAntiOmegaMinusInelasticProcess =
0466 new G4HadronInelasticProcess("anti_omega-Inelastic", G4AntiOmegaMinus::Definition());
0467 fProcessMap.insert(
0468 ProcessPair(G4AntiOmegaMinus::Definition(), theAntiOmegaMinusInelasticProcess));
0469
0470 G4HadronicProcess* theDPlusInelasticProcess =
0471 new G4HadronInelasticProcess("D+Inelastic", G4DMesonPlus::Definition());
0472 fProcessMap.insert(ProcessPair(G4DMesonPlus::Definition(), theDPlusInelasticProcess));
0473 G4HadronicProcess* theDMinusInelasticProcess =
0474 new G4HadronInelasticProcess("D-Inelastic", G4DMesonMinus::Definition());
0475 fProcessMap.insert(ProcessPair(G4DMesonMinus::Definition(), theDMinusInelasticProcess));
0476 G4HadronicProcess* theDZeroInelasticProcess =
0477 new G4HadronInelasticProcess("D0Inelastic", G4DMesonZero::Definition());
0478 fProcessMap.insert(ProcessPair(G4DMesonZero::Definition(), theDZeroInelasticProcess));
0479 G4HadronicProcess* theAntiDZeroInelasticProcess =
0480 new G4HadronInelasticProcess("anti_D0Inelastic", G4AntiDMesonZero::Definition());
0481 fProcessMap.insert(ProcessPair(G4AntiDMesonZero::Definition(), theAntiDZeroInelasticProcess));
0482 G4HadronicProcess* theDsPlusInelasticProcess =
0483 new G4HadronInelasticProcess("Ds+Inelastic", G4DsMesonPlus::Definition());
0484 fProcessMap.insert(ProcessPair(G4DsMesonPlus::Definition(), theDsPlusInelasticProcess));
0485 G4HadronicProcess* theDsMinusInelasticProcess =
0486 new G4HadronInelasticProcess("Ds-Inelastic", G4DsMesonMinus::Definition());
0487 fProcessMap.insert(ProcessPair(G4DsMesonMinus::Definition(), theDsMinusInelasticProcess));
0488 G4HadronicProcess* theBPlusInelasticProcess =
0489 new G4HadronInelasticProcess("B+Inelastic", G4BMesonPlus::Definition());
0490 fProcessMap.insert(ProcessPair(G4BMesonPlus::Definition(), theBPlusInelasticProcess));
0491 G4HadronicProcess* theBMinusInelasticProcess =
0492 new G4HadronInelasticProcess("B-Inelastic", G4BMesonMinus::Definition());
0493 fProcessMap.insert(ProcessPair(G4BMesonMinus::Definition(), theBMinusInelasticProcess));
0494 G4HadronicProcess* theBZeroInelasticProcess =
0495 new G4HadronInelasticProcess("B0Inelastic", G4BMesonZero::Definition());
0496 fProcessMap.insert(ProcessPair(G4BMesonZero::Definition(), theBZeroInelasticProcess));
0497 G4HadronicProcess* theAntiBZeroInelasticProcess =
0498 new G4HadronInelasticProcess("anti_B0Inelastic", G4AntiBMesonZero::Definition());
0499 fProcessMap.insert(ProcessPair(G4AntiBMesonZero::Definition(), theAntiBZeroInelasticProcess));
0500 G4HadronicProcess* theBsZeroInelasticProcess =
0501 new G4HadronInelasticProcess("Bs0Inelastic", G4BsMesonZero::Definition());
0502 fProcessMap.insert(ProcessPair(G4BsMesonZero::Definition(), theBsZeroInelasticProcess));
0503 G4HadronicProcess* theAntiBsZeroInelasticProcess =
0504 new G4HadronInelasticProcess("anti_Bs0Inelastic", G4AntiBsMesonZero::Definition());
0505 fProcessMap.insert(ProcessPair(G4AntiBsMesonZero::Definition(), theAntiBsZeroInelasticProcess));
0506 G4HadronicProcess* theBcPlusInelasticProcess =
0507 new G4HadronInelasticProcess("Bc+Inelastic", G4BcMesonPlus::Definition());
0508 fProcessMap.insert(ProcessPair(G4BcMesonPlus::Definition(), theBcPlusInelasticProcess));
0509 G4HadronicProcess* theBcMinusInelasticProcess =
0510 new G4HadronInelasticProcess("Bc-Inelastic", G4BcMesonMinus::Definition());
0511 fProcessMap.insert(ProcessPair(G4BcMesonMinus::Definition(), theBcMinusInelasticProcess));
0512 G4HadronicProcess* theLambdacPlusInelasticProcess =
0513 new G4HadronInelasticProcess("lambda_c+Inelastic", G4LambdacPlus::Definition());
0514 fProcessMap.insert(ProcessPair(G4LambdacPlus::Definition(), theLambdacPlusInelasticProcess));
0515 G4HadronicProcess* theAntiLambdacPlusInelasticProcess =
0516 new G4HadronInelasticProcess("anti_lambda_c+Inelastic", G4AntiLambdacPlus::Definition());
0517 fProcessMap.insert(
0518 ProcessPair(G4AntiLambdacPlus::Definition(), theAntiLambdacPlusInelasticProcess));
0519 G4HadronicProcess* theXicPlusInelasticProcess =
0520 new G4HadronInelasticProcess("xi_c+Inelastic", G4XicPlus::Definition());
0521 fProcessMap.insert(ProcessPair(G4XicPlus::Definition(), theXicPlusInelasticProcess));
0522 G4HadronicProcess* theAntiXicPlusInelasticProcess =
0523 new G4HadronInelasticProcess("anti_xi_c+Inelastic", G4AntiXicPlus::Definition());
0524 fProcessMap.insert(ProcessPair(G4AntiXicPlus::Definition(), theAntiXicPlusInelasticProcess));
0525 G4HadronicProcess* theXicZeroInelasticProcess =
0526 new G4HadronInelasticProcess("xi_c0Inelastic", G4XicZero::Definition());
0527 fProcessMap.insert(ProcessPair(G4XicZero::Definition(), theXicZeroInelasticProcess));
0528 G4HadronicProcess* theAntiXicZeroInelasticProcess =
0529 new G4HadronInelasticProcess("anti_xi_c0Inelastic", G4AntiXicZero::Definition());
0530 fProcessMap.insert(ProcessPair(G4AntiXicZero::Definition(), theAntiXicZeroInelasticProcess));
0531 G4HadronicProcess* theOmegacZeroInelasticProcess =
0532 new G4HadronInelasticProcess("omega_c0Inelastic", G4OmegacZero::Definition());
0533 fProcessMap.insert(ProcessPair(G4OmegacZero::Definition(), theOmegacZeroInelasticProcess));
0534 G4HadronicProcess* theAntiOmegacZeroInelasticProcess =
0535 new G4HadronInelasticProcess("anti_omega_c0Inelastic", G4AntiOmegacZero::Definition());
0536 fProcessMap.insert(
0537 ProcessPair(G4AntiOmegacZero::Definition(), theAntiOmegacZeroInelasticProcess));
0538 G4HadronicProcess* theLambdabInelasticProcess =
0539 new G4HadronInelasticProcess("lambda_bInelastic", G4Lambdab::Definition());
0540 fProcessMap.insert(ProcessPair(G4Lambdab::Definition(), theLambdabInelasticProcess));
0541 G4HadronicProcess* theAntiLambdabInelasticProcess =
0542 new G4HadronInelasticProcess("anti_lambda_bInelastic", G4AntiLambdab::Definition());
0543 fProcessMap.insert(ProcessPair(G4AntiLambdab::Definition(), theAntiLambdabInelasticProcess));
0544 G4HadronicProcess* theXibZeroInelasticProcess =
0545 new G4HadronInelasticProcess("xi_b0Inelastic", G4XibZero::Definition());
0546 fProcessMap.insert(ProcessPair(G4XibZero::Definition(), theXibZeroInelasticProcess));
0547 G4HadronicProcess* theAntiXibZeroInelasticProcess =
0548 new G4HadronInelasticProcess("anti_xi_b0Inelastic", G4AntiXibZero::Definition());
0549 fProcessMap.insert(ProcessPair(G4AntiXibZero::Definition(), theAntiXibZeroInelasticProcess));
0550 G4HadronicProcess* theXibMinusInelasticProcess =
0551 new G4HadronInelasticProcess("xi_b-Inelastic", G4XibMinus::Definition());
0552 fProcessMap.insert(ProcessPair(G4XibMinus::Definition(), theXibMinusInelasticProcess));
0553 G4HadronicProcess* theAntiXibMinusInelasticProcess =
0554 new G4HadronInelasticProcess("anti_xi_b-Inelastic", G4AntiXibMinus::Definition());
0555 fProcessMap.insert(ProcessPair(G4AntiXibMinus::Definition(), theAntiXibMinusInelasticProcess));
0556 G4HadronicProcess* theOmegabMinusInelasticProcess =
0557 new G4HadronInelasticProcess("omega_b-Inelastic", G4OmegabMinus::Definition());
0558 fProcessMap.insert(ProcessPair(G4OmegabMinus::Definition(), theOmegabMinusInelasticProcess));
0559 G4HadronicProcess* theAntiOmegabMinusInelasticProcess =
0560 new G4HadronInelasticProcess("anti_omega_b-Inelastic", G4AntiOmegabMinus::Definition());
0561 fProcessMap.insert(
0562 ProcessPair(G4AntiOmegabMinus::Definition(), theAntiOmegabMinusInelasticProcess));
0563
0564
0565 if (fPhysicsCase == "CFLUKAHI") {
0566 #ifdef G4_USE_FLUKA
0567 fPhysicsCaseIsSupported = true;
0568
0569
0570 auto flukaInelasticScatteringXS = new FLUKAInelasticScatteringXS();
0571
0572 theProtonInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0573 theAntiProtonInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0574 theNeutronInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0575 theAntiNeutronInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0576 thePionMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0577 thePionPlusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0578 theKaonMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0579 theKaonPlusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0580 theKaonZeroLInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0581 theKaonZeroSInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0582
0583
0584
0585
0586
0587 theLambdaInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0588 theSigmaMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0589 theSigmaPlusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0590 theXiMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0591 theXiZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0592 theOmegaMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0593 theAntiLambdaInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0594 theAntiSigmaMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0595 theAntiSigmaPlusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0596 theAntiXiMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0597 theAntiXiZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0598 theAntiOmegaMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0599 theDPlusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0600 theDMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0601 theDZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0602 theAntiDZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0603 theDsPlusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0604 theDsMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0605 theBPlusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0606 theBMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0607 theBZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0608 theAntiBZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0609 theBsZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0610 theAntiBsZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0611 theBcPlusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0612 theBcMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0613 theLambdacPlusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0614 theAntiLambdacPlusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0615 theXicPlusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0616 theAntiXicPlusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0617 theXicZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0618 theAntiXicZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0619 theOmegacZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0620 theAntiOmegacZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0621 theLambdabInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0622 theAntiLambdabInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0623 theXibZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0624 theAntiXibZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0625 theXibMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0626 theAntiXibMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0627 theOmegabMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0628 theAntiOmegabMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0629
0630
0631
0632 theDeuteronInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0633 theTritonInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0634 theHe3InelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0635 theAlphaInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0636
0637
0638
0639 theAntiDeuteronInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0640 theAntiTritonInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0641 theAntiHe3InelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0642 theAntiAlphaInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0643
0644
0645
0646 theIonInelasticProcess->AddDataSet(flukaInelasticScatteringXS);
0647 #endif
0648 }
0649
0650 else {
0651
0652 thePionMinusInelasticProcess->AddDataSet(thePionMinusXSdata);
0653 thePionPlusInelasticProcess->AddDataSet(thePionPlusXSdata);
0654 theKaonMinusInelasticProcess->AddDataSet(theKaonXSdata);
0655 theKaonPlusInelasticProcess->AddDataSet(theKaonXSdata);
0656 theKaonZeroLInelasticProcess->AddDataSet(theKaonXSdata);
0657 theKaonZeroSInelasticProcess->AddDataSet(theKaonXSdata);
0658 theProtonInelasticProcess->AddDataSet(theProtonXSdata);
0659 theNeutronInelasticProcess->AddDataSet(theNeutronXSdata);
0660 theDeuteronInelasticProcess->AddDataSet(theNuclNuclXSdata);
0661 theTritonInelasticProcess->AddDataSet(theNuclNuclXSdata);
0662 theHe3InelasticProcess->AddDataSet(theNuclNuclXSdata);
0663 theAlphaInelasticProcess->AddDataSet(theNuclNuclXSdata);
0664 theIonInelasticProcess->AddDataSet(theNuclNuclXSdata);
0665 theLambdaInelasticProcess->AddDataSet(theHyperonsXSdata);
0666 theSigmaMinusInelasticProcess->AddDataSet(theHyperonsXSdata);
0667 theSigmaPlusInelasticProcess->AddDataSet(theHyperonsXSdata);
0668 theXiMinusInelasticProcess->AddDataSet(theHyperonsXSdata);
0669 theXiZeroInelasticProcess->AddDataSet(theHyperonsXSdata);
0670 theOmegaMinusInelasticProcess->AddDataSet(theHyperonsXSdata);
0671 theAntiProtonInelasticProcess->AddDataSet(theAntibaryonsXSdata);
0672 theAntiNeutronInelasticProcess->AddDataSet(theAntibaryonsXSdata);
0673 theAntiDeuteronInelasticProcess->AddDataSet(theAntibaryonsXSdata);
0674 theAntiTritonInelasticProcess->AddDataSet(theAntibaryonsXSdata);
0675 theAntiHe3InelasticProcess->AddDataSet(theAntibaryonsXSdata);
0676 theAntiAlphaInelasticProcess->AddDataSet(theHyperonsXSdata);
0677 theAntiLambdaInelasticProcess->AddDataSet(theHyperonsXSdata);
0678 theAntiSigmaMinusInelasticProcess->AddDataSet(theHyperonsXSdata);
0679 theAntiSigmaPlusInelasticProcess->AddDataSet(theHyperonsXSdata);
0680 theAntiXiMinusInelasticProcess->AddDataSet(theHyperonsXSdata);
0681 theAntiXiZeroInelasticProcess->AddDataSet(theHyperonsXSdata);
0682 theAntiOmegaMinusInelasticProcess->AddDataSet(theHyperonsXSdata);
0683
0684 theDPlusInelasticProcess->AddDataSet(theHyperonsXSdata);
0685 theDMinusInelasticProcess->AddDataSet(theHyperonsXSdata);
0686 theDZeroInelasticProcess->AddDataSet(theHyperonsXSdata);
0687 theAntiDZeroInelasticProcess->AddDataSet(theHyperonsXSdata);
0688 theDsPlusInelasticProcess->AddDataSet(theHyperonsXSdata);
0689 theDsMinusInelasticProcess->AddDataSet(theHyperonsXSdata);
0690 theBPlusInelasticProcess->AddDataSet(theHyperonsXSdata);
0691 theBMinusInelasticProcess->AddDataSet(theHyperonsXSdata);
0692 theBZeroInelasticProcess->AddDataSet(theHyperonsXSdata);
0693 theAntiBZeroInelasticProcess->AddDataSet(theHyperonsXSdata);
0694 theBsZeroInelasticProcess->AddDataSet(theHyperonsXSdata);
0695 theAntiBsZeroInelasticProcess->AddDataSet(theHyperonsXSdata);
0696 theBcPlusInelasticProcess->AddDataSet(theHyperonsXSdata);
0697 theBcMinusInelasticProcess->AddDataSet(theHyperonsXSdata);
0698 theLambdacPlusInelasticProcess->AddDataSet(theHyperonsXSdata);
0699 theAntiLambdacPlusInelasticProcess->AddDataSet(theHyperonsXSdata);
0700 theXicPlusInelasticProcess->AddDataSet(theHyperonsXSdata);
0701 theAntiXicPlusInelasticProcess->AddDataSet(theHyperonsXSdata);
0702 theXicZeroInelasticProcess->AddDataSet(theHyperonsXSdata);
0703 theAntiXicZeroInelasticProcess->AddDataSet(theHyperonsXSdata);
0704 theOmegacZeroInelasticProcess->AddDataSet(theHyperonsXSdata);
0705 theAntiOmegacZeroInelasticProcess->AddDataSet(theHyperonsXSdata);
0706 theLambdabInelasticProcess->AddDataSet(theHyperonsXSdata);
0707 theAntiLambdabInelasticProcess->AddDataSet(theHyperonsXSdata);
0708 theXibZeroInelasticProcess->AddDataSet(theHyperonsXSdata);
0709 theAntiXibZeroInelasticProcess->AddDataSet(theHyperonsXSdata);
0710 theXibMinusInelasticProcess->AddDataSet(theHyperonsXSdata);
0711 theAntiXibMinusInelasticProcess->AddDataSet(theHyperonsXSdata);
0712 theOmegabMinusInelasticProcess->AddDataSet(theHyperonsXSdata);
0713 theAntiOmegabMinusInelasticProcess->AddDataSet(theHyperonsXSdata);
0714 }
0715
0716
0717 if (fPhysicsCase == "CFLUKAHI") {
0718 #ifdef G4_USE_FLUKA
0719
0720
0721 FLUKANuclearInelasticModel* flukaModel = new FLUKANuclearInelasticModel();
0722
0723 theProtonInelasticProcess->RegisterMe(flukaModel);
0724 theAntiProtonInelasticProcess->RegisterMe(flukaModel);
0725 theNeutronInelasticProcess->RegisterMe(flukaModel);
0726 theAntiNeutronInelasticProcess->RegisterMe(flukaModel);
0727 thePionMinusInelasticProcess->RegisterMe(flukaModel);
0728 thePionPlusInelasticProcess->RegisterMe(flukaModel);
0729 theKaonMinusInelasticProcess->RegisterMe(flukaModel);
0730 theKaonPlusInelasticProcess->RegisterMe(flukaModel);
0731 theKaonZeroLInelasticProcess->RegisterMe(flukaModel);
0732 theKaonZeroSInelasticProcess->RegisterMe(flukaModel);
0733
0734
0735
0736
0737
0738 theLambdaInelasticProcess->RegisterMe(flukaModel);
0739 theAntiLambdaInelasticProcess->RegisterMe(flukaModel);
0740 theSigmaMinusInelasticProcess->RegisterMe(flukaModel);
0741 theAntiSigmaMinusInelasticProcess->RegisterMe(flukaModel);
0742 theSigmaPlusInelasticProcess->RegisterMe(flukaModel);
0743 theAntiSigmaPlusInelasticProcess->RegisterMe(flukaModel);
0744 theXiMinusInelasticProcess->RegisterMe(flukaModel);
0745 theAntiXiMinusInelasticProcess->RegisterMe(flukaModel);
0746 theXiZeroInelasticProcess->RegisterMe(flukaModel);
0747 theAntiXiZeroInelasticProcess->RegisterMe(flukaModel);
0748 theOmegaMinusInelasticProcess->RegisterMe(flukaModel);
0749 theAntiOmegaMinusInelasticProcess->RegisterMe(flukaModel);
0750 theDPlusInelasticProcess->RegisterMe(flukaModel);
0751 theDMinusInelasticProcess->RegisterMe(flukaModel);
0752 theDZeroInelasticProcess->RegisterMe(flukaModel);
0753 theAntiDZeroInelasticProcess->RegisterMe(flukaModel);
0754 theDsPlusInelasticProcess->RegisterMe(flukaModel);
0755 theDsMinusInelasticProcess->RegisterMe(flukaModel);
0756 theBPlusInelasticProcess->RegisterMe(flukaModel);
0757 theBMinusInelasticProcess->RegisterMe(flukaModel);
0758 theBZeroInelasticProcess->RegisterMe(flukaModel);
0759 theAntiBZeroInelasticProcess->RegisterMe(flukaModel);
0760 theBsZeroInelasticProcess->RegisterMe(flukaModel);
0761 theAntiBsZeroInelasticProcess->RegisterMe(flukaModel);
0762 theBcPlusInelasticProcess->RegisterMe(flukaModel);
0763 theBcMinusInelasticProcess->RegisterMe(flukaModel);
0764 theLambdacPlusInelasticProcess->RegisterMe(flukaModel);
0765 theAntiLambdacPlusInelasticProcess->RegisterMe(flukaModel);
0766 theXicPlusInelasticProcess->RegisterMe(flukaModel);
0767 theAntiXicPlusInelasticProcess->RegisterMe(flukaModel);
0768 theXicZeroInelasticProcess->RegisterMe(flukaModel);
0769 theAntiXicZeroInelasticProcess->RegisterMe(flukaModel);
0770 theOmegacZeroInelasticProcess->RegisterMe(flukaModel);
0771 theAntiOmegacZeroInelasticProcess->RegisterMe(flukaModel);
0772 theLambdabInelasticProcess->RegisterMe(flukaModel);
0773 theAntiLambdabInelasticProcess->RegisterMe(flukaModel);
0774 theXibZeroInelasticProcess->RegisterMe(flukaModel);
0775 theAntiXibZeroInelasticProcess->RegisterMe(flukaModel);
0776 theXibMinusInelasticProcess->RegisterMe(flukaModel);
0777 theAntiXibMinusInelasticProcess->RegisterMe(flukaModel);
0778 theOmegabMinusInelasticProcess->RegisterMe(flukaModel);
0779 theAntiOmegabMinusInelasticProcess->RegisterMe(flukaModel);
0780
0781
0782
0783 theDeuteronInelasticProcess->RegisterMe(flukaModel);
0784 theTritonInelasticProcess->RegisterMe(flukaModel);
0785 theHe3InelasticProcess->RegisterMe(flukaModel);
0786 theAlphaInelasticProcess->RegisterMe(flukaModel);
0787
0788
0789
0790 theAntiDeuteronInelasticProcess->RegisterMe(flukaModel);
0791 theAntiTritonInelasticProcess->RegisterMe(flukaModel);
0792 theAntiHe3InelasticProcess->RegisterMe(flukaModel);
0793 theAntiAlphaInelasticProcess->RegisterMe(flukaModel);
0794
0795
0796
0797 theIonInelasticProcess->RegisterMe(flukaModel);
0798 #endif
0799 }
0800
0801
0802
0803
0804
0805
0806
0807
0808
0809
0810
0811 if (fPhysicsCase == "BIC" || fPhysicsCase == "QGSP_BIC") {
0812
0813
0814 fPhysicsCaseIsSupported = true;
0815 theProtonInelasticProcess->RegisterMe(theBICmodel);
0816 theNeutronInelasticProcess->RegisterMe(theBICmodel);
0817 if (fPhysicsCase == "BIC") {
0818 thePionMinusInelasticProcess->RegisterMe(theBICmodel);
0819 thePionPlusInelasticProcess->RegisterMe(theBICmodel);
0820 }
0821 else {
0822 thePionMinusInelasticProcess->RegisterMe(theBERTmodel);
0823 thePionPlusInelasticProcess->RegisterMe(theBERTmodel);
0824 }
0825 }
0826 else if (fPhysicsCase == "INCL" || fPhysicsCase == "FTFP_INCLXX") {
0827
0828
0829 fPhysicsCaseIsSupported = true;
0830 thePionMinusInelasticProcess->RegisterMe(theINCLmodel);
0831 thePionPlusInelasticProcess->RegisterMe(theINCLmodel);
0832 theProtonInelasticProcess->RegisterMe(theINCLmodel);
0833 theNeutronInelasticProcess->RegisterMe(theINCLmodel);
0834 }
0835 if (fPhysicsCase == "IonBIC" || fPhysicsCase == "FTFP_BERT_ATL" || fPhysicsCase == "FTFP_BERT"
0836 || fPhysicsCase == "FTFP_INCLXX" || fPhysicsCase == "QGSP_BERT" || fPhysicsCase == "QGSP_BIC")
0837 {
0838
0839 fPhysicsCaseIsSupported = true;
0840 theDeuteronInelasticProcess->RegisterMe(theIonBICmodel);
0841 theTritonInelasticProcess->RegisterMe(theIonBICmodel);
0842 theHe3InelasticProcess->RegisterMe(theIonBICmodel);
0843 theAlphaInelasticProcess->RegisterMe(theIonBICmodel);
0844 theIonInelasticProcess->RegisterMe(theIonBICmodel);
0845 }
0846 if (fPhysicsCase == "QGSP" || fPhysicsCase == "QGSP_BERT" || fPhysicsCase == "QGSP_BIC") {
0847 fPhysicsCaseIsSupported = true;
0848 thePionMinusInelasticProcess->RegisterMe(theQGSPmodel);
0849 thePionPlusInelasticProcess->RegisterMe(theQGSPmodel);
0850 theKaonMinusInelasticProcess->RegisterMe(theQGSPmodel);
0851 theKaonPlusInelasticProcess->RegisterMe(theQGSPmodel);
0852 theKaonZeroLInelasticProcess->RegisterMe(theQGSPmodel);
0853 theKaonZeroSInelasticProcess->RegisterMe(theQGSPmodel);
0854 theProtonInelasticProcess->RegisterMe(theQGSPmodel);
0855 theNeutronInelasticProcess->RegisterMe(theQGSPmodel);
0856 theLambdaInelasticProcess->RegisterMe(theQGSPmodel);
0857 theSigmaMinusInelasticProcess->RegisterMe(theQGSPmodel);
0858 theSigmaPlusInelasticProcess->RegisterMe(theQGSPmodel);
0859 theXiMinusInelasticProcess->RegisterMe(theQGSPmodel);
0860 theXiZeroInelasticProcess->RegisterMe(theQGSPmodel);
0861 theOmegaMinusInelasticProcess->RegisterMe(theQGSPmodel);
0862 theAntiProtonInelasticProcess->RegisterMe(theQGSPmodel);
0863 theAntiNeutronInelasticProcess->RegisterMe(theQGSPmodel);
0864 theAntiLambdaInelasticProcess->RegisterMe(theQGSPmodel);
0865 theAntiSigmaMinusInelasticProcess->RegisterMe(theQGSPmodel);
0866 theAntiSigmaPlusInelasticProcess->RegisterMe(theQGSPmodel);
0867 theAntiXiMinusInelasticProcess->RegisterMe(theQGSPmodel);
0868 theAntiXiZeroInelasticProcess->RegisterMe(theQGSPmodel);
0869 theAntiOmegaMinusInelasticProcess->RegisterMe(theQGSPmodel);
0870 theDPlusInelasticProcess->RegisterMe(theQGSPmodel);
0871 theDMinusInelasticProcess->RegisterMe(theQGSPmodel);
0872 theDZeroInelasticProcess->RegisterMe(theQGSPmodel);
0873 theAntiDZeroInelasticProcess->RegisterMe(theQGSPmodel);
0874 theDsPlusInelasticProcess->RegisterMe(theQGSPmodel);
0875 theDsMinusInelasticProcess->RegisterMe(theQGSPmodel);
0876 theBPlusInelasticProcess->RegisterMe(theQGSPmodel);
0877 theBMinusInelasticProcess->RegisterMe(theQGSPmodel);
0878 theBZeroInelasticProcess->RegisterMe(theQGSPmodel);
0879 theAntiBZeroInelasticProcess->RegisterMe(theQGSPmodel);
0880 theBsZeroInelasticProcess->RegisterMe(theQGSPmodel);
0881 theAntiBsZeroInelasticProcess->RegisterMe(theQGSPmodel);
0882 theBcPlusInelasticProcess->RegisterMe(theQGSPmodel);
0883 theBcMinusInelasticProcess->RegisterMe(theQGSPmodel);
0884 theLambdacPlusInelasticProcess->RegisterMe(theQGSPmodel);
0885 theAntiLambdacPlusInelasticProcess->RegisterMe(theQGSPmodel);
0886 theXicPlusInelasticProcess->RegisterMe(theQGSPmodel);
0887 theAntiXicPlusInelasticProcess->RegisterMe(theQGSPmodel);
0888 theXicZeroInelasticProcess->RegisterMe(theQGSPmodel);
0889 theAntiXicZeroInelasticProcess->RegisterMe(theQGSPmodel);
0890 theOmegacZeroInelasticProcess->RegisterMe(theQGSPmodel);
0891 theAntiOmegacZeroInelasticProcess->RegisterMe(theQGSPmodel);
0892 theLambdabInelasticProcess->RegisterMe(theQGSPmodel);
0893 theAntiLambdabInelasticProcess->RegisterMe(theQGSPmodel);
0894 theXibZeroInelasticProcess->RegisterMe(theQGSPmodel);
0895 theAntiXibZeroInelasticProcess->RegisterMe(theQGSPmodel);
0896 theXibMinusInelasticProcess->RegisterMe(theQGSPmodel);
0897 theAntiXibMinusInelasticProcess->RegisterMe(theQGSPmodel);
0898 theOmegabMinusInelasticProcess->RegisterMe(theQGSPmodel);
0899 theAntiOmegabMinusInelasticProcess->RegisterMe(theQGSPmodel);
0900 }
0901 if (fPhysicsCase == "BERT" || fPhysicsCase == "FTFP_BERT_ATL" || fPhysicsCase == "FTFP_BERT"
0902 || fPhysicsCase == "QGSP_BERT")
0903 {
0904
0905 fPhysicsCaseIsSupported = true;
0906 thePionMinusInelasticProcess->RegisterMe(theBERTmodel);
0907 thePionPlusInelasticProcess->RegisterMe(theBERTmodel);
0908 theProtonInelasticProcess->RegisterMe(theBERTmodel);
0909 theNeutronInelasticProcess->RegisterMe(theBERTmodel);
0910 }
0911 if (fPhysicsCase == "BERT" || fPhysicsCase == "FTFP_BERT_ATL" || fPhysicsCase == "FTFP_BERT"
0912 || fPhysicsCase == "FTFP_INCLXX" || fPhysicsCase == "QGSP_BERT" || fPhysicsCase == "QGSP_BIC")
0913 {
0914
0915 fPhysicsCaseIsSupported = true;
0916 theKaonMinusInelasticProcess->RegisterMe(theBERTmodel);
0917 theKaonPlusInelasticProcess->RegisterMe(theBERTmodel);
0918 theKaonZeroLInelasticProcess->RegisterMe(theBERTmodel);
0919 theKaonZeroSInelasticProcess->RegisterMe(theBERTmodel);
0920 theLambdaInelasticProcess->RegisterMe(theBERTmodel);
0921 theSigmaMinusInelasticProcess->RegisterMe(theBERTmodel);
0922 theSigmaPlusInelasticProcess->RegisterMe(theBERTmodel);
0923 theXiMinusInelasticProcess->RegisterMe(theBERTmodel);
0924 theXiZeroInelasticProcess->RegisterMe(theBERTmodel);
0925 theOmegaMinusInelasticProcess->RegisterMe(theBERTmodel);
0926 }
0927 if (fPhysicsCase == "FTFP" || fPhysicsCase == "FTFP_BERT_ATL" || fPhysicsCase == "FTFP_BERT"
0928 || fPhysicsCase == "FTFP_INCLXX" || fPhysicsCase == "QGSP_BERT" || fPhysicsCase == "QGSP_BIC")
0929 {
0930
0931
0932 fPhysicsCaseIsSupported = true;
0933 theAntiDeuteronInelasticProcess->RegisterMe(theFTFPmodel);
0934 theAntiTritonInelasticProcess->RegisterMe(theFTFPmodel);
0935 theAntiHe3InelasticProcess->RegisterMe(theFTFPmodel);
0936 theAntiAlphaInelasticProcess->RegisterMe(theFTFPmodel);
0937 G4TheoFSGenerator* theFTFPmodelToBeUsed = theFTFPmodel_aboveThreshold;
0938 if (fPhysicsCase == "FTFP") {
0939 theFTFPmodelToBeUsed = theFTFPmodel;
0940 }
0941 else if (fPhysicsCase == "QGSP_BERT" || fPhysicsCase == "QGSP_BIC") {
0942 theFTFPmodelToBeUsed = theFTFPmodel_constrained;
0943 }
0944 thePionMinusInelasticProcess->RegisterMe(theFTFPmodelToBeUsed);
0945 thePionPlusInelasticProcess->RegisterMe(theFTFPmodelToBeUsed);
0946 theKaonMinusInelasticProcess->RegisterMe(theFTFPmodelToBeUsed);
0947 theKaonPlusInelasticProcess->RegisterMe(theFTFPmodelToBeUsed);
0948 theKaonZeroLInelasticProcess->RegisterMe(theFTFPmodelToBeUsed);
0949 theKaonZeroSInelasticProcess->RegisterMe(theFTFPmodelToBeUsed);
0950 theProtonInelasticProcess->RegisterMe(theFTFPmodelToBeUsed);
0951 theAntiProtonInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0952 theNeutronInelasticProcess->RegisterMe(theFTFPmodelToBeUsed);
0953 theAntiNeutronInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0954 theLambdaInelasticProcess->RegisterMe(theFTFPmodelToBeUsed);
0955 theAntiLambdaInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0956 theSigmaMinusInelasticProcess->RegisterMe(theFTFPmodelToBeUsed);
0957 theAntiSigmaMinusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0958 theSigmaPlusInelasticProcess->RegisterMe(theFTFPmodelToBeUsed);
0959 theAntiSigmaPlusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0960 theXiMinusInelasticProcess->RegisterMe(theFTFPmodelToBeUsed);
0961 theAntiXiMinusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0962 theXiZeroInelasticProcess->RegisterMe(theFTFPmodelToBeUsed);
0963 theAntiXiZeroInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0964 theOmegaMinusInelasticProcess->RegisterMe(theFTFPmodelToBeUsed);
0965 theAntiOmegaMinusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0966 theDPlusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0967 theDMinusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0968 theDZeroInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0969 theAntiDZeroInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0970 theDsPlusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0971 theDsMinusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0972 theBPlusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0973 theBMinusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0974 theBZeroInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0975 theAntiBZeroInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0976 theBsZeroInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0977 theAntiBsZeroInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0978 theBcPlusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0979 theBcMinusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0980 theLambdacPlusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0981 theAntiLambdacPlusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0982 theXicPlusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0983 theAntiXicPlusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0984 theXicZeroInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0985 theAntiXicZeroInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0986 theOmegacZeroInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0987 theAntiOmegacZeroInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0988 theLambdabInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0989 theAntiLambdabInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0990 theXibZeroInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0991 theAntiXibZeroInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0992 theXibMinusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0993 theAntiXibMinusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0994 theOmegabMinusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0995 theAntiOmegabMinusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold);
0996 theFTFPmodelToBeUsed = theFTFPmodel_aboveThreshold;
0997 if (fPhysicsCase == "FTFP") theFTFPmodelToBeUsed = theFTFPmodel;
0998 theDeuteronInelasticProcess->RegisterMe(theFTFPmodelToBeUsed);
0999 theTritonInelasticProcess->RegisterMe(theFTFPmodelToBeUsed);
1000 theHe3InelasticProcess->RegisterMe(theFTFPmodelToBeUsed);
1001 theAlphaInelasticProcess->RegisterMe(theFTFPmodelToBeUsed);
1002 theIonInelasticProcess->RegisterMe(theFTFPmodelToBeUsed);
1003 }
1004
1005 if (!fPhysicsCaseIsSupported) {
1006 G4cerr << "ERROR: Not supported final-state hadronic inelastic physics case !" << fPhysicsCase
1007 << G4endl << "\t Re-try by choosing one of the following:" << G4endl
1008 << "\t - Hadronic models : BERT, BIC, IonBIC, INCL, FTFP, QGSP" << G4endl
1009 << "\t - \"Physics-list proxies\" : FTFP_BERT (default), FTFP_BERT_ATL, \
1010 QGSP_BERT, QGSP_BIC, FTFP_INCLXX"
1011 << G4endl;
1012 }
1013 }
1014
1015
1016
1017 HadronicGenerator::~HadronicGenerator()
1018 {
1019 fPartTable->DeleteAllParticles();
1020 }
1021
1022
1023
1024 G4bool HadronicGenerator::IsApplicable(const G4String& nameProjectile,
1025 const G4double projectileEnergy) const
1026 {
1027 G4ParticleDefinition* projectileDefinition = fPartTable->FindParticle(nameProjectile);
1028 return IsApplicable(projectileDefinition, projectileEnergy);
1029 }
1030
1031
1032
1033 G4bool HadronicGenerator::IsApplicable(G4ParticleDefinition* projectileDefinition,
1034 const G4double projectileEnergy) const
1035 {
1036 if (projectileDefinition == nullptr) return false;
1037 G4bool isApplicable = true;
1038
1039
1040 if (fPhysicsCase == "BERT") {
1041
1042 if (((projectileDefinition != G4PionMinus::Definition())
1043 && (projectileDefinition != G4PionPlus::Definition())
1044 && (projectileDefinition != G4Proton::Definition())
1045 && (projectileDefinition != G4Neutron::Definition())
1046 && (projectileDefinition != G4Lambda::Definition())
1047 && (projectileDefinition != G4SigmaMinus::Definition())
1048 && (projectileDefinition != G4SigmaPlus::Definition())
1049 && (projectileDefinition != G4XiMinus::Definition())
1050 && (projectileDefinition != G4XiZero::Definition())
1051 && (projectileDefinition != G4OmegaMinus::Definition()))
1052 || (projectileEnergy > 15.0 * CLHEP::GeV))
1053 {
1054 isApplicable = false;
1055 }
1056 }
1057 else if (fPhysicsCase == "QGSP") {
1058
1059 if (projectileEnergy < 2.0 * CLHEP::GeV || projectileDefinition == G4Deuteron::Definition()
1060 || projectileDefinition == G4Triton::Definition()
1061 || projectileDefinition == G4He3::Definition()
1062 || projectileDefinition == G4Alpha::Definition()
1063 || projectileDefinition == G4GenericIon::Definition()
1064 || projectileDefinition == G4AntiDeuteron::Definition()
1065 || projectileDefinition == G4AntiTriton::Definition()
1066 || projectileDefinition == G4AntiHe3::Definition()
1067 || projectileDefinition == G4AntiAlpha::Definition())
1068 {
1069 isApplicable = false;
1070 }
1071 }
1072 else if (fPhysicsCase == "BIC" || fPhysicsCase == "INCL") {
1073
1074
1075
1076 if (((projectileDefinition != G4PionMinus::Definition())
1077 && (projectileDefinition != G4PionPlus::Definition())
1078 && (projectileDefinition != G4Proton::Definition())
1079 && (projectileDefinition != G4Neutron::Definition()))
1080 || (projectileEnergy > 10.0 * CLHEP::GeV))
1081 {
1082 isApplicable = false;
1083 }
1084 }
1085 else if (fPhysicsCase == "IonBIC") {
1086
1087
1088 if (!((projectileDefinition == G4Deuteron::Definition()
1089 && projectileEnergy < 2 * 10.0 * CLHEP::GeV)
1090 || (projectileDefinition == G4Triton::Definition()
1091 && projectileEnergy < 3 * 10.0 * CLHEP::GeV)
1092 || (projectileDefinition == G4He3::Definition()
1093 && projectileEnergy < 3 * 10.0 * CLHEP::GeV)
1094 || (projectileDefinition == G4Alpha::Definition()
1095 && projectileEnergy < 4 * 10.0 * CLHEP::GeV)))
1096 {
1097 isApplicable = false;
1098 }
1099 }
1100 return isApplicable;
1101 }
1102
1103
1104
1105 G4VParticleChange* HadronicGenerator::GenerateInteraction(const G4String& nameProjectile,
1106 const G4double projectileEnergy,
1107 const G4ThreeVector& projectileDirection,
1108 G4Material* targetMaterial)
1109 {
1110 G4ParticleDefinition* projectileDefinition = fPartTable->FindParticle(nameProjectile);
1111 return GenerateInteraction(projectileDefinition, projectileEnergy, projectileDirection,
1112 targetMaterial);
1113 }
1114
1115
1116
1117 G4VParticleChange* HadronicGenerator::GenerateInteraction(
1118 G4ParticleDefinition* projectileDefinition, const G4double projectileEnergy,
1119 const G4ThreeVector& projectileDirection, G4Material* targetMaterial)
1120 {
1121
1122
1123
1124
1125
1126
1127
1128
1129 G4VParticleChange* aChange = nullptr;
1130
1131 if (projectileDefinition == nullptr) {
1132 G4cerr << "ERROR: projectileDefinition is NULL !" << G4endl;
1133 return aChange;
1134 }
1135
1136
1137
1138
1139
1140
1141 if (!IsApplicable(projectileDefinition, projectileEnergy)) {
1142
1143 return aChange;
1144 }
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163 G4DynamicParticle dParticle(projectileDefinition, projectileDirection, projectileEnergy);
1164 const G4double aTime = 0.0;
1165 const G4ThreeVector aPosition = G4ThreeVector(0.0, 0.0, 0.0);
1166 G4Track* gTrack = new G4Track(&dParticle, aTime, aPosition);
1167 G4TouchableHandle fpTouchable(new G4TouchableHistory);
1168 gTrack->SetTouchableHandle(fpTouchable);
1169 G4Step* step = new G4Step;
1170 step->SetTrack(gTrack);
1171 gTrack->SetStep(step);
1172 G4StepPoint* aPoint = new G4StepPoint;
1173 aPoint->SetPosition(aPosition);
1174 aPoint->SetMaterial(targetMaterial);
1175 step->SetPreStepPoint(aPoint);
1176 dParticle.SetKineticEnergy(projectileEnergy);
1177 gTrack->SetStep(step);
1178 gTrack->SetKineticEnergy(projectileEnergy);
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188 G4HadronicProcess* theProcess = nullptr;
1189 G4ParticleDefinition* theProjectileDef = nullptr;
1190 if (projectileDefinition->IsGeneralIon()) {
1191 theProjectileDef = G4GenericIon::Definition();
1192 }
1193 else {
1194 theProjectileDef = projectileDefinition;
1195 }
1196 auto mapIndex = fProcessMap.find(theProjectileDef);
1197 if (mapIndex != fProcessMap.end()) theProcess = mapIndex->second;
1198 if (theProcess != nullptr) {
1199 aChange = theProcess->PostStepDoIt(*gTrack, *step);
1200
1201 }
1202 else {
1203 G4cerr << "ERROR: theProcess is nullptr !" << G4endl;
1204 }
1205 fLastHadronicProcess = theProcess;
1206
1207
1208
1209
1210 return aChange;
1211 }
1212
1213
1214
1215 G4double HadronicGenerator::GetImpactParameter() const
1216 {
1217 G4double impactParameter = -999.0 * fermi;
1218 G4HadronicProcess* hadProcess = GetHadronicProcess();
1219 G4HadronicInteraction* hadInteraction = GetHadronicInteraction();
1220 G4HadronicInteraction* wantedHadInteraction =
1221 const_cast<G4HadronicProcess*>(hadProcess)->GetHadronicModel("FTFP");
1222 if (hadInteraction != nullptr && hadInteraction == wantedHadInteraction) {
1223
1224 G4TheoFSGenerator* theoFSGenerator = dynamic_cast<G4TheoFSGenerator*>(hadInteraction);
1225 if (theoFSGenerator != nullptr) {
1226 const G4FTFModel* ftfModel =
1227 dynamic_cast<const G4FTFModel*>(theoFSGenerator->GetHighEnergyGenerator());
1228 if (ftfModel != nullptr) {
1229
1230
1231 impactParameter = ftfModel->GetImpactParameter();
1232
1233 }
1234 }
1235 }
1236 return impactParameter;
1237 }
1238
1239
1240
1241 G4int HadronicGenerator::GetNumberOfProjectileSpectatorNucleons() const
1242 {
1243 G4double numProjectileSpectatorNucleons = -999;
1244 G4HadronicProcess* hadProcess = GetHadronicProcess();
1245 G4HadronicInteraction* hadInteraction = GetHadronicInteraction();
1246 G4HadronicInteraction* wantedHadInteraction =
1247 const_cast<G4HadronicProcess*>(hadProcess)->GetHadronicModel("FTFP");
1248 if (hadInteraction != nullptr && hadInteraction == wantedHadInteraction) {
1249 G4TheoFSGenerator* theoFSGenerator = dynamic_cast<G4TheoFSGenerator*>(hadInteraction);
1250 if (theoFSGenerator != nullptr) {
1251 const G4FTFModel* ftfModel =
1252 dynamic_cast<const G4FTFModel*>(theoFSGenerator->GetHighEnergyGenerator());
1253 if (ftfModel != nullptr) {
1254 numProjectileSpectatorNucleons = ftfModel->GetNumberOfProjectileSpectatorNucleons();
1255
1256
1257 }
1258 }
1259 }
1260 return numProjectileSpectatorNucleons;
1261 }
1262
1263
1264
1265 G4int HadronicGenerator::GetNumberOfTargetSpectatorNucleons() const
1266 {
1267 G4double numTargetSpectatorNucleons = -999;
1268 G4HadronicProcess* hadProcess = GetHadronicProcess();
1269 G4HadronicInteraction* hadInteraction = GetHadronicInteraction();
1270 G4HadronicInteraction* wantedHadInteraction =
1271 const_cast<G4HadronicProcess*>(hadProcess)->GetHadronicModel("FTFP");
1272 if (hadInteraction != nullptr && hadInteraction == wantedHadInteraction) {
1273 G4TheoFSGenerator* theoFSGenerator = dynamic_cast<G4TheoFSGenerator*>(hadInteraction);
1274 if (theoFSGenerator != nullptr) {
1275 const G4FTFModel* ftfModel =
1276 dynamic_cast<const G4FTFModel*>(theoFSGenerator->GetHighEnergyGenerator());
1277 if (ftfModel != nullptr) {
1278 numTargetSpectatorNucleons = ftfModel->GetNumberOfTargetSpectatorNucleons();
1279
1280 }
1281 }
1282 }
1283 return numTargetSpectatorNucleons;
1284 }
1285
1286
1287
1288 G4int HadronicGenerator::GetNumberOfNNcollisions() const
1289 {
1290 G4double numNNcollisions = -999;
1291 G4HadronicProcess* hadProcess = GetHadronicProcess();
1292 G4HadronicInteraction* hadInteraction = GetHadronicInteraction();
1293 G4HadronicInteraction* wantedHadInteraction =
1294 const_cast<G4HadronicProcess*>(hadProcess)->GetHadronicModel("FTFP");
1295 if (hadInteraction != nullptr && hadInteraction == wantedHadInteraction) {
1296 G4TheoFSGenerator* theoFSGenerator = dynamic_cast<G4TheoFSGenerator*>(hadInteraction);
1297 if (theoFSGenerator != nullptr) {
1298 const G4FTFModel* ftfModel =
1299 dynamic_cast<const G4FTFModel*>(theoFSGenerator->GetHighEnergyGenerator());
1300 if (ftfModel != nullptr) {
1301 numNNcollisions = ftfModel->GetNumberOfNNcollisions();
1302
1303 }
1304 }
1305 }
1306 return numNNcollisions;
1307 }
1308
1309