Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 // $Id: G4HadronHElasticPhysics.cc 106611 2017-10-16 14:51:20Z gcosmo $
27 //
28 //---------------------------------------------------------------------------
29 //
30 // ClassName: G4HadronHElasticPhysics
31 //
32 // Author: 23 November 2006 V. Ivanchenko
33 //
34 // Modified:
35 // 21.03.07 (V.Ivanchenko) Use G4BGGNucleonElasticXS and G4BGGPionElasticXS;
36 // Reduce thresholds for HE and Q-models to zero
37 // 03.06.2010 V.Ivanchenko cleanup constructors and ConstructProcess method
38 // 06.06.2014 A.Ribon Use the current best elastic models.
39 //
40 //----------------------------------------------------------------------------
41 //
42 // CHIPS for sampling scattering for p and n
43 // Glauber model for samplimg of high energy pi+- (E > 1GeV)
44 // LHEP sampling model for the other particle
45 // BBG cross sections for p, n and pi+-
46 // LHEP cross sections for other particles
47 
49 
50 #include "G4SystemOfUnits.hh"
51 #include "G4ParticleDefinition.hh"
52 #include "G4ProcessManager.hh"
53 
54 #include "G4MesonConstructor.hh"
55 #include "G4BaryonConstructor.hh"
56 #include "G4IonConstructor.hh"
57 #include "G4Neutron.hh"
58 
60 #include "G4HadronElastic.hh"
61 #include "G4ChipsElasticModel.hh"
62 #include "G4AntiNuclElastic.hh"
63 #include "G4DiffuseElastic.hh"
65 
66 #include "G4CrossSectionElastic.hh"
67 #include "G4BGGNucleonElasticXS.hh"
68 #include "G4BGGPionElasticXS.hh"
69 #include "G4NeutronElasticXS.hh"
78 
79 #include "G4LMsdGenerator.hh"
80 #include "G4DiffElasticRatio.hh"
81 #include "G4AutoDelete.hh"
82 
83 // factory
85 //
87 
88 G4ThreadLocal G4DiffElasticRatio* G4HadronHElasticPhysics::diffRatio = 0;
89 
91  : G4VPhysicsConstructor( "hElastic_BEST" ), verbose( ver ),
92  fDiffraction(diffraction)
93 {
94  if ( verbose > 1 ) {
95  G4cout << "### G4HadronHElasticPhysics: " << GetPhysicsName()
96  << " low-mass diffraction: " << fDiffraction << G4endl;
97  }
98 }
99 
101 
102 
104  // G4cout << "G4HadronElasticPhysics::ConstructParticle" << G4endl;
105  G4MesonConstructor pMesonConstructor;
106  pMesonConstructor.ConstructParticle();
107 
108  G4BaryonConstructor pBaryonConstructor;
109  pBaryonConstructor.ConstructParticle();
110 
111  // Construct light ions
112  G4IonConstructor pConstructor;
113  pConstructor.ConstructParticle();
114 }
115 
117 
118  const G4double elimitDiffuse = 0.0;
119  const G4double elimitAntiNuc = 100.0*MeV;
120  const G4double delta = 0.1*MeV;
121 
122  if ( verbose > 1 ) {
123  G4cout << "### HadronHElasticPhysics::ConstructProcess: lower energy limit for DiffuseElastic : "
124  << elimitDiffuse/GeV << " GeV" << G4endl
125  << " transition energy for anti-nuclei : "
126  << elimitAntiNuc/GeV << " GeV" << G4endl;
127  }
128 
129  G4AntiNuclElastic* anuc = new G4AntiNuclElastic();
130  anuc->SetMinEnergy( elimitAntiNuc );
131  G4CrossSectionElastic* anucxs =
133 
134  G4HadronElastic* lhep = new G4HadronElastic();
135  lhep->SetMaxEnergy( elimitAntiNuc + delta );
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 = 0;
160  if(fDiffraction) {
161  diffGen = new G4LMsdGenerator("LMsdDiffraction");
162  diffRatio = new G4DiffElasticRatio();
163  G4AutoDelete::Register(diffRatio);
164  }
165 
166  auto myParticleIterator=GetParticleIterator();
167  myParticleIterator->reset();
168  while( (*myParticleIterator)() ) {
169 
170  G4ParticleDefinition* particle = myParticleIterator->value();
171  G4ProcessManager* pmanager = particle->GetProcessManager();
172  G4String pname = particle->GetParticleName();
173 
174  if ( pname == "anti_lambda" ||
175  pname == "anti_sigma-" ||
176  pname == "anti_sigma0" ||
177  pname == "anti_sigma+" ||
178  pname == "anti_xi-" ||
179  pname == "anti_xi0" ||
180  pname == "anti_omega-"
181  ) {
184  hel->RegisterMe( chips1 );
185  pmanager->AddDiscreteProcess( hel );
186  if ( verbose > 1 ) {
187  G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
188  << " added for " << particle->GetParticleName() << G4endl;
189  }
190 
191  } else if ( pname == "lambda" ||
192  pname == "sigma-" ||
193  pname == "sigma0" ||
194  pname == "sigma+" ||
195  pname == "xi-" ||
196  pname == "xi0" ||
197  pname == "omega-"
198  ) {
201  hel->RegisterMe( chips1 );
202  pmanager->AddDiscreteProcess( hel );
203  if ( verbose > 1 ) {
204  G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
205  << " added for " << particle->GetParticleName() << G4endl;
206  }
207 
208  } else if ( pname == "proton" ) {
210  hel->AddDataSet( new G4BGGNucleonElasticXS( particle ) );
211  // To preserve reproducibility, a different instance of
212  // G4DiffuseElastic must be used for each particle type.
213  G4DiffuseElastic* protonDiffuseElastic = new G4DiffuseElastic();
214  protonDiffuseElastic->SetMinEnergy( elimitDiffuse );
215  hel->RegisterMe( chipsH ); // Use Chips only for Hydrogen element
216  hel->RegisterMe( protonDiffuseElastic );
217  pmanager->AddDiscreteProcess( hel );
218  if(fDiffraction) { hel->SetDiffraction(diffGen, diffRatio); }
219  if ( verbose > 1 ) {
220  G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
221  << " added for " << particle->GetParticleName() << G4endl;
222  }
223 
224  } else if ( pname == "neutron" ) {
227  // To preserve reproducibility, a different instance of
228  // G4DiffuseElastic must be used for each particle type.
229  G4DiffuseElastic* neutronDiffuseElastic = new G4DiffuseElastic();
230  neutronDiffuseElastic->SetMinEnergy( elimitDiffuse );
231  hel->RegisterMe( chipsH ); // Use Chips only for Hydrogen element
232  hel->RegisterMe( neutronDiffuseElastic );
233  pmanager->AddDiscreteProcess( hel );
234  if(fDiffraction) { hel->SetDiffraction(diffGen, diffRatio); }
235  if ( verbose > 1 ) {
236  G4cout << "### HadronElasticPhysics: "
237  << hel->GetProcessName()
238  << " added for " << particle->GetParticleName() << G4endl;
239  }
240 
241  } else if ( pname == "pi-" ) {
243  hel->AddDataSet( new G4BGGPionElasticXS( particle ) );
244  // To preserve reproducibility, a different instance of
245  // G4DiffuseElastic must be used for each particle type.
246  G4DiffuseElastic* pionMinusDiffuseElastic = new G4DiffuseElastic();
247  pionMinusDiffuseElastic->SetMinEnergy( elimitDiffuse );
248  hel->RegisterMe( chipsH ); // Use Chips only for Hydrogen element
249  hel->RegisterMe( pionMinusDiffuseElastic );
250  pmanager->AddDiscreteProcess( hel );
251  if(fDiffraction) { hel->SetDiffraction(diffGen, diffRatio); }
252  if ( verbose > 1 ) {
253  G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
254  << " added for " << particle->GetParticleName() << G4endl;
255  }
256 
257  } else if ( pname == "pi+" ) {
259  hel->AddDataSet( new G4BGGPionElasticXS( particle ) );
260  // To preserve reproducibility, a different instance of
261  // G4DiffuseElastic must be used for each particle type.
262  G4DiffuseElastic* pionPlusDiffuseElastic = new G4DiffuseElastic();
263  hel->RegisterMe( chipsH ); // Use Chips only for Hydrogen element
264  pionPlusDiffuseElastic->SetMinEnergy( elimitDiffuse );
265  hel->RegisterMe( pionPlusDiffuseElastic );
266  pmanager->AddDiscreteProcess( hel );
267  if(fDiffraction) { hel->SetDiffraction(diffGen, diffRatio); }
268  if ( verbose > 1 ) {
269  G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
270  << " added for " << particle->GetParticleName() << G4endl;
271  }
272 
273  } else if ( pname == "kaon-" ||
274  pname == "kaon+" ||
275  pname == "kaon0S" ||
276  pname == "kaon0L"
277  ) {
279  //AR-14Aug2017 : Replaced Chips elastic kaon cross sections with
280  // Grichine's Glauber-Gribov ones. In this way, the
281  // total (elastic + inelastic) kaon cross sections
282  // are consistent with the PDG ones.
283  // For the time being, kept Chips elastic as
284  // final-state model.
285  hel->AddDataSet( theComponentGGHadronNucleusData );
286  hel->RegisterMe( chips1 );
287  pmanager->AddDiscreteProcess( hel );
288  if(fDiffraction) { hel->SetDiffraction(diffGen, diffRatio); }
289  if ( verbose > 1 ) {
290  G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
291  << " added for " << particle->GetParticleName() << G4endl;
292  }
293 
294  } else if (
295  pname == "deuteron" ||
296  pname == "triton" ||
297  pname == "He3" ||
298  pname == "alpha"
299  ) {
301  hel->AddDataSet( theComponentGGNuclNuclData );
302  // To preserve reproducibility, replace temporarily
303  // G4NuclNuclDiffuseElastic with the Gheisha elastic model.
304  //hel->RegisterMe( diffuseNuclNuclElastic );
305  G4HadronElastic* lhepLightIon = new G4HadronElastic();
306  hel->RegisterMe( lhepLightIon );
307  pmanager->AddDiscreteProcess( hel );
308  if ( verbose > 1 ) {
309  G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
310  << " added for " << particle->GetParticleName() << G4endl;
311  }
312 
313  } else if ( pname == "anti_proton" || pname == "anti_neutron" ) {
315  hel->AddDataSet( anucxs );
316  hel->RegisterMe( chips2 );
317  hel->RegisterMe( anuc );
318  pmanager->AddDiscreteProcess( hel );
319  if ( verbose > 1 ) {
320  G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
321  << " added for " << particle->GetParticleName() << G4endl;
322  }
323 
324  } else if ( pname == "anti_deuteron" ||
325  pname == "anti_triton" ||
326  pname == "anti_He3" ||
327  pname == "anti_alpha"
328  ) {
330  hel->AddDataSet( anucxs );
331  hel->RegisterMe( lhep );
332  hel->RegisterMe( anuc );
333  pmanager->AddDiscreteProcess( hel );
334  if ( verbose > 1 ) {
335  G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
336  << " added for " << particle->GetParticleName() << G4endl;
337  }
338 
339  } else if ( pname == "GenericIon" ) {
340  // To preserve reproducibility, disable temporarily
341  // G4NuclNuclDiffuseElastic.
342  //G4HadronElasticProcess* hel = new G4HadronElasticProcess();
343  //hel->AddDataSet( theComponentGGNuclNuclData );
344  //hel->RegisterMe( diffuseNuclNuclElastic );
345  //pmanager->AddDiscreteProcess( hel );
346  //if ( verbose > 1 ) {
347  // G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
348  // << " added for " << particle->GetParticleName() << G4endl;
349  //}
350 
351  }
352 
353  }
354 
355 }
G4double GetZ() const
Definition: G4Element.hh:131
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
static void ConstructParticle()
#define G4ThreadLocal
Definition: tls.hh:89
void SetDiffraction(G4HadronicInteraction *, G4VCrossSectionRatio *)
static void ConstructParticle()
int G4int
Definition: G4Types.hh:78
G4ComponentAntiNuclNuclearXS * GetComponentCrossSection()
const G4String & GetParticleName() const
static void ConstructParticle()
void RegisterMe(G4HadronicInteraction *a)
static const char * Default_Name()
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
void SetMinEnergy(G4double anEnergy)
void AddDataSet(G4VCrossSectionDataSet *aDataSet)
void Register(T *inst)
Definition: G4AutoDelete.hh:65
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
static G4CrossSectionDataSetRegistry * Instance()
const G4String & GetPhysicsName() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
string pname
Definition: eplot.py:33
G4ProcessManager * GetProcessManager() const
G4HadronHElasticPhysics(G4int ver=0, G4bool diffraction=false)
static constexpr double GeV
Definition: G4SIunits.hh:217
static const char * Default_Name()
void SetMaxEnergy(const G4double anEnergy)
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
std::vector< G4Element * > G4ElementTable
double G4double
Definition: G4Types.hh:76
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:398
#define G4_DECLARE_PHYSCONSTR_FACTORY(physics_constructor)
void DeActivateFor(const G4Material *aMaterial)