Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups 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 78189 2013-12-04 16:34:08Z 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 fHadron;
90  delete fSAID;
91 }
92 
93 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
94 
95 G4bool
97  const G4Material*)
98 {
99  return true;
100 }
101 
102 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
103 
105  G4int Z, G4int,
106  const G4Element*,
107  const G4Material*)
108 {
109  return (1 == Z);
110 }
111 
112 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
113 
114 G4double
116  G4int ZZ, const G4Material*)
117 {
118  // this method should be called only for Z > 1
119 
120  G4double cross = 0.0;
121  G4double ekin = dp->GetKineticEnergy();
122  G4int Z = ZZ;
123  if(1 == Z) {
124  cross = 1.0115*GetIsoCrossSection(dp,1,1);
125  } else {
126  if(Z > 92) { Z = 92; }
127 
128  if(ekin <= fLowEnergy) {
129  cross = theCoulombFac[Z]*CoulombFactor(ekin, Z);
130 
131  } else if(ekin > fGlauberEnergy) {
132  cross = theGlauberFac[Z]*fGlauber->GetElasticGlauberGribov(dp, Z, theA[Z]);
133  } else {
134  cross = fNucleon->GetElasticCrossSection(dp, Z);
135  }
136  }
137 
138  if(verboseLevel > 1) {
139  G4cout << "G4BGGNucleonElasticXS::GetElementCrossSection for "
140  << dp->GetDefinition()->GetParticleName()
141  << " Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
142  << " in nucleus Z= " << Z << " A= " << theA[Z]
143  << " XS(b)= " << cross/barn
144  << G4endl;
145  }
146  //AR-04Dec2013 if(cross <= fLowestXSection) { cross = 0.0; }
147  return cross;
148 }
149 
150 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
151 
152 G4double
154  G4int Z, G4int A,
155  const G4Isotope*,
156  const G4Element*,
157  const G4Material*)
158 {
159  // this method should be called only for Z = 1
160 
161  G4double cross = 0.0;
162  G4double ekin = dp->GetKineticEnergy();
163 
164  if(ekin <= fSAIDLowEnergyLimit) {
165  cross = theCoulombFac[0]*CoulombFactor(ekin, 1);
166  } else if(ekin <= fSAIDHighEnergyLimit) {
167  cross = fSAID->GetElasticIsotopeCrossSection(particle, ekin, 1, 1);
168  } else if(ekin <= fPDGEnergy) {
169  G4double cross1 =
170  fSAID->GetElasticIsotopeCrossSection(particle, fPDGEnergy, 1, 1);
171  G4double cross2 = theCoulombFac[1];
172  cross = cross1 + (cross2 - cross1)*(ekin - fSAIDHighEnergyLimit)
173  /(fPDGEnergy - fSAIDHighEnergyLimit);
174  } else {
175  fHadron->GetHadronNucleonXscPDG(dp, theProton);
176  cross = fHadron->GetElasticHadronNucleonXsc();
177  }
178  cross *= A;
179 
180  if(verboseLevel > 1) {
181  G4cout << "G4BGGNucleonElasticXS::GetIsoCrossSection for "
182  << dp->GetDefinition()->GetParticleName()
183  << " Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
184  << " in nucleus Z= " << Z << " A= " << A
185  << " XS(b)= " << cross/barn
186  << G4endl;
187  }
188  //AR-04Dec2013 if(cross <= fLowestXSection) { cross = 0.0; }
189  return cross;
190 }
191 
192 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
193 
195 {
196  if(&p == theProton || &p == G4Neutron::Neutron()) {
197  particle = &p;
198 
199  } else {
200  G4cout << "### G4BGGNucleonElasticXS WARNING: is not applicable to "
201  << p.GetParticleName()
202  << G4endl;
203  throw G4HadronicException(__FILE__, __LINE__,
204  "G4BGGNucleonElasticXS::BuildPhysicsTable is used for wrong particle");
205  return;
206  }
207 
208  if(isInitialized) { return; }
209  isInitialized = true;
210 
213 
214  fHadron = new G4HadronNucleonXsc();
215  fSAID = new G4ComponentSAIDTotalXS();
216  fNucleon->BuildPhysicsTable(*particle);
217  fGlauber->BuildPhysicsTable(*particle);
218  if(particle == theProton) {
219  isProton = true;
220  fSAIDHighEnergyLimit = 3*GeV;
221  }
222 
223  G4ParticleDefinition* part = const_cast<G4ParticleDefinition*>(particle);
224  G4ThreeVector mom(0.0,0.0,1.0);
225  G4DynamicParticle dp(part, mom, fGlauberEnergy);
226 
228 
229  G4double csup, csdn;
230  G4int A;
231 
232  if(verboseLevel > 0) {
233  G4cout << "### G4BGGNucleonElasticXS::Initialise for "
234  << particle->GetParticleName() << G4endl;
235  }
236 
237  for(G4int iz=2; iz<93; iz++) {
238 
239  A = G4lrint(nist->GetAtomicMassAmu(iz));
240  theA[iz] = A;
241 
242  csup = fGlauber->GetElasticGlauberGribov(&dp, iz, A);
243  csdn = fNucleon->GetElasticCrossSection(&dp, iz);
244 
245  theGlauberFac[iz] = csdn/csup;
246  if(verboseLevel > 0) {
247  G4cout << "Z= " << iz << " A= " << A
248  << " factor= " << theGlauberFac[iz] << G4endl;
249  }
250  }
251 
252  theCoulombFac[0] =
253  fSAID->GetElasticIsotopeCrossSection(particle,fSAIDLowEnergyLimit,1,1)
254  /CoulombFactor(fSAIDLowEnergyLimit, 1);
255 
256  dp.SetKineticEnergy(fPDGEnergy);
257  fHadron->GetHadronNucleonXscPDG(&dp, theProton);
258  theCoulombFac[1] = fHadron->GetElasticHadronNucleonXsc();
259 
260  if(verboseLevel > 0) {
261  G4cout << "Z=1 A=1 " << " factor0= " << theCoulombFac[0]
262  << " factor1= " << theCoulombFac[1]
263  << G4endl;
264  }
265 
266  dp.SetKineticEnergy(fLowEnergy);
267  for(G4int iz=2; iz<93; iz++) {
268  theCoulombFac[iz] =
269  fNucleon->GetElasticCrossSection(&dp, iz)/CoulombFactor(fLowEnergy, iz);
270  if(verboseLevel > 0) {
271  G4cout << "Z= " << iz << " A= " << theA[iz]
272  << " factor= " << theCoulombFac[iz] << G4endl;
273  }
274  }
275 }
276 
277 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
278 
279 G4double G4BGGNucleonElasticXS::CoulombFactor(G4double kinEnergy, G4int Z)
280 {
281  G4double res= 1.0;
282 
283  // from G4ProtonInelasticCrossSection
284  if(isProton) {
285 
286  if (Z <= 1) { return kinEnergy*kinEnergy; }
287 
288  G4double elog = G4Log(kinEnergy/GeV)/llog10;
289  G4double aa = theA[Z];
290 
291  G4double ff1 = 5.6 - 0.016*aa; // slope of the drop at medium energies.
292  G4double ff2 = 1.37 + 1.37/aa; // start of the slope.
293  G4double ff3 = 0.8 + 18./aa - 0.002*aa; // stephight
294  res = 1.0 + ff3*(1.0 - (1.0/(1+G4Exp(-ff1*(elog + ff2)))));
295 
296  ff1 = 8. - 8./aa - 0.008*aa; // slope of the rise
297  ff2 = 2.34 - 5.4/aa - 0.0028*aa; // start of the rise
298  res /= (1.0 + G4Exp(-ff1*(elog + ff2)));
299 
300  }
301  return res;
302 }
303 
304 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
305 
307 {
308  outFile << "The Barashenkov-Glauber-Gribov cross section handles elastic\n"
309  << "scattering of protons and neutrons from nuclei using the\n"
310  << "Barashenkov parameterization below 91 GeV and the Glauber-Gribov\n"
311  << "parameterization above 91 GeV. n";
312 }
313 
314 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4VCrossSectionDataSet * GetCrossSectionDataSet(const G4String &name, G4bool warning=true)
G4double GetElasticHadronNucleonXsc()
G4double GetKineticEnergy() const
std::ofstream outFile
Definition: GammaRayTel.cc:68
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)
int millibarn
Definition: hepunit.py:40
static G4NistManager * Instance()
const G4String & GetParticleName() const
G4double GetElasticCrossSection(const G4DynamicParticle *aParticle, G4int Z)
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
G4double iz
Definition: TRTMaterials.hh:39
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 G4Log(G4double x)
Definition: G4Log.hh:227
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
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
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
#define G4endl
Definition: G4ios.hh:61
virtual void CrossSectionDescription(std::ostream &) const
double G4double
Definition: G4Types.hh:76
G4double GetElasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)