Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4HadronInelasticQBBC.cc
Go to the documentation of this file.
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
40
41#include "G4SystemOfUnits.hh"
42
43#include "G4HadronicProcess.hh"
47
49
50#include "G4TheoFSGenerator.hh"
51#include "G4FTFModel.hh"
54
57
60#include "G4NeutronCaptureXS.hh"
62
64
65#include "G4CascadeInterface.hh"
66#include "G4BinaryCascade.hh"
68
69#include "G4PreCompoundModel.hh"
71
73#include "G4HadronicBuilder.hh"
74#include "G4HadParticles.hh"
75#include "G4HadProcesses.hh"
76#include "G4BuilderType.hh"
77
78// factory
80//
82
84 : G4VHadronPhysics("hInelasticQBBC")
85{
87 auto param = G4HadronicParameters::Instance();
88 param->SetEnableBCParticles(true);
89 param->SetEnableNeutronGeneralProcess(false);
90 param->SetUseRFilesForXS(true);
91 param->SetVerboseLevel(ver);
92}
93
97
99{
101 G4bool useFactorXS = param->ApplyFactorXS();
103
104 // configure models
105 const G4double eminFtf = param->GetMinEnergyTransitionFTF_Cascade();
106 const G4double eminBert = 1.0*CLHEP::GeV;
107 const G4double emaxBic = 1.5*CLHEP::GeV;
108 const G4double emaxBert = param->GetMaxEnergyTransitionFTF_Cascade();
109 const G4double emaxBertPions = 12.*CLHEP::GeV;
110 const G4double emax = param->GetMaxEnergy();
111
112 if(G4Threading::IsMasterThread() && param->GetVerboseLevel() > 0) {
113 G4cout << "### HadronInelasticQBBC Construct Process:\n"
114 << " Emin(FTFP)= " << eminFtf/CLHEP::GeV
115 << " GeV; Emax(FTFP)= " << emax/CLHEP::GeV << " GeV\n"
116 << " Emin(BERT)= " << eminBert/CLHEP::GeV
117 << " GeV; Emax(BERT)= " << emaxBert/CLHEP::GeV
118 << " GeV; Emax(BERTpions)= " << emaxBertPions/CLHEP::GeV
119 << " GeV;\n" << " Emin(BIC) = 0 GeV; Emax(BIC)= "
120 << emaxBic/CLHEP::GeV << " GeV." << G4endl;
121 }
122
123 // PreCompound and Evaporation models are instantiated here
124 G4PreCompoundModel* thePreCompound = nullptr;
127 thePreCompound = static_cast<G4PreCompoundModel*>(p);
128 if(!thePreCompound) { thePreCompound = new G4PreCompoundModel(); }
129
130 auto theFTFP = new G4TheoFSGenerator("FTFP");
131 auto theStringModel = new G4FTFModel();
132 theStringModel->SetFragmentationModel(new G4ExcitedStringDecay());
133 theFTFP->SetHighEnergyGenerator( theStringModel );
134 theFTFP->SetTransport( new G4GeneratorPrecompoundInterface() );
135 theFTFP->SetMinEnergy( eminFtf );
136 theFTFP->SetMaxEnergy( emax );
137
138 auto theBERT = new G4CascadeInterface();
139 theBERT->SetMinEnergy( eminBert );
140 theBERT->SetMaxEnergy( emaxBert );
141 theBERT->usePreCompoundDeexcitation();
142
143 auto theBERT1 = new G4CascadeInterface();
144 theBERT1->SetMinEnergy( eminBert );
145 theBERT1->SetMaxEnergy( emaxBertPions );
146 theBERT1->usePreCompoundDeexcitation();
147
148 auto theBIC = new G4BinaryCascade(thePreCompound);
149 theBIC->SetMaxEnergy( emaxBic );
150
151 // p
153 G4HadronicProcess* hp =
154 new G4HadronInelasticProcess( particle->GetParticleName()+"Inelastic", particle );
155 hp->AddDataSet( new G4BGGNucleonInelasticXS(particle) );
156 hp->RegisterMe(theFTFP);
157 hp->RegisterMe(theBERT);
158 hp->RegisterMe(theBIC);
159 ph->RegisterProcess(hp, particle);
160 if( useFactorXS ) hp->MultiplyCrossSectionBy( param->XSFactorNucleonInelastic() );
161
162 // n
163 particle = G4Neutron::Neutron();
164 G4HadronicProcess* ni =
165 new G4HadronInelasticProcess( "neutronInelastic", particle );
166 ni->RegisterMe(theFTFP);
167 ni->RegisterMe(theBERT);
168 ni->RegisterMe(theBIC);
170
171 // pi+
172 particle = G4PionPlus::PionPlus();
173 hp = new G4HadronInelasticProcess( particle->GetParticleName()+"Inelastic", particle );
174 hp->AddDataSet(new G4BGGPionInelasticXS(particle));
175 hp->RegisterMe(theFTFP);
176 hp->RegisterMe(theBERT1);
177 hp->RegisterMe(theBIC);
178 ph->RegisterProcess(hp, particle);
179 if( useFactorXS ) hp->MultiplyCrossSectionBy( param->XSFactorPionInelastic() );
180
181 // pi-
182 particle = G4PionMinus::PionMinus();
183 hp = new G4HadronInelasticProcess( particle->GetParticleName()+"Inelastic", particle );
184 hp->AddDataSet(new G4BGGPionInelasticXS(particle));
185 hp->RegisterMe(theFTFP);
186 hp->RegisterMe(theBERT1);
187 hp->RegisterMe(theBIC);
188 ph->RegisterProcess(hp, particle);
189 if( useFactorXS ) hp->MultiplyCrossSectionBy( param->XSFactorPionInelastic() );
190
191 // kaons
193
194 // high energy particles
195 if( emax > param->EnergyThresholdForHeavyHadrons() ) {
196
197 // pbar, nbar, anti light ions
199
200 // hyperons
202
203 // b-, c- baryons and mesons
204 if( param->EnableBCParticles() ) {
206 }
207 }
208}
@ bHadronInelastic
#define G4_DECLARE_PHYSCONSTR_FACTORY(physics_constructor)
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static void BuildNeutronInelasticAndCapture(G4HadronicProcess *)
static void BuildBCHadronsFTFP_BERT()
static void BuildHyperonsFTFP_BERT()
static void BuildKaonsFTFP_BERT()
static void BuildAntiLightIonsFTFP()
G4HadronicInteraction * FindModel(const G4String &name)
static G4HadronicInteractionRegistry * Instance()
static G4HadronicParameters * Instance()
G4double GetMinEnergyTransitionFTF_Cascade() const
G4double GetMaxEnergyTransitionFTF_Cascade() const
G4double EnergyThresholdForHeavyHadrons() const
G4double XSFactorPionInelastic() const
G4double XSFactorNucleonInelastic() const
void AddDataSet(G4VCrossSectionDataSet *aDataSet)
void MultiplyCrossSectionBy(G4double factor)
void RegisterMe(G4HadronicInteraction *a)
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
const G4String & GetParticleName() const
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
static G4PionMinus * PionMinus()
static G4PionPlus * PionPlus()
Definition G4PionPlus.cc:93
static G4Proton * Proton()
Definition G4Proton.cc:90
G4VHadronPhysics(const G4String &name="hInelastic", G4int verbose=0)
G4bool IsMasterThread()