Geant4  10.02.p03
G4BGGNucleonElasticXS Class Reference

#include <G4BGGNucleonElasticXS.hh>

Inheritance diagram for G4BGGNucleonElasticXS:
Collaboration diagram for G4BGGNucleonElasticXS:

Public Member Functions

 G4BGGNucleonElasticXS (const G4ParticleDefinition *)
 
virtual ~G4BGGNucleonElasticXS ()
 
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)
 
G4BGGNucleonElasticXSoperator= (const G4BGGNucleonElasticXS &right)
 
 G4BGGNucleonElasticXS (const G4BGGNucleonElasticXS &)
 

Private Attributes

G4double fGlauberEnergy
 
G4double fPDGEnergy
 
G4double fLowEnergy
 
G4double fSAIDLowEnergyLimit
 
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 65 of file G4BGGNucleonElasticXS.hh.

Constructor & Destructor Documentation

◆ G4BGGNucleonElasticXS() [1/2]

G4BGGNucleonElasticXS::G4BGGNucleonElasticXS ( const G4ParticleDefinition p)

Definition at line 60 of file G4BGGNucleonElasticXS.cc.

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

◆ ~G4BGGNucleonElasticXS()

G4BGGNucleonElasticXS::~G4BGGNucleonElasticXS ( )
virtual

Definition at line 87 of file G4BGGNucleonElasticXS.cc.

88 {
89  delete fSAID;
90  delete fHadron;
91  // The cross section registry will delete fNucleon
92  delete fGlauber;
93 }
G4HadronNucleonXsc * fHadron
G4ComponentGGHadronNucleusXsc * fGlauber
G4ComponentSAIDTotalXS * fSAID

◆ G4BGGNucleonElasticXS() [2/2]

G4BGGNucleonElasticXS::G4BGGNucleonElasticXS ( const G4BGGNucleonElasticXS )
private

Member Function Documentation

◆ BuildPhysicsTable()

void G4BGGNucleonElasticXS::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 196 of file G4BGGNucleonElasticXS.cc.

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 
215  fHadron = new G4HadronNucleonXsc();
217 
220 
221  if(particle == theProton) {
222  isProton = true;
224  }
225 
226  G4ThreeVector mom(0.0,0.0,1.0);
228 
230 
231  G4double csup, csdn;
232  G4int A;
233 
234  if(verboseLevel > 0) {
235  G4cout << "### G4BGGNucleonElasticXS::Initialise for "
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] =
257 
258  dp.SetKineticEnergy(fPDGEnergy);
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] =
272  if(verboseLevel > 0) {
273  G4cout << "Z= " << iz << " A= " << theA[iz]
274  << " factor= " << theCoulombFac[iz] << G4endl;
275  }
276  }
277 }
G4VCrossSectionDataSet * GetCrossSectionDataSet(const G4String &name, G4bool warning=true)
G4double GetElasticHadronNucleonXsc()
G4HadronNucleonXsc * fHadron
G4NucleonNuclearCrossSection * fNucleon
G4ComponentGGHadronNucleusXsc * fGlauber
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
virtual G4double GetElasticIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
G4ComponentSAIDTotalXS * fSAID
const G4ParticleDefinition * theProton
G4double GetElasticCrossSection(const G4DynamicParticle *aParticle, G4int Z)
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
G4double iz
Definition: TRTMaterials.hh:39
static G4CrossSectionDataSetRegistry * Instance()
static const double GeV
Definition: G4SIunits.hh:214
G4double CoulombFactor(G4double kinEnergy, G4int Z)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
const G4ParticleDefinition * particle
G4double GetElasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
int G4lrint(double ad)
Definition: templates.hh:163
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetAtomicMassAmu(const G4String &symb) const
Here is the call graph for this function:

◆ CoulombFactor()

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

Definition at line 281 of file G4BGGNucleonElasticXS.cc.

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 }
const G4double llog10
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
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 308 of file G4BGGNucleonElasticXS.cc.

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 }

◆ GetElementCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 117 of file G4BGGNucleonElasticXS.cc.

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 }
G4NucleonNuclearCrossSection * fNucleon
G4ComponentGGHadronNucleusXsc * fGlauber
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 const double GeV
G4double GetElasticCrossSection(const G4DynamicParticle *aParticle, G4int Z)
G4double GetKineticEnergy() const
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
Float_t Z
G4double CoulombFactor(G4double kinEnergy, G4int Z)
G4double GetElasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
G4ParticleDefinition * GetDefinition() const
#define G4endl
Definition: G4ios.hh:61
static const double barn
Definition: G4SIunits.hh:104
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ GetIsoCrossSection()

G4double G4BGGNucleonElasticXS::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 155 of file G4BGGNucleonElasticXS.cc.

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 =
173  G4double cross2 = theCoulombFac[1];
174  cross = cross1 + (cross2 - cross1)*(ekin - fSAIDHighEnergyLimit)
176  } else {
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 }
G4double GetElasticHadronNucleonXsc()
G4HadronNucleonXsc * fHadron
virtual G4double GetElasticIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)
G4ComponentSAIDTotalXS * fSAID
static const double GeV
const G4ParticleDefinition * theProton
G4double GetKineticEnergy() const
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
Float_t Z
G4double CoulombFactor(G4double kinEnergy, G4int Z)
const G4ParticleDefinition * particle
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
Here is the call graph for this function:
Here is the caller graph for this function:

◆ IsElementApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 98 of file G4BGGNucleonElasticXS.cc.

100 {
101  return true;
102 }

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 106 of file G4BGGNucleonElasticXS.cc.

110 {
111  return (1 == Z);
112 }
Float_t Z

◆ operator=()

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

◆ SetLowestCrossSection()

void G4BGGNucleonElasticXS::SetLowestCrossSection ( G4double  val)
inline

Definition at line 126 of file G4BGGNucleonElasticXS.hh.

127 {
128  fLowestXSection = val;
129 }

Member Data Documentation

◆ fGlauber

G4ComponentGGHadronNucleusXsc* G4BGGNucleonElasticXS::fGlauber
private

Definition at line 118 of file G4BGGNucleonElasticXS.hh.

◆ fGlauberEnergy

G4double G4BGGNucleonElasticXS::fGlauberEnergy
private

Definition at line 106 of file G4BGGNucleonElasticXS.hh.

◆ fHadron

G4HadronNucleonXsc* G4BGGNucleonElasticXS::fHadron
private

Definition at line 120 of file G4BGGNucleonElasticXS.hh.

◆ fLowEnergy

G4double G4BGGNucleonElasticXS::fLowEnergy
private

Definition at line 108 of file G4BGGNucleonElasticXS.hh.

◆ fLowestXSection

G4double G4BGGNucleonElasticXS::fLowestXSection
private

Definition at line 111 of file G4BGGNucleonElasticXS.hh.

◆ fNucleon

G4NucleonNuclearCrossSection* G4BGGNucleonElasticXS::fNucleon
private

Definition at line 119 of file G4BGGNucleonElasticXS.hh.

◆ fPDGEnergy

G4double G4BGGNucleonElasticXS::fPDGEnergy
private

Definition at line 107 of file G4BGGNucleonElasticXS.hh.

◆ fSAID

G4ComponentSAIDTotalXS* G4BGGNucleonElasticXS::fSAID
private

Definition at line 121 of file G4BGGNucleonElasticXS.hh.

◆ fSAIDHighEnergyLimit

G4double G4BGGNucleonElasticXS::fSAIDHighEnergyLimit
private

Definition at line 110 of file G4BGGNucleonElasticXS.hh.

◆ fSAIDLowEnergyLimit

G4double G4BGGNucleonElasticXS::fSAIDLowEnergyLimit
private

Definition at line 109 of file G4BGGNucleonElasticXS.hh.

◆ isInitialized

G4bool G4BGGNucleonElasticXS::isInitialized
private

Definition at line 123 of file G4BGGNucleonElasticXS.hh.

◆ isProton

G4bool G4BGGNucleonElasticXS::isProton
private

Definition at line 122 of file G4BGGNucleonElasticXS.hh.

◆ particle

const G4ParticleDefinition* G4BGGNucleonElasticXS::particle
private

Definition at line 116 of file G4BGGNucleonElasticXS.hh.

◆ theA

G4int G4BGGNucleonElasticXS::theA[93]
private

Definition at line 114 of file G4BGGNucleonElasticXS.hh.

◆ theCoulombFac

G4double G4BGGNucleonElasticXS::theCoulombFac[93]
private

Definition at line 113 of file G4BGGNucleonElasticXS.hh.

◆ theGlauberFac

G4double G4BGGNucleonElasticXS::theGlauberFac[93]
private

Definition at line 112 of file G4BGGNucleonElasticXS.hh.

◆ theProton

const G4ParticleDefinition* G4BGGNucleonElasticXS::theProton
private

Definition at line 117 of file G4BGGNucleonElasticXS.hh.


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