Geant4  10.02.p03
G4BGGNucleonInelasticXS Class Reference

#include <G4BGGNucleonInelasticXS.hh>

Inheritance diagram for G4BGGNucleonInelasticXS:
Collaboration diagram for G4BGGNucleonInelasticXS:

Public Member Functions

 G4BGGNucleonInelasticXS (const G4ParticleDefinition *)
 
virtual ~G4BGGNucleonInelasticXS ()
 
virtual G4bool IsElementApplicable (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)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
 
virtual G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void CrossSectionDescription (std::ostream &) const
 
void SetLowestCrossSection (G4double val)
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
virtual G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual G4int GetVerboseLevel () const
 
virtual void SetVerboseLevel (G4int value)
 
G4double GetMinKinEnergy () const
 
void SetMinKinEnergy (G4double value)
 
G4double GetMaxKinEnergy () const
 
void SetMaxKinEnergy (G4double value)
 
const G4StringGetName () const
 

Private Member Functions

G4double CoulombFactor (G4double kinEnergy, G4int Z)
 
G4BGGNucleonInelasticXSoperator= (const G4BGGNucleonInelasticXS &right)
 
 G4BGGNucleonInelasticXS (const G4BGGNucleonInelasticXS &)
 

Private Attributes

G4double fGlauberEnergy
 
G4double fLowEnergy
 
G4double fHighEnergy
 
G4double fSAIDHighEnergyLimit
 
G4double fLowestXSection
 
G4double theGlauberFac [93]
 
G4double theCoulombFac [93]
 
G4int theA [93]
 
const G4ParticleDefinitionparticle
 
const G4ParticleDefinitiontheProton
 
G4ComponentGGHadronNucleusXscfGlauber
 
G4NucleonNuclearCrossSectionfNucleon
 
G4HadronNucleonXscfHadron
 
G4ComponentSAIDTotalXSfSAID
 
G4bool isProton
 
G4bool isInitialized
 

Additional Inherited Members

- Protected Member Functions inherited from G4VCrossSectionDataSet
void SetName (const G4String &)
 
- Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel
 

Detailed Description

Definition at line 64 of file G4BGGNucleonInelasticXS.hh.

Constructor & Destructor Documentation

◆ G4BGGNucleonInelasticXS() [1/2]

G4BGGNucleonInelasticXS::G4BGGNucleonInelasticXS ( const G4ParticleDefinition p)

Definition at line 65 of file G4BGGNucleonInelasticXS.cc.

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;
88 }
static const double MeV
Definition: G4SIunits.hh:211
G4VCrossSectionDataSet(const G4String &nam="")
G4ComponentSAIDTotalXS * fSAID
int G4int
Definition: G4Types.hh:78
G4ComponentGGHadronNucleusXsc * fGlauber
const G4ParticleDefinition * particle
const G4ParticleDefinition * theProton
static G4Proton * Proton()
Definition: G4Proton.cc:93
static const double GeV
Definition: G4SIunits.hh:214
static const double millibarn
Definition: G4SIunits.hh:105
G4NucleonNuclearCrossSection * fNucleon
Here is the call graph for this function:

◆ ~G4BGGNucleonInelasticXS()

G4BGGNucleonInelasticXS::~G4BGGNucleonInelasticXS ( )
virtual

Definition at line 92 of file G4BGGNucleonInelasticXS.cc.

93 {
94  delete fSAID;
95  delete fHadron;
96  // The cross section registry will delete fNucleon
97  delete fGlauber;
98 }
G4ComponentSAIDTotalXS * fSAID
G4ComponentGGHadronNucleusXsc * fGlauber

◆ G4BGGNucleonInelasticXS() [2/2]

G4BGGNucleonInelasticXS::G4BGGNucleonInelasticXS ( const G4BGGNucleonInelasticXS )
private

Member Function Documentation

◆ BuildPhysicsTable()

void G4BGGNucleonInelasticXS::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 201 of file G4BGGNucleonInelasticXS.cc.

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 
219  fHadron = new G4HadronNucleonXsc();
221 
224 
225  if(particle == theProton) {
226  isProton = true;
228  fHighEnergy = 2*GeV;
229  }
230 
231  G4ThreeVector mom(0.0,0.0,1.0);
233 
235  G4int A;
236 
237  G4double csup, csdn;
238 
239  if(verboseLevel > 0) {
240  G4cout << "### G4BGGNucleonInelasticXS::Initialise for "
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);
264 
265  //G4cout << "Z=1 E(GeV)= " << fSAIDHighEnergyLimit/GeV
266  // << " xsNS(b)= " << fHadron->GetInelasticHadronNucleonXsc()/barn;
268  //G4cout << " xsPDG(b)= " << fHadron->GetInelasticHadronNucleonXsc()/barn;
269  //G4cout << " xsSAID(b)= " << fSAID->GetInelasticIsotopeCrossSection(particle,fSAIDHighEnergyLimit,1,1)/barn << G4endl;
270 
271  dp.SetKineticEnergy(fHighEnergy);
274 
275  //G4cout << "Z=1 E(GeV)= " << fHighEnergy/GeV
276  // << " xsPDG(b)= " << fHadron->GetInelasticHadronNucleonXsc()/barn;
277 
281 
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] =
295 
296  if(verboseLevel > 0) {
297  G4cout << "Z= " << iz << " A= " << theA[iz]
298  << " CoulombFactor= " << theCoulombFac[iz] << G4endl;
299  }
300  }
301 }
G4VCrossSectionDataSet * GetCrossSectionDataSet(const G4String &name, G4bool warning=true)
virtual G4double GetInelasticIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
G4ComponentSAIDTotalXS * fSAID
virtual G4double GetElementCrossSection(const G4DynamicParticle *aParticle, G4int Z, const G4Material *mat=0)
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4ParticleDefinition *)
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
G4ComponentGGHadronNucleusXsc * fGlauber
const G4ParticleDefinition * particle
G4double iz
Definition: TRTMaterials.hh:39
const G4ParticleDefinition * theProton
static G4CrossSectionDataSetRegistry * Instance()
G4double CoulombFactor(G4double kinEnergy, G4int Z)
static const double GeV
Definition: G4SIunits.hh:214
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
int G4lrint(double ad)
Definition: templates.hh:163
G4NucleonNuclearCrossSection * fNucleon
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
#define G4endl
Definition: G4ios.hh:61
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
double G4double
Definition: G4Types.hh:76
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetAtomicMassAmu(const G4String &symb) const
G4double GetInelasticHadronNucleonXsc()
Here is the call graph for this function:

◆ CoulombFactor()

G4double G4BGGNucleonInelasticXS::CoulombFactor ( G4double  kinEnergy,
G4int  Z 
)
private

Definition at line 305 of file G4BGGNucleonInelasticXS.cc.

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 }
Float_t Z
static const double GeV
Definition: G4SIunits.hh:214
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
double G4double
Definition: G4Types.hh:76
const G4double llog10
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CrossSectionDescription()

void G4BGGNucleonInelasticXS::CrossSectionDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 346 of file G4BGGNucleonInelasticXS.cc.

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 }

◆ GetElementCrossSection()

G4double G4BGGNucleonInelasticXS::GetElementCrossSection ( const G4DynamicParticle dp,
G4int  Z,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 121 of file G4BGGNucleonInelasticXS.cc.

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) {
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) {
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 }
virtual G4double GetElementCrossSection(const G4DynamicParticle *aParticle, G4int Z, const G4Material *mat=0)
int G4int
Definition: G4Types.hh:78
static const double GeV
G4double GetKineticEnergy() const
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
G4ComponentGGHadronNucleusXsc * fGlauber
Float_t Z
G4double CoulombFactor(G4double kinEnergy, G4int Z)
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4NucleonNuclearCrossSection * fNucleon
G4ParticleDefinition * GetDefinition() const
#define G4endl
Definition: G4ios.hh:61
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
static const double barn
Definition: G4SIunits.hh:104
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ GetIsoCrossSection()

G4double G4BGGNucleonInelasticXS::GetIsoCrossSection ( const G4DynamicParticle dp,
G4int  Z,
G4int  A,
const G4Isotope iso = 0,
const G4Element elm = 0,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 165 of file G4BGGNucleonInelasticXS.cc.

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) {
180  cross = (theCoulombFac[0]/ekin + 1)*fHadron->GetInelasticHadronNucleonXsc();
181  } else {
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 }
virtual G4double GetInelasticIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)
G4ComponentSAIDTotalXS * fSAID
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4ParticleDefinition *)
static const double GeV
G4double GetKineticEnergy() const
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
Float_t Z
const G4ParticleDefinition * particle
const G4ParticleDefinition * theProton
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
G4ParticleDefinition * GetDefinition() const
#define G4endl
Definition: G4ios.hh:61
static const double barn
Definition: G4SIunits.hh:104
double G4double
Definition: G4Types.hh:76
G4double GetInelasticHadronNucleonXsc()
Here is the call graph for this function:
Here is the caller graph for this function:

◆ IsElementApplicable()

G4bool G4BGGNucleonInelasticXS::IsElementApplicable ( const G4DynamicParticle ,
G4int  Z,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Reimplemented in G4ParticleHPBGGNucleonInelasticXS.

Definition at line 102 of file G4BGGNucleonInelasticXS.cc.

104 {
105  return true;
106 }

◆ IsIsoApplicable()

G4bool G4BGGNucleonInelasticXS::IsIsoApplicable ( const G4DynamicParticle ,
G4int  Z,
G4int  A,
const G4Element elm = 0,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Reimplemented in G4ParticleHPBGGNucleonInelasticXS.

Definition at line 110 of file G4BGGNucleonInelasticXS.cc.

114 {
115  return (1 == Z);
116 }
Float_t Z

◆ operator=()

G4BGGNucleonInelasticXS& G4BGGNucleonInelasticXS::operator= ( const G4BGGNucleonInelasticXS right)
private

◆ SetLowestCrossSection()

void G4BGGNucleonInelasticXS::SetLowestCrossSection ( G4double  val)
inline

Definition at line 125 of file G4BGGNucleonInelasticXS.hh.

126 {
127  fLowestXSection = val;
128 }

Member Data Documentation

◆ fGlauber

G4ComponentGGHadronNucleusXsc* G4BGGNucleonInelasticXS::fGlauber
private

Definition at line 117 of file G4BGGNucleonInelasticXS.hh.

◆ fGlauberEnergy

G4double G4BGGNucleonInelasticXS::fGlauberEnergy
private

Definition at line 105 of file G4BGGNucleonInelasticXS.hh.

◆ fHadron

G4HadronNucleonXsc* G4BGGNucleonInelasticXS::fHadron
private

Definition at line 119 of file G4BGGNucleonInelasticXS.hh.

◆ fHighEnergy

G4double G4BGGNucleonInelasticXS::fHighEnergy
private

Definition at line 107 of file G4BGGNucleonInelasticXS.hh.

◆ fLowEnergy

G4double G4BGGNucleonInelasticXS::fLowEnergy
private

Definition at line 106 of file G4BGGNucleonInelasticXS.hh.

◆ fLowestXSection

G4double G4BGGNucleonInelasticXS::fLowestXSection
private

Definition at line 109 of file G4BGGNucleonInelasticXS.hh.

◆ fNucleon

G4NucleonNuclearCrossSection* G4BGGNucleonInelasticXS::fNucleon
private

Definition at line 118 of file G4BGGNucleonInelasticXS.hh.

◆ fSAID

G4ComponentSAIDTotalXS* G4BGGNucleonInelasticXS::fSAID
private

Definition at line 120 of file G4BGGNucleonInelasticXS.hh.

◆ fSAIDHighEnergyLimit

G4double G4BGGNucleonInelasticXS::fSAIDHighEnergyLimit
private

Definition at line 108 of file G4BGGNucleonInelasticXS.hh.

◆ isInitialized

G4bool G4BGGNucleonInelasticXS::isInitialized
private

Definition at line 122 of file G4BGGNucleonInelasticXS.hh.

◆ isProton

G4bool G4BGGNucleonInelasticXS::isProton
private

Definition at line 121 of file G4BGGNucleonInelasticXS.hh.

◆ particle

const G4ParticleDefinition* G4BGGNucleonInelasticXS::particle
private

Definition at line 114 of file G4BGGNucleonInelasticXS.hh.

◆ theA

G4int G4BGGNucleonInelasticXS::theA[93]
private

Definition at line 112 of file G4BGGNucleonInelasticXS.hh.

◆ theCoulombFac

G4double G4BGGNucleonInelasticXS::theCoulombFac[93]
private

Definition at line 111 of file G4BGGNucleonInelasticXS.hh.

◆ theGlauberFac

G4double G4BGGNucleonInelasticXS::theGlauberFac[93]
private

Definition at line 110 of file G4BGGNucleonInelasticXS.hh.

◆ theProton

const G4ParticleDefinition* G4BGGNucleonInelasticXS::theProton
private

Definition at line 115 of file G4BGGNucleonInelasticXS.hh.


The documentation for this class was generated from the following files: