Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4BGGNucleonElasticXS.cc
이 파일의 문서화 페이지로 가기
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 110543 2018-05-29 13:38:54Z 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;
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;
81  isProton = false;
82  isInitialized = false;
84 }
85 
86 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
87 
89 {
90  delete fSAID;
91  delete fHadron;
92  // The cross section registry will delete fNucleon
93  delete fGlauber;
94 }
95 
96 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
97 
98 G4bool
100  const G4Material*)
101 {
102  return true;
103 }
104 
105 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
106 
108  G4int Z, G4int,
109  const G4Element*,
110  const G4Material*)
111 {
112  return (1 == Z);
113 }
114 
115 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
116 
117 G4double
119  G4int ZZ, const G4Material*)
120 {
121  // this method should be called only for Z > 1
122 
123  G4double cross = 0.0;
124  G4double ekin = dp->GetKineticEnergy();
125  G4int Z = ZZ;
126  if(1 == Z) {
127  cross = 1.0115*GetIsoCrossSection(dp,1,1);
128  } else {
129  if(Z > 92) { Z = 92; }
130 
131  if(ekin <= fLowEnergy) {
132  cross = theCoulombFac[Z]*CoulombFactor(ekin, Z);
133 
134  } else if(ekin > fGlauberEnergy) {
135  cross = theGlauberFac[Z]*fGlauber->GetElasticGlauberGribov(dp, Z, theA[Z]);
136  } else {
137  cross = fNucleon->GetElasticCrossSection(dp, Z);
138  }
139  }
140 
141  if(verboseLevel > 1) {
142  G4cout << "G4BGGNucleonElasticXS::GetElementCrossSection for "
143  << dp->GetDefinition()->GetParticleName()
144  << " Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
145  << " in nucleus Z= " << Z << " A= " << theA[Z]
146  << " XS(b)= " << cross/barn
147  << G4endl;
148  }
149  //AR-04Dec2013 if(cross <= fLowestXSection) { cross = 0.0; }
150  return cross;
151 }
152 
153 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
154 
155 G4double
157  G4int Z, G4int A,
158  const G4Isotope*,
159  const G4Element*,
160  const G4Material*)
161 {
162  // this method should be called only for Z = 1
163 
164  G4double cross = 0.0;
165  G4double ekin = dp->GetKineticEnergy();
166 
167  if(ekin <= fSAIDLowEnergyLimit) {
168  cross = theCoulombFac[0]*CoulombFactor(ekin, 1);
169  } else if(ekin <= fSAIDHighEnergyLimit) {
170  cross = fSAID->GetElasticIsotopeCrossSection(particle, ekin, 1, 1);
171  } else if(ekin <= fPDGEnergy) {
172  G4double cross1 =
174  G4double cross2 = theCoulombFac[1];
175  cross = cross1 + (cross2 - cross1)*(ekin - fSAIDHighEnergyLimit)
177  } else {
180  }
181  cross *= A;
182 
183  if(verboseLevel > 1) {
184  G4cout << "G4BGGNucleonElasticXS::GetIsoCrossSection for "
185  << dp->GetDefinition()->GetParticleName()
186  << " Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
187  << " in nucleus Z= " << Z << " A= " << A
188  << " XS(b)= " << cross/barn
189  << G4endl;
190  }
191  //AR-04Dec2013 if(cross <= fLowestXSection) { cross = 0.0; }
192  return cross;
193 }
194 
195 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
196 
198 {
199  if(&p == theProton || &p == G4Neutron::Neutron()) {
200  particle = &p;
201 
202  } else {
203  G4cout << "### G4BGGNucleonElasticXS WARNING: is not applicable to "
204  << p.GetParticleName()
205  << G4endl;
206  throw G4HadronicException(__FILE__, __LINE__,
207  "G4BGGNucleonElasticXS::BuildPhysicsTable is used for wrong particle");
208  return;
209  }
210 
211  if(isInitialized) { return; }
212  isInitialized = true;
213 
216  fHadron = new G4HadronNucleonXsc();
218 
221 
222  if(particle == theProton) {
223  isProton = true;
225  }
226 
227  G4ThreeVector mom(0.0,0.0,1.0);
229 
231 
232  G4double csup, csdn;
233  G4int A;
234 
235  if(verboseLevel > 0) {
236  G4cout << "### G4BGGNucleonElasticXS::Initialise for "
238  }
239 
240  for(G4int iz=2; iz<93; iz++) {
241 
242  A = G4lrint(nist->GetAtomicMassAmu(iz));
243  theA[iz] = A;
244 
245  csup = fGlauber->GetElasticGlauberGribov(&dp, iz, A);
246  csdn = fNucleon->GetElasticCrossSection(&dp, iz);
247 
248  theGlauberFac[iz] = csdn/csup;
249  if(verboseLevel > 0) {
250  G4cout << "Z= " << iz << " A= " << A
251  << " factor= " << theGlauberFac[iz] << G4endl;
252  }
253  }
254 
255  theCoulombFac[0] =
258 
262 
263  if(verboseLevel > 0) {
264  G4cout << "Z=1 A=1 " << " factor0= " << theCoulombFac[0]
265  << " factor1= " << theCoulombFac[1]
266  << G4endl;
267  }
268 
270  for(G4int iz=2; iz<93; iz++) {
271  theCoulombFac[iz] =
273  if(verboseLevel > 0) {
274  G4cout << "Z= " << iz << " A= " << theA[iz]
275  << " factor= " << theCoulombFac[iz] << G4endl;
276  }
277  }
278 }
279 
280 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
281 
283 {
284  G4double res= 1.0;
285 
286  // from G4ProtonInelasticCrossSection
287  if(isProton) {
288 
289  if (Z <= 1) { return kinEnergy*kinEnergy; }
290 
291  G4double elog = G4Log(kinEnergy/GeV)/llog10;
292  G4double aa = theA[Z];
293 
294  G4double ff1 = 5.6 - 0.016*aa; // slope of the drop at medium energies.
295  G4double ff2 = 1.37 + 1.37/aa; // start of the slope.
296  G4double ff3 = 0.8 + 18./aa - 0.002*aa; // stephight
297  res = 1.0 + ff3*(1.0 - (1.0/(1+G4Exp(-ff1*(elog + ff2)))));
298 
299  ff1 = 8. - 8./aa - 0.008*aa; // slope of the rise
300  ff2 = 2.34 - 5.4/aa - 0.0028*aa; // start of the rise
301  res /= (1.0 + G4Exp(-ff1*(elog + ff2)));
302 
303  }
304  return res;
305 }
306 
307 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
308 
309 void G4BGGNucleonElasticXS::CrossSectionDescription(std::ostream& outFile) const
310 {
311  outFile << "The Barashenkov-Glauber-Gribov cross section handles elastic\n"
312  << "scattering of protons and neutrons from nuclei using the\n"
313  << "Barashenkov parameterization below 91 GeV and the Glauber-Gribov\n"
314  << "parameterization above 91 GeV. n";
315 }
316 
317 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4double CoulombFactor(G4double kinEnergy, G4int Z)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
static constexpr double MeV
Definition: G4SIunits.hh:214
void SetKineticEnergy(G4double aEnergy)
#define G4endl
Definition: G4ios.hh:61
const char * p
Definition: xmltok.h:285
const G4String & GetParticleName() const
G4double GetElasticCrossSection(const G4DynamicParticle *aParticle, G4int Z)
static G4CrossSectionDataSetRegistry * Instance()
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4double GetElasticHadronNucleonXsc()
G4double GetAtomicMassAmu(const G4String &symb) const
G4NucleonNuclearCrossSection * fNucleon
G4double G4Log(G4double x)
Definition: G4Log.hh:230
Float_t Z
G4VCrossSectionDataSet * GetCrossSectionDataSet(const G4String &name, G4bool warning=true)
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
static constexpr double millibarn
Definition: G4SIunits.hh:106
G4ParticleDefinition * GetDefinition() const
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
double A(double temperature)
void SetForAllAtomsAndEnergies(G4bool val)
const G4ParticleDefinition * particle
virtual G4double GetElasticIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)
G4double GetElasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
const G4ParticleDefinition * theProton
G4HadronNucleonXsc * fHadron
G4ComponentSAIDTotalXS * fSAID
G4ComponentGGHadronNucleusXsc * fGlauber
virtual void CrossSectionDescription(std::ostream &) const
int G4lrint(double ad)
Definition: templates.hh:151
int G4int
Definition: G4Types.hh:78
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4BGGNucleonElasticXS(const G4ParticleDefinition *)
static constexpr double barn
Definition: G4SIunits.hh:105
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
virtual G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=0, const G4Material *mat=0)
const G4double llog10
static constexpr double GeV
static constexpr double GeV
Definition: G4SIunits.hh:217
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
static G4NistManager * Instance()