Geant4 Cross Reference |
1 // 2 // ******************************************************************** 3 // * License and Disclaimer * 4 // * * 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. * 10 // * * 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitation of liability. * 17 // * * 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************************************** 25 // 26 /// \file HadronicGenerator.cc 27 /// \brief Implementation of the HadronicGenerator class 28 // 29 //------------------------------------------------------------------------ 30 // Class: HadronicGenerator 31 // Author: Alberto Ribon (CERN EP/SFT), May 2020 32 // Modified: G. Hugo, 8 December 2022 33 //------------------------------------------------------------------------ 34 35 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 36 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 37 38 #include "HadronicGenerator.hh" 39 40 #include "G4AblaInterface.hh" 41 #include "G4Alpha.hh" 42 #include "G4AntiAlpha.hh" 43 #include "G4AntiBMesonZero.hh" 44 #include "G4AntiBsMesonZero.hh" 45 #include "G4AntiDMesonZero.hh" 46 #include "G4AntiDeuteron.hh" 47 #include "G4AntiHe3.hh" 48 #include "G4AntiLambda.hh" 49 #include "G4AntiLambdab.hh" 50 #include "G4AntiLambdacPlus.hh" 51 #include "G4AntiNeutron.hh" 52 #include "G4AntiOmegaMinus.hh" 53 #include "G4AntiOmegabMinus.hh" 54 #include "G4AntiOmegacZero.hh" 55 #include "G4AntiProton.hh" 56 #include "G4AntiSigmaMinus.hh" 57 #include "G4AntiSigmaPlus.hh" 58 #include "G4AntiSigmaZero.hh" 59 #include "G4AntiTriton.hh" 60 #include "G4AntiXiMinus.hh" 61 #include "G4AntiXiZero.hh" 62 #include "G4AntiXibMinus.hh" 63 #include "G4AntiXibZero.hh" 64 #include "G4AntiXicPlus.hh" 65 #include "G4AntiXicZero.hh" 66 #include "G4BGGNucleonInelasticXS.hh" 67 #include "G4BGGPionInelasticXS.hh" 68 #include "G4BMesonMinus.hh" 69 #include "G4BMesonPlus.hh" 70 #include "G4BMesonZero.hh" 71 #include "G4BcMesonMinus.hh" 72 #include "G4BcMesonPlus.hh" 73 #include "G4BinaryCascade.hh" 74 #include "G4BinaryLightIonReaction.hh" 75 #include "G4Box.hh" 76 #include "G4BsMesonZero.hh" 77 #include "G4CascadeInterface.hh" 78 #include "G4ChipsHyperonInelasticXS.hh" 79 #include "G4ComponentAntiNuclNuclearXS.hh" 80 #include "G4ComponentGGHadronNucleusXsc.hh" 81 #include "G4ComponentGGNuclNuclXsc.hh" 82 #include "G4CrossSectionInelastic.hh" 83 #include "G4DMesonMinus.hh" 84 #include "G4DMesonPlus.hh" 85 #include "G4DMesonZero.hh" 86 #include "G4DecayPhysics.hh" 87 #include "G4Deuteron.hh" 88 #include "G4DsMesonMinus.hh" 89 #include "G4DsMesonPlus.hh" 90 #include "G4DynamicParticle.hh" 91 #include "G4ExcitationHandler.hh" 92 #include "G4ExcitedStringDecay.hh" 93 #include "G4FTFModel.hh" 94 #include "G4GeneratorPrecompoundInterface.hh" 95 #include "G4GenericIon.hh" 96 #include "G4HadronInelasticProcess.hh" 97 #include "G4HadronicParameters.hh" 98 #include "G4He3.hh" 99 #include "G4INCLXXInterface.hh" 100 #include "G4IonTable.hh" 101 #include "G4KaonMinus.hh" 102 #include "G4KaonPlus.hh" 103 #include "G4KaonZeroLong.hh" 104 #include "G4KaonZeroShort.hh" 105 #include "G4Lambda.hh" 106 #include "G4Lambdab.hh" 107 #include "G4LambdacPlus.hh" 108 #include "G4LundStringFragmentation.hh" 109 #include "G4Material.hh" 110 #include "G4Neutron.hh" 111 #include "G4NeutronInelasticXS.hh" 112 #include "G4OmegaMinus.hh" 113 #include "G4OmegabMinus.hh" 114 #include "G4OmegacZero.hh" 115 #include "G4PVPlacement.hh" 116 #include "G4ParticleTable.hh" 117 #include "G4PhysicalConstants.hh" 118 #include "G4PionMinus.hh" 119 #include "G4PionPlus.hh" 120 #include "G4PreCompoundModel.hh" 121 #include "G4ProcessManager.hh" 122 #include "G4Proton.hh" 123 #include "G4QGSMFragmentation.hh" 124 #include "G4QGSModel.hh" 125 #include "G4QGSParticipants.hh" 126 #include "G4QuasiElasticChannel.hh" 127 #include "G4SigmaMinus.hh" 128 #include "G4SigmaPlus.hh" 129 #include "G4SigmaZero.hh" 130 #include "G4StateManager.hh" 131 #include "G4Step.hh" 132 #include "G4SystemOfUnits.hh" 133 #include "G4TheoFSGenerator.hh" 134 #include "G4TouchableHistory.hh" 135 #include "G4TransportationManager.hh" 136 #include "G4Triton.hh" 137 #include "G4UnitsTable.hh" 138 #include "G4VCrossSectionDataSet.hh" 139 #include "G4VParticleChange.hh" 140 #include "G4XiMinus.hh" 141 #include "G4XiZero.hh" 142 #include "G4XibMinus.hh" 143 #include "G4XibZero.hh" 144 #include "G4XicPlus.hh" 145 #include "G4XicZero.hh" 146 #include "G4ios.hh" 147 #include "globals.hh" 148 149 #include <iomanip> 150 151 // FLUKA ADAPTATION 152 #ifdef G4_USE_FLUKA 153 # include "FLUKAInelasticScatteringXS.hh" 154 # include "FLUKANuclearInelasticModel.hh" 155 # include "FLUKAParticleTable.hh" 156 # include "fluka_interface.hh" 157 #endif 158 159 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 160 161 HadronicGenerator::HadronicGenerator(const G4String physicsCase) 162 : fPhysicsCase(physicsCase), 163 fPhysicsCaseIsSupported(false), 164 fLastHadronicProcess(nullptr), 165 fPartTable(nullptr) 166 { 167 // The constructor set-ups all the particles, models, cross sections and 168 // hadronic inelastic processes. 169 // This should be done only once for each application. 170 // In the case of a multi-threaded application using this class, 171 // the constructor should be invoked for each thread, 172 // i.e. one instance of the class should be kept per thread. 173 // The particles and processes that are created in this constructor 174 // will then be used by the method GenerateInteraction at each interaction. 175 // Notes: 176 // - Neither the hadronic models nor the cross sections are used directly 177 // by the method GenerateInteraction, but they are associated to the 178 // hadronic processes and used by Geant4 to simulate the collision; 179 // - Although the class generates only final states, but not free mean paths, 180 // inelastic hadron-nuclear cross sections are needed by Geant4 to sample 181 // the target nucleus from the target material. 182 183 // Definition of particles 184 G4GenericIon* gion = G4GenericIon::Definition(); 185 gion->SetProcessManager(new G4ProcessManager(gion)); 186 G4DecayPhysics* decays = new G4DecayPhysics; 187 decays->ConstructParticle(); 188 fPartTable = G4ParticleTable::GetParticleTable(); 189 fPartTable->SetReadiness(); 190 G4IonTable* ions = fPartTable->GetIonTable(); 191 ions->CreateAllIon(); 192 ions->CreateAllIsomer(); 193 194 // FLUKA ADAPTATION 195 if (fPhysicsCase == "CFLUKAHI") { 196 #ifdef G4_USE_FLUKA 197 198 // IMPORTANT: Initialize the FLUKA interface here. 199 // Both activation switches should be set to TRUE to provide the most comprehensive results. 200 // NB: COMPARISON WITH G4 DOES NOT SEEM MEANINGFUL 201 // WHEN COALESCENCE IS ACTIVATED IN BOTH FLUKA AND G4. 202 // Freedom to choose & see the effect of these switches is hence provided here. 203 const G4bool activateCoalescence = true; 204 const G4bool activateHeavyFragmentsEvaporation = true; 205 fluka_interface::initialize(activateCoalescence, activateHeavyFragmentsEvaporation); 206 207 // Initialize FLUKA <-> G4 particles conversions tables. 208 fluka_particle_table::initialize(); 209 #else 210 G4Exception("HadronicGenerator.cc", "Wrong compilation mode.", FatalException, 211 "Requested G4_HP_CernFLUKAHadronInelastic physics list.\n" 212 "This requires COMPILATION in FLUKA mode.\n" 213 "Please fully recompile the example with G4_USE_FLUKA=yes.\n" 214 "For example:\n" 215 "source geant4/examples/extended/hadronic/FlukaCern/FlukaInterface/" 216 "env_FLUKA_G4_interface.sh\n" 217 "cd geant4/examples/extended/hadronic/FlukaCern/ProcessLevel/CrossSection/\n" 218 "mkdir build\n" 219 "cd build\n" 220 "cmake -DGeant4_DIR=your_path_to_geant4 -DG4_USE_FLUKA=1 .. \n" 221 "make -j8 G4_USE_FLUKA=1\n" 222 "NB: First time use of FLUKA interface:\n" 223 "Do not forget to first compile the FLUKA interface itself.\n" 224 "For example: cd geant4/examples/extended/hadronic/FlukaCern/FlukaInterface/ " 225 "&& make interface && make env\n" 226 "FlukaInterface/env_FLUKA_G4_interface.sh then needs to be sourced\n" 227 "in whichever terminal you want to use the FLUKA interface.\n"); 228 #endif 229 } 230 231 // Build BERT model 232 G4CascadeInterface* theBERTmodel = new G4CascadeInterface; 233 234 // Build BIC model 235 G4BinaryCascade* theBICmodel = new G4BinaryCascade; 236 G4PreCompoundModel* thePreEquilib = new G4PreCompoundModel(new G4ExcitationHandler); 237 theBICmodel->SetDeExcitation(thePreEquilib); 238 239 // Build BinaryLightIon model 240 G4PreCompoundModel* thePreEquilibBis = new G4PreCompoundModel(new G4ExcitationHandler); 241 G4BinaryLightIonReaction* theIonBICmodel = new G4BinaryLightIonReaction(thePreEquilibBis); 242 243 // Build the INCL model 244 G4INCLXXInterface* theINCLmodel = new G4INCLXXInterface; 245 const G4bool useAblaDeExcitation = false; // By default INCL uses Preco: set "true" to use 246 // ABLA DeExcitation 247 if (theINCLmodel && useAblaDeExcitation) { 248 G4AblaInterface* theAblaInterface = new G4AblaInterface; 249 theINCLmodel->SetDeExcitation(theAblaInterface); 250 } 251 252 // Build the FTFP model (FTF/Preco) : 4 instances with different kinetic energy intervals. 253 // (Notice that these kinetic energy intervals are applied per nucleons, so they are fine 254 // for all types of hadron and ion projectile). 255 // Model instance without energy constraint. 256 // (Used for the case of FTFP model, and for light anti-ions in all physics lists.) 257 G4TheoFSGenerator* theFTFPmodel = new G4TheoFSGenerator("FTFP"); 258 theFTFPmodel->SetMaxEnergy(G4HadronicParameters::Instance()->GetMaxEnergy()); 259 G4GeneratorPrecompoundInterface* theCascade = new G4GeneratorPrecompoundInterface; 260 theCascade->SetDeExcitation(thePreEquilib); 261 theFTFPmodel->SetTransport(theCascade); 262 G4LundStringFragmentation* theLundFragmentation = new G4LundStringFragmentation; 263 G4ExcitedStringDecay* theStringDecay = new G4ExcitedStringDecay(theLundFragmentation); 264 G4FTFModel* theStringModel = new G4FTFModel; 265 theStringModel->SetFragmentationModel(theStringDecay); 266 267 // If the following line is set, then the square of the impact parameter is sampled 268 // randomly from a flat distribution in the range [ Bmin*Bmin, Bmax*Bmax ] 269 // theStringModel->SetBminBmax( 0.0, 2.0*fermi ); 270 //***LOOKHERE*** CHOOSE IMPACT PARAMETER MIN & MAX 271 272 theFTFPmodel->SetHighEnergyGenerator(theStringModel); 273 // Model instance with constraint to be above a kinetic energy threshold. 274 // (Used for ions in all physics lists, and, in the case of non-QGS-based physics lists, 275 // also for pions, kaons, nucleons and hyperons.) 276 G4TheoFSGenerator* theFTFPmodel_aboveThreshold = new G4TheoFSGenerator("FTFP"); 277 theFTFPmodel_aboveThreshold->SetMaxEnergy(G4HadronicParameters::Instance()->GetMaxEnergy()); 278 theFTFPmodel_aboveThreshold->SetTransport(theCascade); 279 theFTFPmodel_aboveThreshold->SetHighEnergyGenerator(theStringModel); 280 // Model instance with constraint to be within two kinetic energy thresholds. 281 // (Used in the case of QGS-based physics lists for pions, kaons, nucleons and hyperons.) 282 G4TheoFSGenerator* theFTFPmodel_constrained = new G4TheoFSGenerator("FTFP"); 283 theFTFPmodel_constrained->SetMaxEnergy(G4HadronicParameters::Instance()->GetMaxEnergy()); 284 theFTFPmodel_constrained->SetTransport(theCascade); 285 theFTFPmodel_constrained->SetHighEnergyGenerator(theStringModel); 286 // Model instance to be used down to zero kinetic energy, with eventual constraint 287 // - in the case of QGS-based physics lists - to be below a kinetic energy threshold. 288 // (Used for anti-baryons, anti-hyperons, and charmed and bottom hadrons.) 289 G4TheoFSGenerator* theFTFPmodel_belowThreshold = new G4TheoFSGenerator("FTFP"); 290 theFTFPmodel_belowThreshold->SetMaxEnergy(G4HadronicParameters::Instance()->GetMaxEnergy()); 291 theFTFPmodel_belowThreshold->SetTransport(theCascade); 292 theFTFPmodel_belowThreshold->SetHighEnergyGenerator(theStringModel); 293 294 // Build the QGSP model (QGS/Preco) 295 G4TheoFSGenerator* theQGSPmodel = new G4TheoFSGenerator("QGSP"); 296 theQGSPmodel->SetMaxEnergy(G4HadronicParameters::Instance()->GetMaxEnergy()); 297 theQGSPmodel->SetTransport(theCascade); 298 G4QGSMFragmentation* theQgsmFragmentation = new G4QGSMFragmentation; 299 G4ExcitedStringDecay* theQgsmStringDecay = new G4ExcitedStringDecay(theQgsmFragmentation); 300 G4VPartonStringModel* theQgsmStringModel = new G4QGSModel<G4QGSParticipants>; 301 theQgsmStringModel->SetFragmentationModel(theQgsmStringDecay); 302 theQGSPmodel->SetHighEnergyGenerator(theQgsmStringModel); 303 G4QuasiElasticChannel* theQuasiElastic = new G4QuasiElasticChannel; // QGSP uses quasi-elastic 304 theQGSPmodel->SetQuasiElasticChannel(theQuasiElastic); 305 306 // For the case of "physics-list proxies", select the energy range for each hadronic model. 307 // Note: the transition energy between hadronic models vary between physics lists, 308 // type of hadrons, and version of Geant4. Here, for simplicity, we use an uniform 309 // energy transition for all types of hadrons and regardless of the Geant4 version; 310 // moreover, for "FTFP_INCLXX" we use a different energy transition range 311 // between FTFP and INCL than in the real physics list. 312 if (fPhysicsCase == "FTFP_BERT_ATL" || fPhysicsCase == "FTFP_BERT" 313 || fPhysicsCase == "FTFP_INCLXX" || fPhysicsCase == "QGSP_BERT" || fPhysicsCase == "QGSP_BIC") 314 { 315 const G4double ftfpMinE = G4HadronicParameters::Instance()->GetMinEnergyTransitionFTF_Cascade(); 316 const G4double bertMaxE = G4HadronicParameters::Instance()->GetMaxEnergyTransitionFTF_Cascade(); 317 const G4double ftfpMinE_ATL = 9.0 * CLHEP::GeV; 318 const G4double bertMaxE_ATL = 12.0 * CLHEP::GeV; 319 const G4double ftfpMaxE = G4HadronicParameters::Instance()->GetMaxEnergyTransitionQGS_FTF(); 320 const G4double qgspMinE = G4HadronicParameters::Instance()->GetMinEnergyTransitionQGS_FTF(); 321 theFTFPmodel->SetMinEnergy(0.0); 322 theFTFPmodel_belowThreshold->SetMinEnergy(0.0); 323 if (fPhysicsCase == "FTFP_BERT_ATL") { 324 theBERTmodel->SetMaxEnergy(bertMaxE_ATL); 325 theIonBICmodel->SetMaxEnergy(bertMaxE_ATL); 326 theFTFPmodel_aboveThreshold->SetMinEnergy(ftfpMinE_ATL); 327 theFTFPmodel_constrained->SetMinEnergy(ftfpMinE_ATL); 328 } 329 else { 330 theBERTmodel->SetMaxEnergy(bertMaxE); 331 theIonBICmodel->SetMaxEnergy(bertMaxE); 332 theFTFPmodel_aboveThreshold->SetMinEnergy(ftfpMinE); 333 theFTFPmodel_constrained->SetMinEnergy(ftfpMinE); 334 } 335 if (fPhysicsCase == "FTFP_INCLXX") { 336 theINCLmodel->SetMaxEnergy(bertMaxE); 337 } 338 if (fPhysicsCase == "QGSP_BERT" || fPhysicsCase == "QGSP_BIC") { 339 theFTFPmodel_constrained->SetMaxEnergy(ftfpMaxE); 340 theFTFPmodel_belowThreshold->SetMaxEnergy(ftfpMaxE); 341 theQGSPmodel->SetMinEnergy(qgspMinE); 342 theBICmodel->SetMaxEnergy(bertMaxE); 343 } 344 } 345 346 // Cross sections (needed by Geant4 to sample the target nucleus from the target material) 347 G4VCrossSectionDataSet* thePionMinusXSdata = new G4BGGPionInelasticXS(G4PionMinus::Definition()); 348 thePionMinusXSdata->BuildPhysicsTable(*(G4PionMinus::Definition())); 349 G4VCrossSectionDataSet* thePionPlusXSdata = new G4BGGPionInelasticXS(G4PionPlus::Definition()); 350 thePionPlusXSdata->BuildPhysicsTable(*(G4PionPlus::Definition())); 351 G4VCrossSectionDataSet* theKaonXSdata = 352 new G4CrossSectionInelastic(new G4ComponentGGHadronNucleusXsc); 353 theKaonXSdata->BuildPhysicsTable(*(G4KaonMinus::Definition())); 354 theKaonXSdata->BuildPhysicsTable(*(G4KaonPlus::Definition())); 355 theKaonXSdata->BuildPhysicsTable(*(G4KaonZeroLong::Definition())); 356 theKaonXSdata->BuildPhysicsTable(*(G4KaonZeroShort::Definition())); 357 G4VCrossSectionDataSet* theProtonXSdata = new G4BGGNucleonInelasticXS(G4Proton::Proton()); 358 theProtonXSdata->BuildPhysicsTable(*(G4Proton::Definition())); 359 G4VCrossSectionDataSet* theNeutronXSdata = new G4NeutronInelasticXS; 360 theNeutronXSdata->BuildPhysicsTable(*(G4Neutron::Definition())); 361 // For hyperon and anti-hyperons we can use either Chips or, for G4 >= 10.5, 362 // Glauber-Gribov cross sections 363 // G4VCrossSectionDataSet* theHyperonsXSdata = new G4ChipsHyperonInelasticXS; 364 G4VCrossSectionDataSet* theHyperonsXSdata = 365 new G4CrossSectionInelastic(new G4ComponentGGHadronNucleusXsc); 366 G4VCrossSectionDataSet* theAntibaryonsXSdata = 367 new G4CrossSectionInelastic(new G4ComponentAntiNuclNuclearXS); 368 G4VCrossSectionDataSet* theNuclNuclXSdata = 369 new G4CrossSectionInelastic(new G4ComponentGGNuclNuclXsc); 370 371 // Set up inelastic processes : store them in a map (with particle definition as key) 372 // for convenience 373 typedef std::pair<G4ParticleDefinition*, G4HadronicProcess*> ProcessPair; 374 G4HadronicProcess* thePionMinusInelasticProcess = 375 new G4HadronInelasticProcess("pi-Inelastic", G4PionMinus::Definition()); 376 fProcessMap.insert(ProcessPair(G4PionMinus::Definition(), thePionMinusInelasticProcess)); 377 G4HadronicProcess* thePionPlusInelasticProcess = 378 new G4HadronInelasticProcess("pi+Inelastic", G4PionPlus::Definition()); 379 fProcessMap.insert(ProcessPair(G4PionPlus::Definition(), thePionPlusInelasticProcess)); 380 G4HadronicProcess* theKaonMinusInelasticProcess = 381 new G4HadronInelasticProcess("kaon-Inelastic", G4KaonMinus::Definition()); 382 fProcessMap.insert(ProcessPair(G4KaonMinus::Definition(), theKaonMinusInelasticProcess)); 383 G4HadronicProcess* theKaonPlusInelasticProcess = 384 new G4HadronInelasticProcess("kaon+Inelastic", G4KaonPlus::Definition()); 385 fProcessMap.insert(ProcessPair(G4KaonPlus::Definition(), theKaonPlusInelasticProcess)); 386 G4HadronicProcess* theKaonZeroLInelasticProcess = 387 new G4HadronInelasticProcess("kaon0LInelastic", G4KaonZeroLong::Definition()); 388 fProcessMap.insert(ProcessPair(G4KaonZeroLong::Definition(), theKaonZeroLInelasticProcess)); 389 G4HadronicProcess* theKaonZeroSInelasticProcess = 390 new G4HadronInelasticProcess("kaon0SInelastic", G4KaonZeroShort::Definition()); 391 fProcessMap.insert(ProcessPair(G4KaonZeroShort::Definition(), theKaonZeroSInelasticProcess)); 392 G4HadronicProcess* theProtonInelasticProcess = 393 new G4HadronInelasticProcess("protonInelastic", G4Proton::Definition()); 394 fProcessMap.insert(ProcessPair(G4Proton::Definition(), theProtonInelasticProcess)); 395 G4HadronicProcess* theNeutronInelasticProcess = 396 new G4HadronInelasticProcess("neutronInelastic", G4Neutron::Definition()); 397 fProcessMap.insert(ProcessPair(G4Neutron::Definition(), theNeutronInelasticProcess)); 398 G4HadronicProcess* theDeuteronInelasticProcess = 399 new G4HadronInelasticProcess("dInelastic", G4Deuteron::Definition()); 400 fProcessMap.insert(ProcessPair(G4Deuteron::Definition(), theDeuteronInelasticProcess)); 401 G4HadronicProcess* theTritonInelasticProcess = 402 new G4HadronInelasticProcess("tInelastic", G4Triton::Definition()); 403 fProcessMap.insert(ProcessPair(G4Triton::Definition(), theTritonInelasticProcess)); 404 G4HadronicProcess* theHe3InelasticProcess = 405 new G4HadronInelasticProcess("he3Inelastic", G4He3::Definition()); 406 fProcessMap.insert(ProcessPair(G4He3::Definition(), theHe3InelasticProcess)); 407 G4HadronicProcess* theAlphaInelasticProcess = 408 new G4HadronInelasticProcess("alphaInelastic", G4Alpha::Definition()); 409 fProcessMap.insert(ProcessPair(G4Alpha::Definition(), theAlphaInelasticProcess)); 410 G4HadronicProcess* theIonInelasticProcess = 411 new G4HadronInelasticProcess("ionInelastic", G4GenericIon::Definition()); 412 fProcessMap.insert(ProcessPair(G4GenericIon::Definition(), theIonInelasticProcess)); 413 G4HadronicProcess* theLambdaInelasticProcess = 414 new G4HadronInelasticProcess("lambdaInelastic", G4Lambda::Definition()); 415 fProcessMap.insert(ProcessPair(G4Lambda::Definition(), theLambdaInelasticProcess)); 416 G4HadronicProcess* theSigmaMinusInelasticProcess = 417 new G4HadronInelasticProcess("sigma-Inelastic", G4SigmaMinus::Definition()); 418 fProcessMap.insert(ProcessPair(G4SigmaMinus::Definition(), theSigmaMinusInelasticProcess)); 419 G4HadronicProcess* theSigmaPlusInelasticProcess = 420 new G4HadronInelasticProcess("sigma+Inelastic", G4SigmaPlus::Definition()); 421 fProcessMap.insert(ProcessPair(G4SigmaPlus::Definition(), theSigmaPlusInelasticProcess)); 422 G4HadronicProcess* theXiMinusInelasticProcess = 423 new G4HadronInelasticProcess("xi-Inelastic", G4XiMinus::Definition()); 424 fProcessMap.insert(ProcessPair(G4XiMinus::Definition(), theXiMinusInelasticProcess)); 425 G4HadronicProcess* theXiZeroInelasticProcess = 426 new G4HadronInelasticProcess("xi0Inelastic", G4XiZero::Definition()); 427 fProcessMap.insert(ProcessPair(G4XiZero::Definition(), theXiZeroInelasticProcess)); 428 G4HadronicProcess* theOmegaMinusInelasticProcess = 429 new G4HadronInelasticProcess("omega-Inelastic", G4OmegaMinus::Definition()); 430 fProcessMap.insert(ProcessPair(G4OmegaMinus::Definition(), theOmegaMinusInelasticProcess)); 431 G4HadronicProcess* theAntiProtonInelasticProcess = 432 new G4HadronInelasticProcess("anti_protonInelastic", G4AntiProton::Definition()); 433 fProcessMap.insert(ProcessPair(G4AntiProton::Definition(), theAntiProtonInelasticProcess)); 434 G4HadronicProcess* theAntiNeutronInelasticProcess = 435 new G4HadronInelasticProcess("anti_neutronInelastic", G4AntiNeutron::Definition()); 436 fProcessMap.insert(ProcessPair(G4AntiNeutron::Definition(), theAntiNeutronInelasticProcess)); 437 G4HadronicProcess* theAntiDeuteronInelasticProcess = 438 new G4HadronInelasticProcess("anti_deuteronInelastic", G4AntiDeuteron::Definition()); 439 fProcessMap.insert(ProcessPair(G4AntiDeuteron::Definition(), theAntiDeuteronInelasticProcess)); 440 G4HadronicProcess* theAntiTritonInelasticProcess = 441 new G4HadronInelasticProcess("anti_tritonInelastic", G4AntiTriton::Definition()); 442 fProcessMap.insert(ProcessPair(G4AntiTriton::Definition(), theAntiTritonInelasticProcess)); 443 G4HadronicProcess* theAntiHe3InelasticProcess = 444 new G4HadronInelasticProcess("anti_He3Inelastic", G4AntiHe3::Definition()); 445 fProcessMap.insert(ProcessPair(G4AntiHe3::Definition(), theAntiHe3InelasticProcess)); 446 G4HadronicProcess* theAntiAlphaInelasticProcess = 447 new G4HadronInelasticProcess("anti_alphaInelastic", G4AntiAlpha::Definition()); 448 fProcessMap.insert(ProcessPair(G4AntiAlpha::Definition(), theAntiAlphaInelasticProcess)); 449 G4HadronicProcess* theAntiLambdaInelasticProcess = 450 new G4HadronInelasticProcess("anti-lambdaInelastic", G4AntiLambda::Definition()); 451 fProcessMap.insert(ProcessPair(G4AntiLambda::Definition(), theAntiLambdaInelasticProcess)); 452 G4HadronicProcess* theAntiSigmaMinusInelasticProcess = 453 new G4HadronInelasticProcess("anti_sigma-Inelastic", G4AntiSigmaMinus::Definition()); 454 fProcessMap.insert( 455 ProcessPair(G4AntiSigmaMinus::Definition(), theAntiSigmaMinusInelasticProcess)); 456 G4HadronicProcess* theAntiSigmaPlusInelasticProcess = 457 new G4HadronInelasticProcess("anti_sigma+Inelastic", G4AntiSigmaPlus::Definition()); 458 fProcessMap.insert(ProcessPair(G4AntiSigmaPlus::Definition(), theAntiSigmaPlusInelasticProcess)); 459 G4HadronicProcess* theAntiXiMinusInelasticProcess = 460 new G4HadronInelasticProcess("anti_xi-Inelastic", G4AntiXiMinus::Definition()); 461 fProcessMap.insert(ProcessPair(G4AntiXiMinus::Definition(), theAntiXiMinusInelasticProcess)); 462 G4HadronicProcess* theAntiXiZeroInelasticProcess = 463 new G4HadronInelasticProcess("anti_xi0Inelastic", G4AntiXiZero::Definition()); 464 fProcessMap.insert(ProcessPair(G4AntiXiZero::Definition(), theAntiXiZeroInelasticProcess)); 465 G4HadronicProcess* theAntiOmegaMinusInelasticProcess = 466 new G4HadronInelasticProcess("anti_omega-Inelastic", G4AntiOmegaMinus::Definition()); 467 fProcessMap.insert( 468 ProcessPair(G4AntiOmegaMinus::Definition(), theAntiOmegaMinusInelasticProcess)); 469 470 G4HadronicProcess* theDPlusInelasticProcess = 471 new G4HadronInelasticProcess("D+Inelastic", G4DMesonPlus::Definition()); 472 fProcessMap.insert(ProcessPair(G4DMesonPlus::Definition(), theDPlusInelasticProcess)); 473 G4HadronicProcess* theDMinusInelasticProcess = 474 new G4HadronInelasticProcess("D-Inelastic", G4DMesonMinus::Definition()); 475 fProcessMap.insert(ProcessPair(G4DMesonMinus::Definition(), theDMinusInelasticProcess)); 476 G4HadronicProcess* theDZeroInelasticProcess = 477 new G4HadronInelasticProcess("D0Inelastic", G4DMesonZero::Definition()); 478 fProcessMap.insert(ProcessPair(G4DMesonZero::Definition(), theDZeroInelasticProcess)); 479 G4HadronicProcess* theAntiDZeroInelasticProcess = 480 new G4HadronInelasticProcess("anti_D0Inelastic", G4AntiDMesonZero::Definition()); 481 fProcessMap.insert(ProcessPair(G4AntiDMesonZero::Definition(), theAntiDZeroInelasticProcess)); 482 G4HadronicProcess* theDsPlusInelasticProcess = 483 new G4HadronInelasticProcess("Ds+Inelastic", G4DsMesonPlus::Definition()); 484 fProcessMap.insert(ProcessPair(G4DsMesonPlus::Definition(), theDsPlusInelasticProcess)); 485 G4HadronicProcess* theDsMinusInelasticProcess = 486 new G4HadronInelasticProcess("Ds-Inelastic", G4DsMesonMinus::Definition()); 487 fProcessMap.insert(ProcessPair(G4DsMesonMinus::Definition(), theDsMinusInelasticProcess)); 488 G4HadronicProcess* theBPlusInelasticProcess = 489 new G4HadronInelasticProcess("B+Inelastic", G4BMesonPlus::Definition()); 490 fProcessMap.insert(ProcessPair(G4BMesonPlus::Definition(), theBPlusInelasticProcess)); 491 G4HadronicProcess* theBMinusInelasticProcess = 492 new G4HadronInelasticProcess("B-Inelastic", G4BMesonMinus::Definition()); 493 fProcessMap.insert(ProcessPair(G4BMesonMinus::Definition(), theBMinusInelasticProcess)); 494 G4HadronicProcess* theBZeroInelasticProcess = 495 new G4HadronInelasticProcess("B0Inelastic", G4BMesonZero::Definition()); 496 fProcessMap.insert(ProcessPair(G4BMesonZero::Definition(), theBZeroInelasticProcess)); 497 G4HadronicProcess* theAntiBZeroInelasticProcess = 498 new G4HadronInelasticProcess("anti_B0Inelastic", G4AntiBMesonZero::Definition()); 499 fProcessMap.insert(ProcessPair(G4AntiBMesonZero::Definition(), theAntiBZeroInelasticProcess)); 500 G4HadronicProcess* theBsZeroInelasticProcess = 501 new G4HadronInelasticProcess("Bs0Inelastic", G4BsMesonZero::Definition()); 502 fProcessMap.insert(ProcessPair(G4BsMesonZero::Definition(), theBsZeroInelasticProcess)); 503 G4HadronicProcess* theAntiBsZeroInelasticProcess = 504 new G4HadronInelasticProcess("anti_Bs0Inelastic", G4AntiBsMesonZero::Definition()); 505 fProcessMap.insert(ProcessPair(G4AntiBsMesonZero::Definition(), theAntiBsZeroInelasticProcess)); 506 G4HadronicProcess* theBcPlusInelasticProcess = 507 new G4HadronInelasticProcess("Bc+Inelastic", G4BcMesonPlus::Definition()); 508 fProcessMap.insert(ProcessPair(G4BcMesonPlus::Definition(), theBcPlusInelasticProcess)); 509 G4HadronicProcess* theBcMinusInelasticProcess = 510 new G4HadronInelasticProcess("Bc-Inelastic", G4BcMesonMinus::Definition()); 511 fProcessMap.insert(ProcessPair(G4BcMesonMinus::Definition(), theBcMinusInelasticProcess)); 512 G4HadronicProcess* theLambdacPlusInelasticProcess = 513 new G4HadronInelasticProcess("lambda_c+Inelastic", G4LambdacPlus::Definition()); 514 fProcessMap.insert(ProcessPair(G4LambdacPlus::Definition(), theLambdacPlusInelasticProcess)); 515 G4HadronicProcess* theAntiLambdacPlusInelasticProcess = 516 new G4HadronInelasticProcess("anti_lambda_c+Inelastic", G4AntiLambdacPlus::Definition()); 517 fProcessMap.insert( 518 ProcessPair(G4AntiLambdacPlus::Definition(), theAntiLambdacPlusInelasticProcess)); 519 G4HadronicProcess* theXicPlusInelasticProcess = 520 new G4HadronInelasticProcess("xi_c+Inelastic", G4XicPlus::Definition()); 521 fProcessMap.insert(ProcessPair(G4XicPlus::Definition(), theXicPlusInelasticProcess)); 522 G4HadronicProcess* theAntiXicPlusInelasticProcess = 523 new G4HadronInelasticProcess("anti_xi_c+Inelastic", G4AntiXicPlus::Definition()); 524 fProcessMap.insert(ProcessPair(G4AntiXicPlus::Definition(), theAntiXicPlusInelasticProcess)); 525 G4HadronicProcess* theXicZeroInelasticProcess = 526 new G4HadronInelasticProcess("xi_c0Inelastic", G4XicZero::Definition()); 527 fProcessMap.insert(ProcessPair(G4XicZero::Definition(), theXicZeroInelasticProcess)); 528 G4HadronicProcess* theAntiXicZeroInelasticProcess = 529 new G4HadronInelasticProcess("anti_xi_c0Inelastic", G4AntiXicZero::Definition()); 530 fProcessMap.insert(ProcessPair(G4AntiXicZero::Definition(), theAntiXicZeroInelasticProcess)); 531 G4HadronicProcess* theOmegacZeroInelasticProcess = 532 new G4HadronInelasticProcess("omega_c0Inelastic", G4OmegacZero::Definition()); 533 fProcessMap.insert(ProcessPair(G4OmegacZero::Definition(), theOmegacZeroInelasticProcess)); 534 G4HadronicProcess* theAntiOmegacZeroInelasticProcess = 535 new G4HadronInelasticProcess("anti_omega_c0Inelastic", G4AntiOmegacZero::Definition()); 536 fProcessMap.insert( 537 ProcessPair(G4AntiOmegacZero::Definition(), theAntiOmegacZeroInelasticProcess)); 538 G4HadronicProcess* theLambdabInelasticProcess = 539 new G4HadronInelasticProcess("lambda_bInelastic", G4Lambdab::Definition()); 540 fProcessMap.insert(ProcessPair(G4Lambdab::Definition(), theLambdabInelasticProcess)); 541 G4HadronicProcess* theAntiLambdabInelasticProcess = 542 new G4HadronInelasticProcess("anti_lambda_bInelastic", G4AntiLambdab::Definition()); 543 fProcessMap.insert(ProcessPair(G4AntiLambdab::Definition(), theAntiLambdabInelasticProcess)); 544 G4HadronicProcess* theXibZeroInelasticProcess = 545 new G4HadronInelasticProcess("xi_b0Inelastic", G4XibZero::Definition()); 546 fProcessMap.insert(ProcessPair(G4XibZero::Definition(), theXibZeroInelasticProcess)); 547 G4HadronicProcess* theAntiXibZeroInelasticProcess = 548 new G4HadronInelasticProcess("anti_xi_b0Inelastic", G4AntiXibZero::Definition()); 549 fProcessMap.insert(ProcessPair(G4AntiXibZero::Definition(), theAntiXibZeroInelasticProcess)); 550 G4HadronicProcess* theXibMinusInelasticProcess = 551 new G4HadronInelasticProcess("xi_b-Inelastic", G4XibMinus::Definition()); 552 fProcessMap.insert(ProcessPair(G4XibMinus::Definition(), theXibMinusInelasticProcess)); 553 G4HadronicProcess* theAntiXibMinusInelasticProcess = 554 new G4HadronInelasticProcess("anti_xi_b-Inelastic", G4AntiXibMinus::Definition()); 555 fProcessMap.insert(ProcessPair(G4AntiXibMinus::Definition(), theAntiXibMinusInelasticProcess)); 556 G4HadronicProcess* theOmegabMinusInelasticProcess = 557 new G4HadronInelasticProcess("omega_b-Inelastic", G4OmegabMinus::Definition()); 558 fProcessMap.insert(ProcessPair(G4OmegabMinus::Definition(), theOmegabMinusInelasticProcess)); 559 G4HadronicProcess* theAntiOmegabMinusInelasticProcess = 560 new G4HadronInelasticProcess("anti_omega_b-Inelastic", G4AntiOmegabMinus::Definition()); 561 fProcessMap.insert( 562 ProcessPair(G4AntiOmegabMinus::Definition(), theAntiOmegabMinusInelasticProcess)); 563 564 // FLUKA ADAPTATION 565 if (fPhysicsCase == "CFLUKAHI") { 566 #ifdef G4_USE_FLUKA 567 fPhysicsCaseIsSupported = true; 568 569 // FLUKA HADRON - NUCLEUS INELASTIC XS 570 auto flukaInelasticScatteringXS = new FLUKAInelasticScatteringXS(); 571 572 theProtonInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 573 theAntiProtonInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 574 theNeutronInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 575 theAntiNeutronInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 576 thePionMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 577 thePionPlusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 578 theKaonMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 579 theKaonPlusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 580 theKaonZeroLInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 581 theKaonZeroSInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 582 583 // * Hyperons and anti-hyperons 584 // NB: Could have used G4HadParticles::GetHyperons() and G4HadParticles::GetAntiHyperons(). 585 // * BC hadrons 586 // NB: Could have used G4HadParticles::GetBCHadrons. 587 theLambdaInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 588 theSigmaMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 589 theSigmaPlusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 590 theXiMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 591 theXiZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 592 theOmegaMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 593 theAntiLambdaInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 594 theAntiSigmaMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 595 theAntiSigmaPlusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 596 theAntiXiMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 597 theAntiXiZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 598 theAntiOmegaMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 599 theDPlusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 600 theDMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 601 theDZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 602 theAntiDZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 603 theDsPlusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 604 theDsMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 605 theBPlusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 606 theBMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 607 theBZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 608 theAntiBZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 609 theBsZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 610 theAntiBsZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 611 theBcPlusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 612 theBcMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 613 theLambdacPlusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 614 theAntiLambdacPlusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 615 theXicPlusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 616 theAntiXicPlusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 617 theXicZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 618 theAntiXicZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 619 theOmegacZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 620 theAntiOmegacZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 621 theLambdabInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 622 theAntiLambdabInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 623 theXibZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 624 theAntiXibZeroInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 625 theXibMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 626 theAntiXibMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 627 theOmegabMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 628 theAntiOmegabMinusInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 629 630 // * Light ions 631 // NB: Could have used G4HadParticles::GetLightIons(). 632 theDeuteronInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 633 theTritonInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 634 theHe3InelasticProcess->AddDataSet(flukaInelasticScatteringXS); 635 theAlphaInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 636 637 // * Light anti-ions 638 // NB: Could have used G4HadParticles::GetLightAntiIons(). 639 theAntiDeuteronInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 640 theAntiTritonInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 641 theAntiHe3InelasticProcess->AddDataSet(flukaInelasticScatteringXS); 642 theAntiAlphaInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 643 644 // * Ions 645 // NB: Could have used G4GenericIon::GenericIon(). 646 theIonInelasticProcess->AddDataSet(flukaInelasticScatteringXS); 647 #endif 648 } 649 // G4 XS 650 else { 651 // Add the cross sections to the corresponding hadronic processes 652 thePionMinusInelasticProcess->AddDataSet(thePionMinusXSdata); 653 thePionPlusInelasticProcess->AddDataSet(thePionPlusXSdata); 654 theKaonMinusInelasticProcess->AddDataSet(theKaonXSdata); 655 theKaonPlusInelasticProcess->AddDataSet(theKaonXSdata); 656 theKaonZeroLInelasticProcess->AddDataSet(theKaonXSdata); 657 theKaonZeroSInelasticProcess->AddDataSet(theKaonXSdata); 658 theProtonInelasticProcess->AddDataSet(theProtonXSdata); 659 theNeutronInelasticProcess->AddDataSet(theNeutronXSdata); 660 theDeuteronInelasticProcess->AddDataSet(theNuclNuclXSdata); 661 theTritonInelasticProcess->AddDataSet(theNuclNuclXSdata); 662 theHe3InelasticProcess->AddDataSet(theNuclNuclXSdata); 663 theAlphaInelasticProcess->AddDataSet(theNuclNuclXSdata); 664 theIonInelasticProcess->AddDataSet(theNuclNuclXSdata); 665 theLambdaInelasticProcess->AddDataSet(theHyperonsXSdata); 666 theSigmaMinusInelasticProcess->AddDataSet(theHyperonsXSdata); 667 theSigmaPlusInelasticProcess->AddDataSet(theHyperonsXSdata); 668 theXiMinusInelasticProcess->AddDataSet(theHyperonsXSdata); 669 theXiZeroInelasticProcess->AddDataSet(theHyperonsXSdata); 670 theOmegaMinusInelasticProcess->AddDataSet(theHyperonsXSdata); 671 theAntiProtonInelasticProcess->AddDataSet(theAntibaryonsXSdata); 672 theAntiNeutronInelasticProcess->AddDataSet(theAntibaryonsXSdata); 673 theAntiDeuteronInelasticProcess->AddDataSet(theAntibaryonsXSdata); 674 theAntiTritonInelasticProcess->AddDataSet(theAntibaryonsXSdata); 675 theAntiHe3InelasticProcess->AddDataSet(theAntibaryonsXSdata); 676 theAntiAlphaInelasticProcess->AddDataSet(theHyperonsXSdata); 677 theAntiLambdaInelasticProcess->AddDataSet(theHyperonsXSdata); 678 theAntiSigmaMinusInelasticProcess->AddDataSet(theHyperonsXSdata); 679 theAntiSigmaPlusInelasticProcess->AddDataSet(theHyperonsXSdata); 680 theAntiXiMinusInelasticProcess->AddDataSet(theHyperonsXSdata); 681 theAntiXiZeroInelasticProcess->AddDataSet(theHyperonsXSdata); 682 theAntiOmegaMinusInelasticProcess->AddDataSet(theHyperonsXSdata); 683 684 theDPlusInelasticProcess->AddDataSet(theHyperonsXSdata); 685 theDMinusInelasticProcess->AddDataSet(theHyperonsXSdata); 686 theDZeroInelasticProcess->AddDataSet(theHyperonsXSdata); 687 theAntiDZeroInelasticProcess->AddDataSet(theHyperonsXSdata); 688 theDsPlusInelasticProcess->AddDataSet(theHyperonsXSdata); 689 theDsMinusInelasticProcess->AddDataSet(theHyperonsXSdata); 690 theBPlusInelasticProcess->AddDataSet(theHyperonsXSdata); 691 theBMinusInelasticProcess->AddDataSet(theHyperonsXSdata); 692 theBZeroInelasticProcess->AddDataSet(theHyperonsXSdata); 693 theAntiBZeroInelasticProcess->AddDataSet(theHyperonsXSdata); 694 theBsZeroInelasticProcess->AddDataSet(theHyperonsXSdata); 695 theAntiBsZeroInelasticProcess->AddDataSet(theHyperonsXSdata); 696 theBcPlusInelasticProcess->AddDataSet(theHyperonsXSdata); 697 theBcMinusInelasticProcess->AddDataSet(theHyperonsXSdata); 698 theLambdacPlusInelasticProcess->AddDataSet(theHyperonsXSdata); 699 theAntiLambdacPlusInelasticProcess->AddDataSet(theHyperonsXSdata); 700 theXicPlusInelasticProcess->AddDataSet(theHyperonsXSdata); 701 theAntiXicPlusInelasticProcess->AddDataSet(theHyperonsXSdata); 702 theXicZeroInelasticProcess->AddDataSet(theHyperonsXSdata); 703 theAntiXicZeroInelasticProcess->AddDataSet(theHyperonsXSdata); 704 theOmegacZeroInelasticProcess->AddDataSet(theHyperonsXSdata); 705 theAntiOmegacZeroInelasticProcess->AddDataSet(theHyperonsXSdata); 706 theLambdabInelasticProcess->AddDataSet(theHyperonsXSdata); 707 theAntiLambdabInelasticProcess->AddDataSet(theHyperonsXSdata); 708 theXibZeroInelasticProcess->AddDataSet(theHyperonsXSdata); 709 theAntiXibZeroInelasticProcess->AddDataSet(theHyperonsXSdata); 710 theXibMinusInelasticProcess->AddDataSet(theHyperonsXSdata); 711 theAntiXibMinusInelasticProcess->AddDataSet(theHyperonsXSdata); 712 theOmegabMinusInelasticProcess->AddDataSet(theHyperonsXSdata); 713 theAntiOmegabMinusInelasticProcess->AddDataSet(theHyperonsXSdata); 714 } 715 716 // FLUKA HADRON - NUCLEUS INELASTIC MODEL 717 if (fPhysicsCase == "CFLUKAHI") { 718 #ifdef G4_USE_FLUKA 719 720 // FLUKA HADRON - NUCLEUS INELASTIC MODEL 721 FLUKANuclearInelasticModel* flukaModel = new FLUKANuclearInelasticModel(); 722 723 theProtonInelasticProcess->RegisterMe(flukaModel); 724 theAntiProtonInelasticProcess->RegisterMe(flukaModel); 725 theNeutronInelasticProcess->RegisterMe(flukaModel); 726 theAntiNeutronInelasticProcess->RegisterMe(flukaModel); 727 thePionMinusInelasticProcess->RegisterMe(flukaModel); 728 thePionPlusInelasticProcess->RegisterMe(flukaModel); 729 theKaonMinusInelasticProcess->RegisterMe(flukaModel); 730 theKaonPlusInelasticProcess->RegisterMe(flukaModel); 731 theKaonZeroLInelasticProcess->RegisterMe(flukaModel); 732 theKaonZeroSInelasticProcess->RegisterMe(flukaModel); 733 734 // * Hyperons and anti-hyperons 735 // NB: Could have used G4HadParticles::GetHyperons() and G4HadParticles::GetAntiHyperons(). 736 // * BC hadrons 737 // NB: Could have used G4HadParticles::GetBCHadrons. 738 theLambdaInelasticProcess->RegisterMe(flukaModel); 739 theAntiLambdaInelasticProcess->RegisterMe(flukaModel); 740 theSigmaMinusInelasticProcess->RegisterMe(flukaModel); 741 theAntiSigmaMinusInelasticProcess->RegisterMe(flukaModel); 742 theSigmaPlusInelasticProcess->RegisterMe(flukaModel); 743 theAntiSigmaPlusInelasticProcess->RegisterMe(flukaModel); 744 theXiMinusInelasticProcess->RegisterMe(flukaModel); 745 theAntiXiMinusInelasticProcess->RegisterMe(flukaModel); 746 theXiZeroInelasticProcess->RegisterMe(flukaModel); 747 theAntiXiZeroInelasticProcess->RegisterMe(flukaModel); 748 theOmegaMinusInelasticProcess->RegisterMe(flukaModel); 749 theAntiOmegaMinusInelasticProcess->RegisterMe(flukaModel); 750 theDPlusInelasticProcess->RegisterMe(flukaModel); 751 theDMinusInelasticProcess->RegisterMe(flukaModel); 752 theDZeroInelasticProcess->RegisterMe(flukaModel); 753 theAntiDZeroInelasticProcess->RegisterMe(flukaModel); 754 theDsPlusInelasticProcess->RegisterMe(flukaModel); 755 theDsMinusInelasticProcess->RegisterMe(flukaModel); 756 theBPlusInelasticProcess->RegisterMe(flukaModel); 757 theBMinusInelasticProcess->RegisterMe(flukaModel); 758 theBZeroInelasticProcess->RegisterMe(flukaModel); 759 theAntiBZeroInelasticProcess->RegisterMe(flukaModel); 760 theBsZeroInelasticProcess->RegisterMe(flukaModel); 761 theAntiBsZeroInelasticProcess->RegisterMe(flukaModel); 762 theBcPlusInelasticProcess->RegisterMe(flukaModel); 763 theBcMinusInelasticProcess->RegisterMe(flukaModel); 764 theLambdacPlusInelasticProcess->RegisterMe(flukaModel); 765 theAntiLambdacPlusInelasticProcess->RegisterMe(flukaModel); 766 theXicPlusInelasticProcess->RegisterMe(flukaModel); 767 theAntiXicPlusInelasticProcess->RegisterMe(flukaModel); 768 theXicZeroInelasticProcess->RegisterMe(flukaModel); 769 theAntiXicZeroInelasticProcess->RegisterMe(flukaModel); 770 theOmegacZeroInelasticProcess->RegisterMe(flukaModel); 771 theAntiOmegacZeroInelasticProcess->RegisterMe(flukaModel); 772 theLambdabInelasticProcess->RegisterMe(flukaModel); 773 theAntiLambdabInelasticProcess->RegisterMe(flukaModel); 774 theXibZeroInelasticProcess->RegisterMe(flukaModel); 775 theAntiXibZeroInelasticProcess->RegisterMe(flukaModel); 776 theXibMinusInelasticProcess->RegisterMe(flukaModel); 777 theAntiXibMinusInelasticProcess->RegisterMe(flukaModel); 778 theOmegabMinusInelasticProcess->RegisterMe(flukaModel); 779 theAntiOmegabMinusInelasticProcess->RegisterMe(flukaModel); 780 781 // * Light ions 782 // NB: Could have used G4HadParticles::GetLightIons(). 783 theDeuteronInelasticProcess->RegisterMe(flukaModel); 784 theTritonInelasticProcess->RegisterMe(flukaModel); 785 theHe3InelasticProcess->RegisterMe(flukaModel); 786 theAlphaInelasticProcess->RegisterMe(flukaModel); 787 788 // * Light anti-ions 789 // NB: Could have used G4HadParticles::GetLightAntiIons(). 790 theAntiDeuteronInelasticProcess->RegisterMe(flukaModel); 791 theAntiTritonInelasticProcess->RegisterMe(flukaModel); 792 theAntiHe3InelasticProcess->RegisterMe(flukaModel); 793 theAntiAlphaInelasticProcess->RegisterMe(flukaModel); 794 795 // * Ions 796 // NB: Could have used G4GenericIon::GenericIon(). 797 theIonInelasticProcess->RegisterMe(flukaModel); 798 #endif 799 } 800 801 // G4 MODELS 802 // Register the proper hadronic model(s) to the corresponding hadronic processes. 803 // Note: hadronic models ("BERT", "BIC", "IonBIC", "INCL", "FTFP", "QGSP") are 804 // used for the hadrons and energies they are applicable 805 // (exception for INCL, which in recent versions of Geant4 can handle 806 // more hadron types and higher energies than considered here). 807 // For "physics-list proxies" ("FTFP_BERT", "FTFP_BERT_ATL", "QGSP_BERT", 808 // "QGSP_BIC", "FTFP_INCLXX"), all hadron types and all energies are covered 809 // by combining different hadronic models - similarly (but not identically) 810 // to the corresponding physics lists. 811 if (fPhysicsCase == "BIC" || fPhysicsCase == "QGSP_BIC") { 812 // The BIC model is applicable to nucleons and pions, 813 // whereas in the physics list QGSP_BIC it is used only for nucleons 814 fPhysicsCaseIsSupported = true; 815 theProtonInelasticProcess->RegisterMe(theBICmodel); 816 theNeutronInelasticProcess->RegisterMe(theBICmodel); 817 if (fPhysicsCase == "BIC") { 818 thePionMinusInelasticProcess->RegisterMe(theBICmodel); 819 thePionPlusInelasticProcess->RegisterMe(theBICmodel); 820 } 821 else { 822 thePionMinusInelasticProcess->RegisterMe(theBERTmodel); 823 thePionPlusInelasticProcess->RegisterMe(theBERTmodel); 824 } 825 } 826 else if (fPhysicsCase == "INCL" || fPhysicsCase == "FTFP_INCLXX") { 827 // We consider here for simplicity only nucleons and pions 828 // (although recent versions of INCL can handle others particles as well) 829 fPhysicsCaseIsSupported = true; 830 thePionMinusInelasticProcess->RegisterMe(theINCLmodel); 831 thePionPlusInelasticProcess->RegisterMe(theINCLmodel); 832 theProtonInelasticProcess->RegisterMe(theINCLmodel); 833 theNeutronInelasticProcess->RegisterMe(theINCLmodel); 834 } 835 if (fPhysicsCase == "IonBIC" || fPhysicsCase == "FTFP_BERT_ATL" || fPhysicsCase == "FTFP_BERT" 836 || fPhysicsCase == "FTFP_INCLXX" || fPhysicsCase == "QGSP_BERT" || fPhysicsCase == "QGSP_BIC") 837 { 838 // The Binary Light Ion model is used for ions in all physics lists 839 fPhysicsCaseIsSupported = true; 840 theDeuteronInelasticProcess->RegisterMe(theIonBICmodel); 841 theTritonInelasticProcess->RegisterMe(theIonBICmodel); 842 theHe3InelasticProcess->RegisterMe(theIonBICmodel); 843 theAlphaInelasticProcess->RegisterMe(theIonBICmodel); 844 theIonInelasticProcess->RegisterMe(theIonBICmodel); 845 } 846 if (fPhysicsCase == "QGSP" || fPhysicsCase == "QGSP_BERT" || fPhysicsCase == "QGSP_BIC") { 847 fPhysicsCaseIsSupported = true; 848 thePionMinusInelasticProcess->RegisterMe(theQGSPmodel); 849 thePionPlusInelasticProcess->RegisterMe(theQGSPmodel); 850 theKaonMinusInelasticProcess->RegisterMe(theQGSPmodel); 851 theKaonPlusInelasticProcess->RegisterMe(theQGSPmodel); 852 theKaonZeroLInelasticProcess->RegisterMe(theQGSPmodel); 853 theKaonZeroSInelasticProcess->RegisterMe(theQGSPmodel); 854 theProtonInelasticProcess->RegisterMe(theQGSPmodel); 855 theNeutronInelasticProcess->RegisterMe(theQGSPmodel); 856 theLambdaInelasticProcess->RegisterMe(theQGSPmodel); 857 theSigmaMinusInelasticProcess->RegisterMe(theQGSPmodel); 858 theSigmaPlusInelasticProcess->RegisterMe(theQGSPmodel); 859 theXiMinusInelasticProcess->RegisterMe(theQGSPmodel); 860 theXiZeroInelasticProcess->RegisterMe(theQGSPmodel); 861 theOmegaMinusInelasticProcess->RegisterMe(theQGSPmodel); 862 theAntiProtonInelasticProcess->RegisterMe(theQGSPmodel); 863 theAntiNeutronInelasticProcess->RegisterMe(theQGSPmodel); 864 theAntiLambdaInelasticProcess->RegisterMe(theQGSPmodel); 865 theAntiSigmaMinusInelasticProcess->RegisterMe(theQGSPmodel); 866 theAntiSigmaPlusInelasticProcess->RegisterMe(theQGSPmodel); 867 theAntiXiMinusInelasticProcess->RegisterMe(theQGSPmodel); 868 theAntiXiZeroInelasticProcess->RegisterMe(theQGSPmodel); 869 theAntiOmegaMinusInelasticProcess->RegisterMe(theQGSPmodel); 870 theDPlusInelasticProcess->RegisterMe(theQGSPmodel); 871 theDMinusInelasticProcess->RegisterMe(theQGSPmodel); 872 theDZeroInelasticProcess->RegisterMe(theQGSPmodel); 873 theAntiDZeroInelasticProcess->RegisterMe(theQGSPmodel); 874 theDsPlusInelasticProcess->RegisterMe(theQGSPmodel); 875 theDsMinusInelasticProcess->RegisterMe(theQGSPmodel); 876 theBPlusInelasticProcess->RegisterMe(theQGSPmodel); 877 theBMinusInelasticProcess->RegisterMe(theQGSPmodel); 878 theBZeroInelasticProcess->RegisterMe(theQGSPmodel); 879 theAntiBZeroInelasticProcess->RegisterMe(theQGSPmodel); 880 theBsZeroInelasticProcess->RegisterMe(theQGSPmodel); 881 theAntiBsZeroInelasticProcess->RegisterMe(theQGSPmodel); 882 theBcPlusInelasticProcess->RegisterMe(theQGSPmodel); 883 theBcMinusInelasticProcess->RegisterMe(theQGSPmodel); 884 theLambdacPlusInelasticProcess->RegisterMe(theQGSPmodel); 885 theAntiLambdacPlusInelasticProcess->RegisterMe(theQGSPmodel); 886 theXicPlusInelasticProcess->RegisterMe(theQGSPmodel); 887 theAntiXicPlusInelasticProcess->RegisterMe(theQGSPmodel); 888 theXicZeroInelasticProcess->RegisterMe(theQGSPmodel); 889 theAntiXicZeroInelasticProcess->RegisterMe(theQGSPmodel); 890 theOmegacZeroInelasticProcess->RegisterMe(theQGSPmodel); 891 theAntiOmegacZeroInelasticProcess->RegisterMe(theQGSPmodel); 892 theLambdabInelasticProcess->RegisterMe(theQGSPmodel); 893 theAntiLambdabInelasticProcess->RegisterMe(theQGSPmodel); 894 theXibZeroInelasticProcess->RegisterMe(theQGSPmodel); 895 theAntiXibZeroInelasticProcess->RegisterMe(theQGSPmodel); 896 theXibMinusInelasticProcess->RegisterMe(theQGSPmodel); 897 theAntiXibMinusInelasticProcess->RegisterMe(theQGSPmodel); 898 theOmegabMinusInelasticProcess->RegisterMe(theQGSPmodel); 899 theAntiOmegabMinusInelasticProcess->RegisterMe(theQGSPmodel); 900 } 901 if (fPhysicsCase == "BERT" || fPhysicsCase == "FTFP_BERT_ATL" || fPhysicsCase == "FTFP_BERT" 902 || fPhysicsCase == "QGSP_BERT") 903 { 904 // The BERT model is used for pions and nucleons in all Bertini-based physics lists 905 fPhysicsCaseIsSupported = true; 906 thePionMinusInelasticProcess->RegisterMe(theBERTmodel); 907 thePionPlusInelasticProcess->RegisterMe(theBERTmodel); 908 theProtonInelasticProcess->RegisterMe(theBERTmodel); 909 theNeutronInelasticProcess->RegisterMe(theBERTmodel); 910 } 911 if (fPhysicsCase == "BERT" || fPhysicsCase == "FTFP_BERT_ATL" || fPhysicsCase == "FTFP_BERT" 912 || fPhysicsCase == "FTFP_INCLXX" || fPhysicsCase == "QGSP_BERT" || fPhysicsCase == "QGSP_BIC") 913 { 914 // The BERT model is used for kaons and hyperons in all physics lists 915 fPhysicsCaseIsSupported = true; 916 theKaonMinusInelasticProcess->RegisterMe(theBERTmodel); 917 theKaonPlusInelasticProcess->RegisterMe(theBERTmodel); 918 theKaonZeroLInelasticProcess->RegisterMe(theBERTmodel); 919 theKaonZeroSInelasticProcess->RegisterMe(theBERTmodel); 920 theLambdaInelasticProcess->RegisterMe(theBERTmodel); 921 theSigmaMinusInelasticProcess->RegisterMe(theBERTmodel); 922 theSigmaPlusInelasticProcess->RegisterMe(theBERTmodel); 923 theXiMinusInelasticProcess->RegisterMe(theBERTmodel); 924 theXiZeroInelasticProcess->RegisterMe(theBERTmodel); 925 theOmegaMinusInelasticProcess->RegisterMe(theBERTmodel); 926 } 927 if (fPhysicsCase == "FTFP" || fPhysicsCase == "FTFP_BERT_ATL" || fPhysicsCase == "FTFP_BERT" 928 || fPhysicsCase == "FTFP_INCLXX" || fPhysicsCase == "QGSP_BERT" || fPhysicsCase == "QGSP_BIC") 929 { 930 // The FTFP model is applied for all hadrons, but in different energy intervals according 931 // whether it is consider as a stand-alone hadronic model, or within physics lists 932 fPhysicsCaseIsSupported = true; 933 theAntiDeuteronInelasticProcess->RegisterMe(theFTFPmodel); 934 theAntiTritonInelasticProcess->RegisterMe(theFTFPmodel); 935 theAntiHe3InelasticProcess->RegisterMe(theFTFPmodel); 936 theAntiAlphaInelasticProcess->RegisterMe(theFTFPmodel); 937 G4TheoFSGenerator* theFTFPmodelToBeUsed = theFTFPmodel_aboveThreshold; 938 if (fPhysicsCase == "FTFP") { 939 theFTFPmodelToBeUsed = theFTFPmodel; 940 } 941 else if (fPhysicsCase == "QGSP_BERT" || fPhysicsCase == "QGSP_BIC") { 942 theFTFPmodelToBeUsed = theFTFPmodel_constrained; 943 } 944 thePionMinusInelasticProcess->RegisterMe(theFTFPmodelToBeUsed); 945 thePionPlusInelasticProcess->RegisterMe(theFTFPmodelToBeUsed); 946 theKaonMinusInelasticProcess->RegisterMe(theFTFPmodelToBeUsed); 947 theKaonPlusInelasticProcess->RegisterMe(theFTFPmodelToBeUsed); 948 theKaonZeroLInelasticProcess->RegisterMe(theFTFPmodelToBeUsed); 949 theKaonZeroSInelasticProcess->RegisterMe(theFTFPmodelToBeUsed); 950 theProtonInelasticProcess->RegisterMe(theFTFPmodelToBeUsed); 951 theAntiProtonInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 952 theNeutronInelasticProcess->RegisterMe(theFTFPmodelToBeUsed); 953 theAntiNeutronInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 954 theLambdaInelasticProcess->RegisterMe(theFTFPmodelToBeUsed); 955 theAntiLambdaInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 956 theSigmaMinusInelasticProcess->RegisterMe(theFTFPmodelToBeUsed); 957 theAntiSigmaMinusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 958 theSigmaPlusInelasticProcess->RegisterMe(theFTFPmodelToBeUsed); 959 theAntiSigmaPlusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 960 theXiMinusInelasticProcess->RegisterMe(theFTFPmodelToBeUsed); 961 theAntiXiMinusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 962 theXiZeroInelasticProcess->RegisterMe(theFTFPmodelToBeUsed); 963 theAntiXiZeroInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 964 theOmegaMinusInelasticProcess->RegisterMe(theFTFPmodelToBeUsed); 965 theAntiOmegaMinusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 966 theDPlusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 967 theDMinusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 968 theDZeroInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 969 theAntiDZeroInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 970 theDsPlusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 971 theDsMinusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 972 theBPlusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 973 theBMinusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 974 theBZeroInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 975 theAntiBZeroInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 976 theBsZeroInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 977 theAntiBsZeroInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 978 theBcPlusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 979 theBcMinusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 980 theLambdacPlusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 981 theAntiLambdacPlusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 982 theXicPlusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 983 theAntiXicPlusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 984 theXicZeroInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 985 theAntiXicZeroInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 986 theOmegacZeroInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 987 theAntiOmegacZeroInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 988 theLambdabInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 989 theAntiLambdabInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 990 theXibZeroInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 991 theAntiXibZeroInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 992 theXibMinusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 993 theAntiXibMinusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 994 theOmegabMinusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 995 theAntiOmegabMinusInelasticProcess->RegisterMe(theFTFPmodel_belowThreshold); 996 theFTFPmodelToBeUsed = theFTFPmodel_aboveThreshold; 997 if (fPhysicsCase == "FTFP") theFTFPmodelToBeUsed = theFTFPmodel; 998 theDeuteronInelasticProcess->RegisterMe(theFTFPmodelToBeUsed); 999 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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 1016 1017 HadronicGenerator::~HadronicGenerator() 1018 { 1019 fPartTable->DeleteAllParticles(); 1020 } 1021 1022 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 1032 1033 G4bool HadronicGenerator::IsApplicable(G4ParticleDefinition* projectileDefinition, 1034 const G4double projectileEnergy) const 1035 { 1036 if (projectileDefinition == nullptr) return false; 1037 G4bool isApplicable = true; 1038 // No restrictions for "physics list proxies" because they cover all hadron types and energies. 1039 // For the individual models, instead, we need to consider their limitations. 1040 if (fPhysicsCase == "BERT") { 1041 // We consider BERT model below 15 GeV 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 // We consider QGSP above 2 GeV and not for ions or anti-ions 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 // We consider BIC and INCL models only for pions and nucleons below 10 GeV 1074 // (although in recent versions INCL is capable of handling more hadrons 1075 // and up to higher energies) 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 // We consider IonBIC models only for deuteron, triton, He3, alpha 1087 // with energies below 10 GeV / nucleon 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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 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 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 1116 1117 G4VParticleChange* HadronicGenerator::GenerateInteraction( 1118 G4ParticleDefinition* projectileDefinition, const G4double projectileEnergy, 1119 const G4ThreeVector& projectileDirection, G4Material* targetMaterial) 1120 { 1121 // This is the most important method of the HadronicGenerator class: 1122 // the method performs the specified hadronic interaction 1123 // (by invoking the "PostStepDoIt" method of the corresponding hadronic process) 1124 // and returns the final state, i.e. the secondaries produced by the collision. 1125 // It is a relatively short method because the heavy load of setting up all 1126 // possible hadronic processes - with their hadronic models, transition regions, 1127 // and cross sections (the latter is needed for sampling the target nucleus from 1128 // the target material) - was already done by the constructor of the class. 1129 G4VParticleChange* aChange = nullptr; 1130 1131 if (projectileDefinition == nullptr) { 1132 G4cerr << "ERROR: projectileDefinition is NULL !" << G4endl; 1133 return aChange; 1134 } 1135 1136 // Debugging print-out 1137 // G4cout << "\t" << projectileDefinition->GetParticleName() 1138 // << "\t" << projectileEnergy/CLHEP::GeV 1139 // << " GeV \t" << projectileDirection 1140 // << "\t" << ( targetMaterial ? targetMaterial->GetName() : "NULL" ); 1141 if (!IsApplicable(projectileDefinition, projectileEnergy)) { 1142 // G4cout << " -> NOT applicable !" ; //<< G4endl; // Debugging print-out 1143 return aChange; 1144 } 1145 // G4cout << G4endl; 1146 1147 // Check Geant4 state (not strictly needed) 1148 // if ( ! G4StateManager::GetStateManager()->SetNewState( G4State_PreInit ) ) { 1149 // G4cerr << "ERROR: No possible to set G4State_PreInit !" << G4endl; 1150 // return aChange; 1151 //} 1152 1153 // Geometry definition (not strictly needed) 1154 // const G4double dimX = 1.0*mm; 1155 // const G4double dimY = 1.0*mm; 1156 // const G4double dimZ = 1.0*mm; 1157 // G4Box* sFrame = new G4Box( "Box", dimX, dimY, dimZ ); 1158 // G4LogicalVolume* lFrame = new G4LogicalVolume( sFrame, targetMaterial, "Box", 0, 0, 0 ); 1159 // G4PVPlacement* pFrame = new G4PVPlacement( 0, G4ThreeVector(), "Box", lFrame, 0, false, 0 ); 1160 // G4TransportationManager::GetTransportationManager()->SetWorldForTracking( pFrame ); 1161 1162 // Projectile track & step 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); // Not strictly needed 1168 gTrack->SetTouchableHandle(fpTouchable); // Not strictly needed 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 // Change Geant4 state: from "PreInit" to "Idle" (not strictly needed) 1181 // if ( ! G4StateManager::GetStateManager()->SetNewState( G4State_Idle ) ) { 1182 // G4cerr << "ERROR: No possible to set G4State_Idle !" << G4endl; 1183 // return aChange; 1184 //} 1185 1186 // Finally, the hadronic interaction: hadron projectile and ion projectile 1187 // need to be treated slightly differently 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 // delete pFrame; 1207 // delete lFrame; 1208 // delete sFrame; 1209 1210 return aChange; 1211 } 1212 1213 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 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 // FTFP has handled the inelastic hadronic interaction. 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 // ftfModel points to the G4FTFModel object instance that handled the 1230 // inelastic hadronic interaction. 1231 impactParameter = ftfModel->GetImpactParameter(); 1232 // G4cout << "\t impactParameter = " << impactParameter/fermi << " fm" << G4endl; 1233 } 1234 } 1235 } 1236 return impactParameter; 1237 } 1238 1239 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 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 // G4cout << "\t numProjectileSpectatorNucleons = " 1256 // << numProjectileSpectatorNucleons << G4endl; 1257 } 1258 } 1259 } 1260 return numProjectileSpectatorNucleons; 1261 } 1262 1263 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 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 // G4cout << "\t numTargetSpectatorNucleons = " << numTargetSpectatorNucleons << G4endl; 1280 } 1281 } 1282 } 1283 return numTargetSpectatorNucleons; 1284 } 1285 1286 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 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 // G4cout << "\t numNNcollisions = " << numNNcollisions << G4endl; 1303 } 1304 } 1305 } 1306 return numNNcollisions; 1307 } 1308 1309 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 1310