Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4BGGNucleonInelasticXS.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: G4BGGNucleonInelasticXS.cc 93682 2015-10-28 10:09:49Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4BGGNucleonInelasticXS
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 13.03.2007
38 // Modifications:
39 //
40 //
41 // -------------------------------------------------------------------
42 //
43 
45 #include "G4SystemOfUnits.hh"
48 #include "G4HadronNucleonXsc.hh"
50 #include "G4Proton.hh"
51 #include "G4Neutron.hh"
52 #include "G4NistManager.hh"
53 #include "G4Material.hh"
54 #include "G4Element.hh"
55 #include "G4Isotope.hh"
56 #include "G4Log.hh"
57 #include "G4Exp.hh"
58 
60 
61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
62 
63 const G4double llog10 = G4Log(10.);
64 
66  : G4VCrossSectionDataSet("Barashenkov-Glauber")
67 {
68  verboseLevel = 0;
69  fGlauberEnergy = 91.*GeV;
70  fLowEnergy = 14.*MeV;
71  fHighEnergy = 5.*GeV;
72  fSAIDHighEnergyLimit = 1.3*GeV;
73  fLowestXSection = millibarn;
74  for (G4int i = 0; i < 93; ++i) {
75  theGlauberFac[i] = 0.0;
76  theCoulombFac[i] = 0.0;
77  theA[i] = 1;
78  }
79  fNucleon = 0;
80  fGlauber = 0;
81  fHadron = 0;
82  fSAID = 0;
83 
84  particle = p;
85  theProton= G4Proton::Proton();
86  isProton = false;
87  isInitialized = false;
88 }
89 
90 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
91 
93 {
94  delete fSAID;
95  delete fHadron;
96  // The cross section registry will delete fNucleon
97  delete fGlauber;
98 }
99 
100 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
101 
103  G4int, const G4Material*)
104 {
105  return true;
106 }
107 
108 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
109 
111  G4int Z, G4int,
112  const G4Element*,
113  const G4Material*)
114 {
115  return (1 == Z);
116 }
117 
118 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
119 
120 G4double
122  G4int ZZ, const G4Material*)
123 {
124  // this method should be called only for Z > 1
125 
126  G4double cross = 0.0;
127  G4double ekin = dp->GetKineticEnergy();
128  G4int Z = ZZ;
129  if(1 == Z) {
130  cross = 1.0115*GetIsoCrossSection(dp,1,1);
131  } else if(2 == Z) {
132  if(ekin > fGlauberEnergy) {
133  cross = theGlauberFac[Z]*fGlauber->GetInelasticGlauberGribov(dp, Z, theA[Z]);
134  } else {
135  cross = fNucleon->GetElementCrossSection(dp, Z);
136  }
137 
138  } else {
139  if(Z > 92) { Z = 92; }
140 
141  if(ekin <= fLowEnergy) {
142  cross = theCoulombFac[Z]*CoulombFactor(ekin, Z);
143  } else if(ekin > fGlauberEnergy) {
144  cross = theGlauberFac[Z]*fGlauber->GetInelasticGlauberGribov(dp, Z, theA[Z]);
145  } else {
146  cross = fNucleon->GetElementCrossSection(dp, Z);
147  }
148  }
149 
150  if(verboseLevel > 1) {
151  G4cout << "G4BGGNucleonInelasticXS::GetCrossSection for "
152  << dp->GetDefinition()->GetParticleName()
153  << " Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
154  << " in nucleus Z= " << Z << " A= " << theA[Z]
155  << " XS(b)= " << cross/barn
156  << G4endl;
157  }
158  //AR-18Dec2013 if(cross <= fLowestXSection) { cross = 0.0; }
159  return cross;
160 }
161 
162 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
163 
164 G4double
166  G4int Z, G4int A,
167  const G4Isotope*,
168  const G4Element*,
169  const G4Material*)
170 {
171  // this method should be called only for Z = 1
172 
173  G4double cross = 0.0;
174  G4double ekin = dp->GetKineticEnergy();
175 
176  if(ekin <= fSAIDHighEnergyLimit) {
177  cross = fSAID->GetInelasticIsotopeCrossSection(particle, ekin, 1, 1);
178  } else if(ekin < fHighEnergy) {
179  fHadron->GetHadronNucleonXscNS(dp, theProton);
180  cross = (theCoulombFac[0]/ekin + 1)*fHadron->GetInelasticHadronNucleonXsc();
181  } else {
182  fHadron->GetHadronNucleonXscPDG(dp, theProton);
183  cross = (theCoulombFac[1]/ekin + 1)*fHadron->GetInelasticHadronNucleonXsc();
184  }
185  cross *= A;
186 
187  if(verboseLevel > 1) {
188  G4cout << "G4BGGNucleonInelasticXS::GetCrossSection for "
189  << dp->GetDefinition()->GetParticleName()
190  << " Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
191  << " in nucleus Z= " << Z << " A= " << A
192  << " XS(b)= " << cross/barn
193  << G4endl;
194  }
195  //AR-18Dec2013 if(cross <= fLowestXSection) { cross = 0.0; }
196  return cross;
197 }
198 
199 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
200 
202 {
203  if(&p == theProton || &p == G4Neutron::Neutron()) {
204  particle = &p;
205  } else {
206  G4cout << "### G4BGGNucleonInelasticXS WARNING: is not applicable to "
207  << p.GetParticleName()
208  << G4endl;
209  throw G4HadronicException(__FILE__, __LINE__,
210  "G4BGGNucleonElasticXS::BuildPhysicsTable is used for wrong particle");
211  return;
212  }
213 
214  if(isInitialized) { return; }
215  isInitialized = true;
216 
218  fGlauber = new G4ComponentGGHadronNucleusXsc();
219  fHadron = new G4HadronNucleonXsc();
220  fSAID = new G4ComponentSAIDTotalXS();
221 
222  fNucleon->BuildPhysicsTable(*particle);
223  fGlauber->BuildPhysicsTable(*particle);
224 
225  if(particle == theProton) {
226  isProton = true;
227  fSAIDHighEnergyLimit = 2*GeV;
228  fHighEnergy = 2*GeV;
229  }
230 
231  G4ThreeVector mom(0.0,0.0,1.0);
232  G4DynamicParticle dp(particle, mom, fGlauberEnergy);
233 
235  G4int A;
236 
237  G4double csup, csdn;
238 
239  if(verboseLevel > 0) {
240  G4cout << "### G4BGGNucleonInelasticXS::Initialise for "
241  << particle->GetParticleName() << G4endl;
242  }
243  for(G4int iz=2; iz<93; iz++) {
244 
245  A = G4lrint(nist->GetAtomicMassAmu(iz));
246  theA[iz] = A;
247 
248  csup = fGlauber->GetInelasticGlauberGribov(&dp, iz, A);
249  csdn = fNucleon->GetElementCrossSection(&dp, iz);
250 
251  theGlauberFac[iz] = csdn/csup;
252  if(verboseLevel > 0) {
253  G4cout << "Z= " << iz << " A= " << A
254  << " GlauberFactor= " << theGlauberFac[iz] << G4endl;
255  }
256  }
257  //const G4Material* mat = 0;
258 
259  dp.SetKineticEnergy(fSAIDHighEnergyLimit);
260  fHadron->GetHadronNucleonXscNS(&dp, theProton);
261  theCoulombFac[0] = fSAIDHighEnergyLimit*
262  (fSAID->GetInelasticIsotopeCrossSection(particle,fSAIDHighEnergyLimit,1,1)
263  /fHadron->GetInelasticHadronNucleonXsc() - 1);
264 
265  //G4cout << "Z=1 E(GeV)= " << fSAIDHighEnergyLimit/GeV
266  // << " xsNS(b)= " << fHadron->GetInelasticHadronNucleonXsc()/barn;
267  fHadron->GetHadronNucleonXscPDG(&dp, theProton);
268  //G4cout << " xsPDG(b)= " << fHadron->GetInelasticHadronNucleonXsc()/barn;
269  //G4cout << " xsSAID(b)= " << fSAID->GetInelasticIsotopeCrossSection(particle,fSAIDHighEnergyLimit,1,1)/barn << G4endl;
270 
271  dp.SetKineticEnergy(fHighEnergy);
272  fHadron->GetHadronNucleonXscPDG(&dp, theProton);
273  G4double x = fHadron->GetInelasticHadronNucleonXsc();
274 
275  //G4cout << "Z=1 E(GeV)= " << fHighEnergy/GeV
276  // << " xsPDG(b)= " << fHadron->GetInelasticHadronNucleonXsc()/barn;
277 
278  fHadron->GetHadronNucleonXscNS(&dp, theProton);
279  theCoulombFac[1] = fHighEnergy*((theCoulombFac[0]/fHighEnergy + 1)
280  *fHadron->GetInelasticHadronNucleonXsc()/x - 1);
281 
282  fHadron->GetHadronNucleonXscNS(&dp, theProton);
283  //G4cout <<" xsNS(b)= "<<fHadron->GetInelasticHadronNucleonXsc()/barn<<G4endl;
284 
285  if(verboseLevel > 0) {
286  G4cout << "Z=1 A=1" << " CoulombFactor[0]= " << theCoulombFac[0]
287  << " CoulombFactor[1]= " << theCoulombFac[1] << G4endl;
288  }
289  theCoulombFac[2] = 1.0;
290 
291  dp.SetKineticEnergy(fLowEnergy);
292  for(G4int iz=3; iz<93; iz++) {
293  theCoulombFac[iz] =
294  fNucleon->GetElementCrossSection(&dp, iz)/CoulombFactor(fLowEnergy, iz);
295 
296  if(verboseLevel > 0) {
297  G4cout << "Z= " << iz << " A= " << theA[iz]
298  << " CoulombFactor= " << theCoulombFac[iz] << G4endl;
299  }
300  }
301 }
302 
303 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
304 
305 G4double G4BGGNucleonInelasticXS::CoulombFactor(G4double kinEnergy, G4int Z)
306 {
307  G4double res= 0.0;
308  if(kinEnergy <= 0.0) { return res; }
309  else if (Z <= 1) { return kinEnergy*kinEnergy; }
310 
311  G4double elog = G4Log(kinEnergy/GeV)/llog10;
312  G4double aa = theA[Z];
313 
314  // from G4ProtonInelasticCrossSection
315  if(isProton) {
316 
317  G4double ff1 = 5.6 - 0.016*aa; // slope of the drop at medium energies.
318  G4double ff2 = 1.37 + 1.37/aa; // start of the slope.
319  G4double ff3 = 0.8 + 18./aa - 0.002*aa; // stephight
320  res = 1.0 + ff3*(1.0 - (1.0/(1+G4Exp(-ff1*(elog + ff2)))));
321 
322  ff1 = 8. - 8./aa - 0.008*aa; // slope of the rise
323  ff2 = 2.34 - 5.4/aa - 0.0028*aa; // start of the rise
324  res /= (1.0 + G4Exp(-ff1*(elog + ff2)));
325 
326  } else {
327 
328  // from G4NeutronInelasticCrossSection
329  G4double p3 = 0.6 + 13./aa - 0.0005*aa;
330  G4double p4 = 7.2449 - 0.018242*aa;
331  G4double p5 = 1.36 + 1.8/aa + 0.0005*aa;
332  G4double p6 = 1. + 200./aa + 0.02*aa;
333  G4double p7 = 3.0 - (aa-70.)*(aa-200.)/11000.;
334 
335  G4double firstexp = G4Exp(-p4*(elog + p5));
336  G4double secondexp = G4Exp(-p6*(elog + p7));
337 
338  res = (1.+p3*firstexp/(1. + firstexp))/(1. + secondexp);
339 
340  }
341  return res;
342 }
343 
344 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
345 
346 void G4BGGNucleonInelasticXS::CrossSectionDescription(std::ostream& outFile) const
347 {
348  outFile << "The Barashenkov-Glauber-Gribov cross section calculates inelastic\n"
349  << "scattering of protons and neutrons from nuclei using the\n"
350  << "Barashenkov parameterization below 91 GeV and the Glauber-Gribov\n"
351  << "parameterization above 91 GeV. It uses the G4HadronNucleonXsc\n"
352  << "cross section component for hydrogen targets, and the\n"
353  << "G4ComponentGGHadronNucleusXsc component for other targets.\n";
354 }
355 
356 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4VCrossSectionDataSet * GetCrossSectionDataSet(const G4String &name, G4bool warning=true)
virtual G4double GetInelasticIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)
G4double GetKineticEnergy() const
const char * p
Definition: xmltok.h:285
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
virtual G4double GetElementCrossSection(const G4DynamicParticle *aParticle, G4int Z, const G4Material *mat=0)
G4ParticleDefinition * GetDefinition() const
const G4double llog10
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
const G4String & GetParticleName() const
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4ParticleDefinition *)
G4BGGNucleonInelasticXS(const G4ParticleDefinition *)
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
virtual void CrossSectionDescription(std::ostream &) const
bool G4bool
Definition: G4Types.hh:79
static G4CrossSectionDataSetRegistry * Instance()
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
void SetKineticEnergy(G4double aEnergy)
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
static constexpr double GeV
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
int G4lrint(double ad)
Definition: templates.hh:163
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4double GetAtomicMassAmu(const G4String &symb) const
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
static constexpr double GeV
Definition: G4SIunits.hh:217
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
double G4double
Definition: G4Types.hh:76
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
static constexpr double barn
Definition: G4SIunits.hh:105
virtual G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=0, const G4Material *mat=0)
static constexpr double millibarn
Definition: G4SIunits.hh:106
G4double GetInelasticHadronNucleonXsc()