#include <G4BGGPionInelasticXS.hh>
|
| G4BGGPionInelasticXS (const G4ParticleDefinition *) |
|
virtual | ~G4BGGPionInelasticXS () |
|
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 |
|
| 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 G4Isotope * | SelectIsotope (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 G4String & | GetName () const |
|
Definition at line 68 of file G4BGGPionInelasticXS.hh.
◆ G4BGGPionInelasticXS() [1/2]
Definition at line 56 of file G4BGGPionInelasticXS.cc.
66 for (
G4int i = 0; i < 93; i++) {
static G4Pow * GetInstance()
G4VCrossSectionDataSet(const G4String &nam="")
const G4ParticleDefinition * particle
G4UPiNuclearCrossSection * fPion
G4double theGlauberFac[93]
G4double theCoulombFac[93]
G4HadronNucleonXsc * fHadron
void SetMinKinEnergy(G4double value)
G4ComponentSAIDTotalXS * fSAID
static G4Proton * Proton()
const G4ParticleDefinition * theProton
G4ComponentGGHadronNucleusXsc * fGlauber
void SetMaxKinEnergy(G4double value)
G4double fSAIDHighEnergyLimit
◆ ~G4BGGPionInelasticXS()
G4BGGPionInelasticXS::~G4BGGPionInelasticXS |
( |
| ) |
|
|
virtual |
Definition at line 85 of file G4BGGPionInelasticXS.cc.
G4UPiNuclearCrossSection * fPion
G4HadronNucleonXsc * fHadron
G4ComponentSAIDTotalXS * fSAID
G4ComponentGGHadronNucleusXsc * fGlauber
◆ G4BGGPionInelasticXS() [2/2]
◆ BuildPhysicsTable()
Reimplemented from G4VCrossSectionDataSet.
Definition at line 184 of file G4BGGPionInelasticXS.cc.
189 G4cout <<
"### G4BGGPionInelasticXS WARNING: is not applicable to " 193 "G4BGGPionInelasticXS::BuildPhysicsTable is used for wrong particle");
219 G4cout <<
"### G4BGGPionInelasticXS::Initialise for " 235 G4cout <<
"Z= " << iz <<
" A= " << A
246 dp.SetKineticEnergy(2*
MeV);
247 for(
G4int iz=2; iz<93; iz++) {
256 for(
G4int iz=2; iz<93; iz++) {
virtual G4double GetInelasticIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)
G4double GetInelasticCrossSection(const G4DynamicParticle *aParticle, G4int Z, G4int A)
const G4ParticleDefinition * particle
static G4NistManager * Instance()
G4UPiNuclearCrossSection * fPion
G4double theGlauberFac[93]
G4double theCoulombFac[93]
G4HadronNucleonXsc * fHadron
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
G4double CoulombFactor(G4double kinEnergy, G4int Z)
G4ComponentSAIDTotalXS * fSAID
static G4PionPlus * PionPlus()
const G4ParticleDefinition * theProton
G4ComponentGGHadronNucleusXsc * fGlauber
static G4PionMinus * PionMinus()
G4double fSAIDHighEnergyLimit
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetAtomicMassAmu(const G4String &symb) const
void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetInelasticHadronNucleonXsc()
◆ CoulombFactor()
Definition at line 264 of file G4BGGPionInelasticXS.cc.
268 if(kinEnergy <=
DBL_MIN) {
return res; }
269 else if(A < 2) {
return kinEnergy*kinEnergy; }
277 G4double ff3 = 0.8 + 18/aa - 0.002*aa;
278 res = 1.0 + ff3*(1.0 - (1.0/(1+
fG4pow->
expA(-8*ff1*(elog + 1.37*ff2)))));
280 ff1 = 1. - 1./aa - 0.001*aa;
281 ff2 = 1.17 - 2.7/aa-0.0014*aa;
282 res /= (1 +
fG4pow->
expA(-8.*ff1*(elog + 2*ff2)));
G4double log10A(G4double A) const
double A(double temperature)
G4double expA(G4double A) const
◆ CrossSectionDescription()
void G4BGGPionInelasticXS::CrossSectionDescription |
( |
std::ostream & |
outFile | ) |
const |
|
virtual |
Reimplemented from G4VCrossSectionDataSet.
Definition at line 289 of file G4BGGPionInelasticXS.cc.
291 outFile <<
"The Barashenkov-Glauber-Gribov cross section handles inelastic\n" 292 <<
"pion scattering from nuclei at all energies. The Barashenkov\n" 293 <<
"parameterization is used below 91 GeV and the Glauber-Gribov\n" 294 <<
"parameterization is used above 91 GeV.\n";
◆ GetElementCrossSection()
Reimplemented from G4VCrossSectionDataSet.
Definition at line 115 of file G4BGGPionInelasticXS.cc.
126 if(Z > 92) { Z = 92; }
139 G4cout <<
"G4BGGPionInelasticXS::GetCrossSection for " 142 <<
" in nucleus Z= " << Z <<
" A= " <<
theA[
Z]
143 <<
" XS(b)= " << cross/
barn
G4double GetInelasticCrossSection(const G4DynamicParticle *aParticle, G4int Z, G4int A)
G4UPiNuclearCrossSection * fPion
G4double theGlauberFac[93]
G4double theCoulombFac[93]
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4double GetKineticEnergy() const
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
G4double CoulombFactor(G4double kinEnergy, G4int Z)
G4ComponentGGHadronNucleusXsc * fGlauber
G4ParticleDefinition * GetDefinition() const
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
◆ GetIsoCrossSection()
Reimplemented from G4VCrossSectionDataSet.
Definition at line 152 of file G4BGGPionInelasticXS.cc.
172 G4cout <<
"G4BGGPionInelasticXS::GetCrossSection for " 175 <<
" in nucleus Z= " <<
Z <<
" A= " <<
A 176 <<
" XS(b)= " << cross/
barn virtual G4double GetInelasticIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)
const G4ParticleDefinition * particle
G4double theCoulombFac[93]
G4HadronNucleonXsc * fHadron
G4double GetKineticEnergy() const
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
G4ComponentSAIDTotalXS * fSAID
const G4ParticleDefinition * theProton
G4double fSAIDHighEnergyLimit
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
G4ParticleDefinition * GetDefinition() const
G4double GetInelasticHadronNucleonXsc()
◆ IsElementApplicable()
◆ IsIsoApplicable()
◆ operator=()
◆ fG4pow
G4Pow* G4BGGPionInelasticXS::fG4pow |
|
private |
◆ fGlauber
◆ fGlauberEnergy
G4double G4BGGPionInelasticXS::fGlauberEnergy |
|
private |
◆ fHadron
◆ fLowEnergy
G4double G4BGGPionInelasticXS::fLowEnergy |
|
private |
◆ fPion
◆ fSAID
◆ fSAIDHighEnergyLimit
G4double G4BGGPionInelasticXS::fSAIDHighEnergyLimit |
|
private |
◆ isInitialized
G4bool G4BGGPionInelasticXS::isInitialized |
|
private |
◆ isPiplus
G4bool G4BGGPionInelasticXS::isPiplus |
|
private |
◆ particle
◆ theA
G4int G4BGGPionInelasticXS::theA[93] |
|
private |
◆ theCoulombFac
G4double G4BGGPionInelasticXS::theCoulombFac[93] |
|
private |
◆ theGlauberFac
G4double G4BGGPionInelasticXS::theGlauberFac[93] |
|
private |
◆ theProton
The documentation for this class was generated from the following files: