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