Geant4  10.01.p02
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 87072 2014-11-24 14:06:09Z 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"
77 
78 #include "G4LMsdGenerator.hh"
79 #include "G4DiffElasticRatio.hh"
80 
81 // factory
83 //
85 
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 
114 
116 
117  if ( wasActivated ) return;
118  wasActivated = true;
119 
120  const G4double elimitDiffuse = 0.0;
121  const G4double elimitAntiNuc = 100.0*MeV;
122  const G4double delta = 0.1*MeV;
123 
124  if ( verbose > 1 ) {
125  G4cout << "### HadronHElasticPhysics::ConstructProcess: lower energy limit for DiffuseElastic : "
126  << elimitDiffuse/GeV << " GeV" << G4endl
127  << " transition energy for anti-nuclei : "
128  << elimitAntiNuc/GeV << " GeV" << G4endl;
129  }
130 
131  G4AntiNuclElastic* anuc = new G4AntiNuclElastic();
132  anuc->SetMinEnergy( elimitAntiNuc );
133  G4CrossSectionElastic* anucxs =
135 
136  G4HadronElastic* lhep = new G4HadronElastic();
137  lhep->SetMaxEnergy( elimitAntiNuc + delta );
138 
139  // Three instances of Chips elastic model: one used everywhere,
140  // one used below a energy threshold, and one used only for the
141  // hydrogen element.
144  chips2->SetMaxEnergy( elimitAntiNuc + delta );
146  const G4ElementTable* theElementTable = G4Element::GetElementTable();
147  for ( size_t i_ele = 0; i_ele < theElementTable->size(); i_ele++ ) {
148  G4Element* element = (*theElementTable)[ i_ele ];
149  if ( element->GetZ() > 1.0 ) chipsH->DeActivateFor( element );
150  }
151 
152  G4NuclNuclDiffuseElastic* diffuseNuclNuclElastic = new G4NuclNuclDiffuseElastic();
153  diffuseNuclNuclElastic->SetMinEnergy( elimitDiffuse );
154 
155  G4VCrossSectionDataSet* theComponentGGNuclNuclData =
157 
158  G4LMsdGenerator* diffGen = 0;
159  G4DiffElasticRatio* diffRatio = 0;
160  if(fDiffraction) {
161  diffGen = new G4LMsdGenerator("LMsdDiffraction");
162  diffRatio = new G4DiffElasticRatio();
163  }
164 
165  aParticleIterator->reset();
166  while( (*aParticleIterator)() ) {
167 
168  G4ParticleDefinition* particle = aParticleIterator->value();
169  G4ProcessManager* pmanager = particle->GetProcessManager();
170  G4String pname = particle->GetParticleName();
171 
172  if ( pname == "anti_lambda" ||
173  pname == "anti_sigma-" ||
174  pname == "anti_sigma0" ||
175  pname == "anti_sigma+" ||
176  pname == "anti_xi-" ||
177  pname == "anti_xi0" ||
178  pname == "anti_omega-"
179  ) {
182  hel->RegisterMe( chips1 );
183  pmanager->AddDiscreteProcess( hel );
184  if ( verbose > 1 ) {
185  G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
186  << " added for " << particle->GetParticleName() << G4endl;
187  }
188 
189  } else if ( pname == "lambda" ||
190  pname == "sigma-" ||
191  pname == "sigma0" ||
192  pname == "sigma+" ||
193  pname == "xi-" ||
194  pname == "xi0" ||
195  pname == "omega-"
196  ) {
199  hel->RegisterMe( chips1 );
200  pmanager->AddDiscreteProcess( hel );
201  if ( verbose > 1 ) {
202  G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
203  << " added for " << particle->GetParticleName() << G4endl;
204  }
205 
206  } else if ( pname == "proton" ) {
208  hel->AddDataSet( new G4BGGNucleonElasticXS( particle ) );
209  // To preserve reproducibility, a different instance of
210  // G4DiffuseElastic must be used for each particle type.
211  G4DiffuseElastic* protonDiffuseElastic = new G4DiffuseElastic();
212  protonDiffuseElastic->SetMinEnergy( elimitDiffuse );
213  hel->RegisterMe( chipsH ); // Use Chips only for Hydrogen element
214  hel->RegisterMe( protonDiffuseElastic );
215  pmanager->AddDiscreteProcess( hel );
216  if(fDiffraction) { hel->SetDiffraction(diffGen, diffRatio); }
217  if ( verbose > 1 ) {
218  G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
219  << " added for " << particle->GetParticleName() << G4endl;
220  }
221 
222  } else if ( pname == "neutron" ) {
225  // To preserve reproducibility, a different instance of
226  // G4DiffuseElastic must be used for each particle type.
227  G4DiffuseElastic* neutronDiffuseElastic = new G4DiffuseElastic();
228  neutronDiffuseElastic->SetMinEnergy( elimitDiffuse );
229  hel->RegisterMe( chipsH ); // Use Chips only for Hydrogen element
230  hel->RegisterMe( neutronDiffuseElastic );
231  pmanager->AddDiscreteProcess( hel );
232  if(fDiffraction) { hel->SetDiffraction(diffGen, diffRatio); }
233  if ( verbose > 1 ) {
234  G4cout << "### HadronElasticPhysics: "
235  << hel->GetProcessName()
236  << " added for " << particle->GetParticleName() << G4endl;
237  }
238 
239  } else if ( pname == "pi-" ) {
241  hel->AddDataSet( new G4BGGPionElasticXS( particle ) );
242  // To preserve reproducibility, a different instance of
243  // G4DiffuseElastic must be used for each particle type.
244  G4DiffuseElastic* pionMinusDiffuseElastic = new G4DiffuseElastic();
245  pionMinusDiffuseElastic->SetMinEnergy( elimitDiffuse );
246  hel->RegisterMe( chipsH ); // Use Chips only for Hydrogen element
247  hel->RegisterMe( pionMinusDiffuseElastic );
248  pmanager->AddDiscreteProcess( hel );
249  if(fDiffraction) { hel->SetDiffraction(diffGen, diffRatio); }
250  if ( verbose > 1 ) {
251  G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
252  << " added for " << particle->GetParticleName() << G4endl;
253  }
254 
255  } else if ( pname == "pi+" ) {
257  hel->AddDataSet( new G4BGGPionElasticXS( particle ) );
258  // To preserve reproducibility, a different instance of
259  // G4DiffuseElastic must be used for each particle type.
260  G4DiffuseElastic* pionPlusDiffuseElastic = new G4DiffuseElastic();
261  hel->RegisterMe( chipsH ); // Use Chips only for Hydrogen element
262  pionPlusDiffuseElastic->SetMinEnergy( elimitDiffuse );
263  hel->RegisterMe( pionPlusDiffuseElastic );
264  pmanager->AddDiscreteProcess( hel );
265  if(fDiffraction) { hel->SetDiffraction(diffGen, diffRatio); }
266  if ( verbose > 1 ) {
267  G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
268  << " added for " << particle->GetParticleName() << G4endl;
269  }
270 
271  } else if ( pname == "kaon-" ||
272  pname == "kaon+" ||
273  pname == "kaon0S" ||
274  pname == "kaon0L"
275  ) {
277  if ( pname == "kaon-" ) {
279  } else if ( pname == "kaon+" ) {
281  } else {
283  }
284  hel->RegisterMe( chips1 );
285  pmanager->AddDiscreteProcess( hel );
286  if(fDiffraction) { hel->SetDiffraction(diffGen, diffRatio); }
287  if ( verbose > 1 ) {
288  G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
289  << " added for " << particle->GetParticleName() << G4endl;
290  }
291 
292  } else if (
293  pname == "deuteron" ||
294  pname == "triton" ||
295  pname == "He3" ||
296  pname == "alpha"
297  ) {
299  hel->AddDataSet( theComponentGGNuclNuclData );
300  // To preserve reproducibility, replace temporarily
301  // G4NuclNuclDiffuseElastic with the Gheisha elastic model.
302  //hel->RegisterMe( diffuseNuclNuclElastic );
303  G4HadronElastic* lhepLightIon = new G4HadronElastic();
304  hel->RegisterMe( lhepLightIon );
305  pmanager->AddDiscreteProcess( hel );
306  if ( verbose > 1 ) {
307  G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
308  << " added for " << particle->GetParticleName() << G4endl;
309  }
310 
311  } else if ( pname == "anti_proton" || pname == "anti_neutron" ) {
313  hel->AddDataSet( anucxs );
314  hel->RegisterMe( chips2 );
315  hel->RegisterMe( anuc );
316  pmanager->AddDiscreteProcess( hel );
317  if ( verbose > 1 ) {
318  G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
319  << " added for " << particle->GetParticleName() << G4endl;
320  }
321 
322  } else if ( pname == "anti_deuteron" ||
323  pname == "anti_triton" ||
324  pname == "anti_He3" ||
325  pname == "anti_alpha"
326  ) {
328  hel->AddDataSet( anucxs );
329  hel->RegisterMe( lhep );
330  hel->RegisterMe( anuc );
331  pmanager->AddDiscreteProcess( hel );
332  if ( verbose > 1 ) {
333  G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
334  << " added for " << particle->GetParticleName() << G4endl;
335  }
336 
337  } else if ( pname == "GenericIon" ) {
338  // To preserve reproducibility, disable temporarily
339  // G4NuclNuclDiffuseElastic.
340  //G4HadronElasticProcess* hel = new G4HadronElasticProcess();
341  //hel->AddDataSet( theComponentGGNuclNuclData );
342  //hel->RegisterMe( diffuseNuclNuclElastic );
343  //pmanager->AddDiscreteProcess( hel );
344  //if ( verbose > 1 ) {
345  // G4cout << "### HadronElasticPhysics: " << hel->GetProcessName()
346  // << " added for " << particle->GetParticleName() << G4endl;
347  //}
348 
349  }
350 
351  }
352 
353 }
static const double MeV
Definition: G4SIunits.hh:193
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()
G4ProcessManager * GetProcessManager() const
int G4int
Definition: G4Types.hh:78
G4ComponentAntiNuclNuclearXS * GetComponentCrossSection()
const G4String & GetParticleName() const
static void ConstructParticle()
void RegisterMe(G4HadronicInteraction *a)
static const char * Default_Name()
void SetMinEnergy(G4double anEnergy)
void AddDataSet(G4VCrossSectionDataSet *aDataSet)
G4GLOB_DLL std::ostream G4cout
static G4ThreadLocal G4bool wasActivated
bool G4bool
Definition: G4Types.hh:79
#define aParticleIterator
static G4CrossSectionDataSetRegistry * Instance()
G4_DECLARE_PHYSCONSTR_FACTORY(G4HadronHElasticPhysics)
static const char * Default_Name()
static const double GeV
Definition: G4SIunits.hh:196
const G4String & GetPhysicsName() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4HadronHElasticPhysics(G4int ver=0, G4bool diffraction=false)
static const char * Default_Name()
void SetMaxEnergy(const G4double anEnergy)
#define G4endl
Definition: G4ios.hh:61
static const char * Default_Name()
std::vector< G4Element * > G4ElementTable
double G4double
Definition: G4Types.hh:76
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
void DeActivateFor(const G4Material *aMaterial)