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