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