Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4BGGNucleonInelasticXS.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: G4BGGNucleonInelasticXS.cc 110543 2018-05-29 13:38:54Z 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;
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;
86  isProton = false;
87  isInitialized = false;
89 }
90 
91 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
92 
94 {
95  delete fSAID;
96  delete fHadron;
97  // The cross section registry will delete fNucleon
98  delete fGlauber;
99 }
100 
101 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
102 
104  G4int, const G4Material*)
105 {
106  return true;
107 }
108 
109 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
110 
112  G4int Z, G4int,
113  const G4Element*,
114  const G4Material*)
115 {
116  return (1 == Z);
117 }
118 
119 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
120 
121 G4double
123  G4int ZZ, const G4Material*)
124 {
125  // this method should be called only for Z > 1
126 
127  G4double cross = 0.0;
128  G4double ekin = dp->GetKineticEnergy();
129  G4int Z = ZZ;
130  if(1 == Z) {
131  cross = 1.0115*GetIsoCrossSection(dp,1,1);
132  } else if(2 == Z) {
133  if(ekin > fGlauberEnergy) {
135  } else {
136  cross = fNucleon->GetElementCrossSection(dp, Z);
137  }
138 
139  } else {
140  if(Z > 92) { Z = 92; }
141 
142  if(ekin <= fLowEnergy) {
143  cross = theCoulombFac[Z]*CoulombFactor(ekin, Z);
144  } else if(ekin > fGlauberEnergy) {
146  } else {
147  cross = fNucleon->GetElementCrossSection(dp, Z);
148  }
149  }
150 
151  if(verboseLevel > 1) {
152  G4cout << "G4BGGNucleonInelasticXS::GetCrossSection for "
153  << dp->GetDefinition()->GetParticleName()
154  << " Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
155  << " in nucleus Z= " << Z << " A= " << theA[Z]
156  << " XS(b)= " << cross/barn
157  << G4endl;
158  }
159  //AR-18Dec2013 if(cross <= fLowestXSection) { cross = 0.0; }
160  return cross;
161 }
162 
163 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
164 
165 G4double
167  G4int Z, G4int A,
168  const G4Isotope*,
169  const G4Element*,
170  const G4Material*)
171 {
172  // this method should be called only for Z = 1
173 
174  G4double cross = 0.0;
175  G4double ekin = dp->GetKineticEnergy();
176 
177  if(ekin <= fSAIDHighEnergyLimit) {
178  cross = fSAID->GetInelasticIsotopeCrossSection(particle, ekin, 1, 1);
179  } else if(ekin < fHighEnergy) {
181  cross = (theCoulombFac[0]/ekin + 1)*fHadron->GetInelasticHadronNucleonXsc();
182  } else {
184  cross = (theCoulombFac[1]/ekin + 1)*fHadron->GetInelasticHadronNucleonXsc();
185  }
186  cross *= A;
187 
188  if(verboseLevel > 1) {
189  G4cout << "G4BGGNucleonInelasticXS::GetCrossSection for "
190  << dp->GetDefinition()->GetParticleName()
191  << " Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
192  << " in nucleus Z= " << Z << " A= " << A
193  << " XS(b)= " << cross/barn
194  << G4endl;
195  }
196  //AR-18Dec2013 if(cross <= fLowestXSection) { cross = 0.0; }
197  return cross;
198 }
199 
200 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
201 
203 {
204  if(&p == theProton || &p == G4Neutron::Neutron()) {
205  particle = &p;
206  } else {
207  G4cout << "### G4BGGNucleonInelasticXS WARNING: is not applicable to "
208  << p.GetParticleName()
209  << G4endl;
210  throw G4HadronicException(__FILE__, __LINE__,
211  "G4BGGNucleonElasticXS::BuildPhysicsTable is used for wrong particle");
212  return;
213  }
214 
215  if(isInitialized) { return; }
216  isInitialized = true;
217 
220  fHadron = new G4HadronNucleonXsc();
222 
225 
226  if(particle == theProton) {
227  isProton = true;
229  fHighEnergy = 2*GeV;
230  }
231 
232  G4ThreeVector mom(0.0,0.0,1.0);
234 
236  G4int A;
237 
238  G4double csup, csdn;
239 
240  if(verboseLevel > 0) {
241  G4cout << "### G4BGGNucleonInelasticXS::Initialise for "
243  }
244  for(G4int iz=2; iz<93; iz++) {
245 
246  A = G4lrint(nist->GetAtomicMassAmu(iz));
247  theA[iz] = A;
248 
249  csup = fGlauber->GetInelasticGlauberGribov(&dp, iz, A);
250  csdn = fNucleon->GetElementCrossSection(&dp, iz);
251 
252  theGlauberFac[iz] = csdn/csup;
253  if(verboseLevel > 0) {
254  G4cout << "Z= " << iz << " A= " << A
255  << " GlauberFactor= " << theGlauberFac[iz] << G4endl;
256  }
257  }
258  //const G4Material* mat = 0;
259 
265 
266  //G4cout << "Z=1 E(GeV)= " << fSAIDHighEnergyLimit/GeV
267  // << " xsNS(b)= " << fHadron->GetInelasticHadronNucleonXsc()/barn;
269  //G4cout << " xsPDG(b)= " << fHadron->GetInelasticHadronNucleonXsc()/barn;
270  //G4cout << " xsSAID(b)= " << fSAID->GetInelasticIsotopeCrossSection(particle,fSAIDHighEnergyLimit,1,1)/barn << G4endl;
271 
275 
276  //G4cout << "Z=1 E(GeV)= " << fHighEnergy/GeV
277  // << " xsPDG(b)= " << fHadron->GetInelasticHadronNucleonXsc()/barn;
278 
282 
284  //G4cout <<" xsNS(b)= "<<fHadron->GetInelasticHadronNucleonXsc()/barn<<G4endl;
285 
286  if(verboseLevel > 0) {
287  G4cout << "Z=1 A=1" << " CoulombFactor[0]= " << theCoulombFac[0]
288  << " CoulombFactor[1]= " << theCoulombFac[1] << G4endl;
289  }
290  theCoulombFac[2] = 1.0;
291 
293  for(G4int iz=3; iz<93; iz++) {
294  theCoulombFac[iz] =
296 
297  if(verboseLevel > 0) {
298  G4cout << "Z= " << iz << " A= " << theA[iz]
299  << " CoulombFactor= " << theCoulombFac[iz] << G4endl;
300  }
301  }
302 }
303 
304 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
305 
307 {
308  G4double res= 0.0;
309  if(kinEnergy <= 0.0) { return res; }
310  else if (Z <= 1) { return kinEnergy*kinEnergy; }
311 
312  G4double elog = G4Log(kinEnergy/GeV)/llog10;
313  G4double aa = theA[Z];
314 
315  // from G4ProtonInelasticCrossSection
316  if(isProton) {
317 
318  G4double ff1 = 5.6 - 0.016*aa; // slope of the drop at medium energies.
319  G4double ff2 = 1.37 + 1.37/aa; // start of the slope.
320  G4double ff3 = 0.8 + 18./aa - 0.002*aa; // stephight
321  res = 1.0 + ff3*(1.0 - (1.0/(1+G4Exp(-ff1*(elog + ff2)))));
322 
323  ff1 = 8. - 8./aa - 0.008*aa; // slope of the rise
324  ff2 = 2.34 - 5.4/aa - 0.0028*aa; // start of the rise
325  res /= (1.0 + G4Exp(-ff1*(elog + ff2)));
326 
327  } else {
328 
329  // from G4NeutronInelasticCrossSection
330  G4double p3 = 0.6 + 13./aa - 0.0005*aa;
331  G4double p4 = 7.2449 - 0.018242*aa;
332  G4double p5 = 1.36 + 1.8/aa + 0.0005*aa;
333  G4double p6 = 1. + 200./aa + 0.02*aa;
334  G4double p7 = 3.0 - (aa-70.)*(aa-200.)/11000.;
335 
336  G4double firstexp = G4Exp(-p4*(elog + p5));
337  G4double secondexp = G4Exp(-p6*(elog + p7));
338 
339  res = (1.+p3*firstexp/(1. + firstexp))/(1. + secondexp);
340 
341  }
342  return res;
343 }
344 
345 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
346 
347 void G4BGGNucleonInelasticXS::CrossSectionDescription(std::ostream& outFile) const
348 {
349  outFile << "The Barashenkov-Glauber-Gribov cross section calculates inelastic\n"
350  << "scattering of protons and neutrons from nuclei using the\n"
351  << "Barashenkov parameterization below 91 GeV and the Glauber-Gribov\n"
352  << "parameterization above 91 GeV. It uses the G4HadronNucleonXsc\n"
353  << "cross section component for hydrogen targets, and the\n"
354  << "G4ComponentGGHadronNucleusXsc component for other targets.\n";
355 }
356 
357 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Float_t x
Definition: compare.C:6
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
const G4ParticleDefinition * particle
static constexpr double MeV
Definition: G4SIunits.hh:214
void SetKineticEnergy(G4double aEnergy)
G4ComponentSAIDTotalXS * fSAID
#define G4endl
Definition: G4ios.hh:61
G4NucleonNuclearCrossSection * fNucleon
const char * p
Definition: xmltok.h:285
G4double CoulombFactor(G4double kinEnergy, G4int Z)
const G4ParticleDefinition * theProton
G4BGGNucleonInelasticXS(const G4ParticleDefinition *)
const G4String & GetParticleName() const
static G4CrossSectionDataSetRegistry * Instance()
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
G4double GetAtomicMassAmu(const G4String &symb) const
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4double G4Log(G4double x)
Definition: G4Log.hh:230
Float_t Z
G4VCrossSectionDataSet * GetCrossSectionDataSet(const G4String &name, G4bool warning=true)
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)
G4double GetInelasticHadronNucleonXsc()
void SetForAllAtomsAndEnergies(G4bool val)
virtual G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=0, const G4Material *mat=0)
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
virtual G4double GetInelasticIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)
virtual void CrossSectionDescription(std::ostream &) const
int G4lrint(double ad)
Definition: templates.hh:151
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
int G4int
Definition: G4Types.hh:78
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4ComponentGGHadronNucleusXsc * fGlauber
static constexpr double barn
Definition: G4SIunits.hh:105
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
virtual G4double GetElementCrossSection(const G4DynamicParticle *aParticle, G4int Z, const G4Material *mat=0)
const G4double llog10
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4ParticleDefinition *)
static constexpr double GeV
static constexpr double GeV
Definition: G4SIunits.hh:217
static G4NistManager * Instance()