Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4HadronHElasticPhysics.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// ClassName: G4HadronHElasticPhysics
29//
30// Author: 23 November 2006 V. Ivanchenko
31//
32// Modified:
33// 21.03.07 (V.Ivanchenko) Use G4BGGNucleonElasticXS and G4BGGPionElasticXS;
34// Reduce thresholds for HE and Q-models to zero
35// 03.06.2010 V.Ivanchenko cleanup constructors and ConstructProcess method
36// 06.06.2014 A.Ribon Use the current best elastic models.
37//
38//----------------------------------------------------------------------------
39//
40
42
43#include "G4SystemOfUnits.hh"
45#include "G4ProcessManager.hh"
46
47#include "G4MesonConstructor.hh"
49#include "G4IonConstructor.hh"
50#include "G4Neutron.hh"
51
53#include "G4HadronElastic.hh"
55#include "G4AntiNuclElastic.hh"
56#include "G4DiffuseElastic.hh"
58
61#include "G4BGGPionElasticXS.hh"
62#include "G4NeutronElasticXS.hh"
71
72#include "G4LMsdGenerator.hh"
73#include "G4DiffElasticRatio.hh"
74#include "G4AutoDelete.hh"
75
77#include "G4HadronicBuilder.hh"
78#include "G4HadParticles.hh"
79
80
81// factory
83//
85
86G4ThreadLocal G4DiffElasticRatio* G4HadronHElasticPhysics::diffRatio = nullptr;
87
89 : G4VPhysicsConstructor( "hElastic_BEST" ), verbose( ver ),
90 fDiffraction(diffraction)
91{
92 if ( verbose > 1 ) {
93 G4cout << "### G4HadronHElasticPhysics: " << GetPhysicsName()
94 << " low-mass diffraction: " << fDiffraction << G4endl;
95 }
96}
97
99
100
102 // G4cout << "G4HadronElasticPhysics::ConstructParticle" << G4endl;
103 G4MesonConstructor pMesonConstructor;
104 pMesonConstructor.ConstructParticle();
105
106 G4BaryonConstructor pBaryonConstructor;
107 pBaryonConstructor.ConstructParticle();
108
109 // Construct light ions
110 G4IonConstructor pConstructor;
111 pConstructor.ConstructParticle();
112}
113
115
116 const G4double elimitDiffuse = 0.0;
117 const G4double elimitAntiNuc = 100.0*MeV;
118 const G4double delta = 0.1*MeV;
119
120 if ( verbose > 1 ) {
121 G4cout << "### HadronHElasticPhysics::ConstructProcess: lower energy limit for DiffuseElastic : "
122 << elimitDiffuse/GeV << " GeV" << G4endl
123 << " transition energy for anti-nuclei : "
124 << elimitAntiNuc/GeV << " GeV" << G4endl;
125 }
126
128 anuc->SetMinEnergy( elimitAntiNuc );
129 G4CrossSectionElastic* anucxs =
131
132 G4HadronElastic* lhep = new G4HadronElastic();
133 lhep->SetMaxEnergy( elimitAntiNuc + delta );
134
135 G4HadronElastic* lhep0 = new G4HadronElastic();
136
137 // Three instances of Chips elastic model: one used everywhere,
138 // one used below a energy threshold, and one used only for the
139 // hydrogen element.
142 chips2->SetMaxEnergy( elimitAntiNuc + delta );
144 const G4ElementTable* theElementTable = G4Element::GetElementTable();
145 for ( size_t i_ele = 0; i_ele < theElementTable->size(); i_ele++ ) {
146 G4Element* element = (*theElementTable)[ i_ele ];
147 if ( element->GetZ() > 1.0 ) chipsH->DeActivateFor( element );
148 }
149
150 G4NuclNuclDiffuseElastic* diffuseNuclNuclElastic = new G4NuclNuclDiffuseElastic();
151 diffuseNuclNuclElastic->SetMinEnergy( elimitDiffuse );
152
153 G4VCrossSectionDataSet* theComponentGGHadronNucleusData =
155
156 G4VCrossSectionDataSet* theComponentGGNuclNuclData =
158
159 G4LMsdGenerator* diffGen = nullptr;
160 if(fDiffraction) {
161 diffGen = new G4LMsdGenerator("LMsdDiffraction");
162 diffRatio = new G4DiffElasticRatio();
163 G4AutoDelete::Register(diffRatio);
164 }
165 G4HadronElasticProcess* hel = nullptr;
166
167 auto myParticleIterator=GetParticleIterator();
168 myParticleIterator->reset();
169 while( (*myParticleIterator)() ) {
170
171 G4ParticleDefinition* particle = myParticleIterator->value();
172 G4ProcessManager* pmanager = particle->GetProcessManager();
173 G4String pname = particle->GetParticleName();
174
175 if ( pname == "anti_lambda" ||
176 pname == "anti_sigma-" ||
177 pname == "anti_sigma0" ||
178 pname == "anti_sigma+" ||
179 pname == "anti_xi-" ||
180 pname == "anti_xi0" ||
181 pname == "anti_omega-" ||
182 pname == "lambda" ||
183 pname == "sigma-" ||
184 pname == "sigma0" ||
185 pname == "sigma+" ||
186 pname == "xi-" ||
187 pname == "xi0" ||
188 pname == "omega-"
189 ) {
190 hel = new G4HadronElasticProcess();
191 hel->AddDataSet( theComponentGGHadronNucleusData );
192 hel->RegisterMe( chips1 );
193 pmanager->AddDiscreteProcess( hel );
194 if ( verbose > 1 ) {
195 G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
196 << " added for " << particle->GetParticleName() << G4endl;
197 }
198
199 } else if ( pname == "proton" ) {
200 hel = new G4HadronElasticProcess();
201 hel->AddDataSet( new G4BGGNucleonElasticXS( particle ) );
202 // To preserve reproducibility, a different instance of
203 // G4DiffuseElastic must be used for each particle type.
204 G4DiffuseElastic* protonDiffuseElastic = new G4DiffuseElastic();
205 protonDiffuseElastic->SetMinEnergy( elimitDiffuse );
206 hel->RegisterMe( chipsH ); // Use Chips only for Hydrogen element
207 hel->RegisterMe( protonDiffuseElastic );
208 pmanager->AddDiscreteProcess( hel );
209 if(fDiffraction) { hel->SetDiffraction(diffGen, diffRatio); }
210 if ( verbose > 1 ) {
211 G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
212 << " added for " << particle->GetParticleName() << G4endl;
213 }
214
215 } else if ( pname == "neutron" ) {
216 hel = new G4HadronElasticProcess();
217 hel->AddDataSet(new G4NeutronElasticXS() );
218 // To preserve reproducibility, a different instance of
219 // G4DiffuseElastic must be used for each particle type.
220 G4DiffuseElastic* neutronDiffuseElastic = new G4DiffuseElastic();
221 neutronDiffuseElastic->SetMinEnergy( elimitDiffuse );
222 hel->RegisterMe( chipsH ); // Use Chips only for Hydrogen element
223 hel->RegisterMe( neutronDiffuseElastic );
224 pmanager->AddDiscreteProcess( hel );
225 if(fDiffraction) { hel->SetDiffraction(diffGen, diffRatio); }
226 if ( verbose > 1 ) {
227 G4cout << "### HadronElasticPhysics: "
228 << hel->GetProcessName()
229 << " added for " << particle->GetParticleName() << G4endl;
230 }
231
232 } else if ( pname == "pi-" || pname == "pi+" ) {
233 hel = new G4HadronElasticProcess();
234 hel->AddDataSet( new G4BGGPionElasticXS( particle ) );
235 // To preserve reproducibility, a different instance of
236 // G4DiffuseElastic must be used for each particle type.
237 G4DiffuseElastic* dElastic = new G4DiffuseElastic();
238 dElastic->SetMinEnergy( elimitDiffuse );
239 hel->RegisterMe( chipsH ); // Use Chips only for Hydrogen element
240 hel->RegisterMe( dElastic );
241 pmanager->AddDiscreteProcess( hel );
242 if(fDiffraction) { hel->SetDiffraction(diffGen, diffRatio); }
243 if ( verbose > 1 ) {
244 G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
245 << " added for " << particle->GetParticleName() << G4endl;
246 }
247
248 } else if ( pname == "kaon-" ||
249 pname == "kaon+" ||
250 pname == "kaon0S" ||
251 pname == "kaon0L"
252 ) {
253 hel = new G4HadronElasticProcess();
254 //AR-14Aug2017 : Replaced Chips elastic kaon cross sections with
255 // Grichine's Glauber-Gribov ones. In this way, the
256 // total (elastic + inelastic) kaon cross sections
257 // are consistent with the PDG ones.
258 // For the time being, kept Chips elastic as
259 // final-state model.
260 hel->AddDataSet( theComponentGGHadronNucleusData );
261 hel->RegisterMe( chips1 );
262 pmanager->AddDiscreteProcess( hel );
263 if(fDiffraction) { hel->SetDiffraction(diffGen, diffRatio); }
264 if ( verbose > 1 ) {
265 G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
266 << " added for " << particle->GetParticleName() << G4endl;
267 }
268
269 } else if(pname == "alpha" ||
270 pname == "deuteron" ||
271 pname == "triton" ||
272 pname == "He3"
273 ) {
274 hel = new G4HadronElasticProcess();
275 hel->AddDataSet(theComponentGGNuclNuclData);
276 hel->RegisterMe(lhep0);
277 pmanager->AddDiscreteProcess(hel);
278 if(verbose > 1) {
279 G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
280 << " added for " << particle->GetParticleName() << G4endl;
281 }
282
283 } else if ( pname == "anti_proton" || pname == "anti_neutron" ) {
284 hel = new G4HadronElasticProcess();
285 hel->AddDataSet( anucxs );
286 hel->RegisterMe( chips2 );
287 hel->RegisterMe( anuc );
288 pmanager->AddDiscreteProcess( hel );
289 if ( verbose > 1 ) {
290 G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
291 << " added for " << particle->GetParticleName() << G4endl;
292 }
293
294 } else if ( pname == "anti_deuteron" ||
295 pname == "anti_triton" ||
296 pname == "anti_He3" ||
297 pname == "anti_alpha"
298 ) {
299 hel = new G4HadronElasticProcess();
300 hel->AddDataSet( anucxs );
301 hel->RegisterMe( lhep );
302 hel->RegisterMe( anuc );
303 pmanager->AddDiscreteProcess( hel );
304 if ( verbose > 1 ) {
305 G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
306 << " added for " << particle->GetParticleName() << G4endl;
307 }
308 }
309 }
310
311 // Charm and bottom hadrons
312 if ( G4HadronicParameters::Instance()->EnableBCParticles() ) {
314 }
315
316}
std::vector< G4Element * > G4ElementTable
#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:57
G4GLOB_DLL std::ostream G4cout
G4ComponentAntiNuclNuclearXS * GetComponentCrossSection()
static void ConstructParticle()
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:397
G4double GetZ() const
Definition: G4Element.hh:130
static const std::vector< G4int > & GetBCHadrons()
void SetDiffraction(G4HadronicInteraction *, G4VCrossSectionRatio *)
G4HadronHElasticPhysics(G4int ver=0, G4bool diffraction=false)
static void BuildElastic(const std::vector< G4int > &particleList)
void DeActivateFor(const G4Material *aMaterial)
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
static G4HadronicParameters * Instance()
void AddDataSet(G4VCrossSectionDataSet *aDataSet)
void RegisterMe(G4HadronicInteraction *a)
static void ConstructParticle()
static void ConstructParticle()
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleName() const
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
const G4String & GetPhysicsName() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382
void Register(T *inst)
Definition: G4AutoDelete.hh:65
#define G4ThreadLocal
Definition: tls.hh:77