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 // 27 //--------------------------------------------------------------------------- 28 // 29 // ClassName: G4HadronInelasticQBBC 30 // 31 // Author: 2 October 2009 V. Ivanchenko 32 // 33 // Modified: 34 // 12.10.2023 V.Ivanchenko added usage of the neutron general process 35 // 36 //---------------------------------------------------------------------------- 37 // 38 39 #include "G4HadronInelasticQBBC.hh" 40 41 #include "G4SystemOfUnits.hh" 42 43 #include "G4HadronicProcess.hh" 44 #include "G4HadronInelasticProcess.hh" 45 #include "G4NeutronCaptureProcess.hh" 46 #include "G4HadronicInteraction.hh" 47 48 #include "G4ParticleDefinition.hh" 49 50 #include "G4TheoFSGenerator.hh" 51 #include "G4FTFModel.hh" 52 #include "G4ExcitedStringDecay.hh" 53 #include "G4GeneratorPrecompoundInterface.hh" 54 55 #include "G4BGGNucleonInelasticXS.hh" 56 #include "G4BGGPionInelasticXS.hh" 57 58 #include "G4ParticleInelasticXS.hh" 59 #include "G4NeutronInelasticXS.hh" 60 #include "G4NeutronCaptureXS.hh" 61 #include "G4NeutronGeneralProcess.hh" 62 63 #include "G4CrossSectionInelastic.hh" 64 65 #include "G4CascadeInterface.hh" 66 #include "G4BinaryCascade.hh" 67 #include "G4NeutronRadCapture.hh" 68 69 #include "G4PreCompoundModel.hh" 70 #include "G4HadronicInteractionRegistry.hh" 71 72 #include "G4HadronicParameters.hh" 73 #include "G4HadronicBuilder.hh" 74 #include "G4HadParticles.hh" 75 #include "G4HadProcesses.hh" 76 #include "G4BuilderType.hh" 77 78 // factory 79 #include "G4PhysicsConstructorFactory.hh" 80 // 81 G4_DECLARE_PHYSCONSTR_FACTORY(G4HadronInelasticQBBC); 82 83 G4HadronInelasticQBBC::G4HadronInelasticQBBC(G4int ver) 84 : G4VHadronPhysics("hInelasticQBBC") 85 { 86 SetPhysicsType(bHadronInelastic); 87 auto param = G4HadronicParameters::Instance(); 88 param->SetEnableBCParticles(true); 89 param->SetEnableNeutronGeneralProcess(true); 90 param->SetVerboseLevel(ver); 91 } 92 93 G4HadronInelasticQBBC::G4HadronInelasticQBBC(const G4String&, G4int ver, 94 G4bool, G4bool,G4bool, G4bool, G4bool) : G4HadronInelasticQBBC(ver) 95 {} 96 97 void G4HadronInelasticQBBC::ConstructProcess() 98 { 99 G4HadronicParameters* param = G4HadronicParameters::Instance(); 100 G4bool useFactorXS = param->ApplyFactorXS(); 101 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper(); 102 103 // configure models 104 const G4double eminFtf = param->GetMinEnergyTransitionFTF_Cascade(); 105 const G4double eminBert = 1.0*CLHEP::GeV; 106 const G4double emaxBic = 1.5*CLHEP::GeV; 107 const G4double emaxBert = param->GetMaxEnergyTransitionFTF_Cascade(); 108 const G4double emaxBertPions = 12.*CLHEP::GeV; 109 const G4double emax = param->GetMaxEnergy(); 110 111 if(G4Threading::IsMasterThread() && param->GetVerboseLevel() > 0) { 112 G4cout << "### HadronInelasticQBBC Construct Process:\n" 113 << " Emin(FTFP)= " << eminFtf/CLHEP::GeV 114 << " GeV; Emax(FTFP)= " << emax/CLHEP::GeV << " GeV\n" 115 << " Emin(BERT)= " << eminBert/CLHEP::GeV 116 << " GeV; Emax(BERT)= " << emaxBert/CLHEP::GeV 117 << " GeV; Emax(BERTpions)= " << emaxBertPions/CLHEP::GeV 118 << " GeV;\n" << " Emin(BIC) = 0 GeV; Emax(BIC)= " 119 << emaxBic/CLHEP::GeV << " GeV." << G4endl; 120 } 121 122 // PreCompound and Evaporation models are instantiated here 123 G4PreCompoundModel* thePreCompound = nullptr; 124 G4HadronicInteraction* p = 125 G4HadronicInteractionRegistry::Instance()->FindModel("PRECO"); 126 thePreCompound = static_cast<G4PreCompoundModel*>(p); 127 if(!thePreCompound) { thePreCompound = new G4PreCompoundModel(); } 128 129 auto theFTFP = new G4TheoFSGenerator("FTFP"); 130 auto theStringModel = new G4FTFModel(); 131 theStringModel->SetFragmentationModel(new G4ExcitedStringDecay()); 132 theFTFP->SetHighEnergyGenerator( theStringModel ); 133 theFTFP->SetTransport( new G4GeneratorPrecompoundInterface() ); 134 theFTFP->SetMinEnergy( eminFtf ); 135 theFTFP->SetMaxEnergy( emax ); 136 137 auto theBERT = new G4CascadeInterface(); 138 theBERT->SetMinEnergy( eminBert ); 139 theBERT->SetMaxEnergy( emaxBert ); 140 theBERT->usePreCompoundDeexcitation(); 141 142 auto theBERT1 = new G4CascadeInterface(); 143 theBERT1->SetMinEnergy( eminBert ); 144 theBERT1->SetMaxEnergy( emaxBertPions ); 145 theBERT1->usePreCompoundDeexcitation(); 146 147 auto theBIC = new G4BinaryCascade(thePreCompound); 148 theBIC->SetMaxEnergy( emaxBic ); 149 150 // p 151 G4ParticleDefinition* particle = G4Proton::Proton(); 152 G4HadronicProcess* hp = 153 new G4HadronInelasticProcess( particle->GetParticleName()+"Inelastic", particle ); 154 hp->AddDataSet(new G4ParticleInelasticXS(particle)); 155 hp->RegisterMe(theFTFP); 156 hp->RegisterMe(theBERT); 157 hp->RegisterMe(theBIC); 158 ph->RegisterProcess(hp, particle); 159 if( useFactorXS ) hp->MultiplyCrossSectionBy( param->XSFactorNucleonInelastic() ); 160 161 // n 162 particle = G4Neutron::Neutron(); 163 G4HadronicProcess* ni = 164 new G4HadronInelasticProcess( "neutronInelastic", particle ); 165 ni->RegisterMe(theFTFP); 166 ni->RegisterMe(theBERT); 167 ni->RegisterMe(theBIC); 168 G4HadProcesses::BuildNeutronInelasticAndCapture(ni); 169 170 // pi+ 171 particle = G4PionPlus::PionPlus(); 172 hp = new G4HadronInelasticProcess( particle->GetParticleName()+"Inelastic", particle ); 173 hp->AddDataSet(new G4BGGPionInelasticXS(particle)); 174 hp->RegisterMe(theFTFP); 175 hp->RegisterMe(theBERT1); 176 hp->RegisterMe(theBIC); 177 ph->RegisterProcess(hp, particle); 178 if( useFactorXS ) hp->MultiplyCrossSectionBy( param->XSFactorPionInelastic() ); 179 180 // pi- 181 particle = G4PionMinus::PionMinus(); 182 hp = new G4HadronInelasticProcess( particle->GetParticleName()+"Inelastic", particle ); 183 hp->AddDataSet(new G4BGGPionInelasticXS(particle)); 184 hp->RegisterMe(theFTFP); 185 hp->RegisterMe(theBERT1); 186 hp->RegisterMe(theBIC); 187 ph->RegisterProcess(hp, particle); 188 if( useFactorXS ) hp->MultiplyCrossSectionBy( param->XSFactorPionInelastic() ); 189 190 // kaons 191 G4HadronicBuilder::BuildKaonsFTFP_BERT(); 192 193 // high energy particles 194 if( emax > param->EnergyThresholdForHeavyHadrons() ) { 195 196 // pbar, nbar, anti light ions 197 G4HadronicBuilder::BuildAntiLightIonsFTFP(); 198 199 // hyperons 200 G4HadronicBuilder::BuildHyperonsFTFP_BERT(); 201 202 // b-, c- baryons and mesons 203 if( param->EnableBCParticles() ) { 204 G4HadronicBuilder::BuildBCHadronsFTFP_BERT(); 205 } 206 } 207 } 208