Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BGGNucleonInelasticXS Class Reference

#include <G4BGGNucleonInelasticXS.hh>

+ Inheritance diagram for G4BGGNucleonInelasticXS:

Public Member Functions

 G4BGGNucleonInelasticXS (const G4ParticleDefinition *)
 
 ~G4BGGNucleonInelasticXS () override
 
G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *mat) override
 
G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm, const G4Material *mat) override
 
G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat) override
 
G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr) override
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
void CrossSectionDescription (std::ostream &) const override
 
G4BGGNucleonInelasticXSoperator= (const G4BGGNucleonInelasticXS &right)=delete
 
 G4BGGNucleonInelasticXS (const G4BGGNucleonInelasticXS &)=delete
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
virtual G4double ComputeCrossSectionPerElement (G4double kinEnergy, G4double loge, const G4ParticleDefinition *, const G4Element *, const G4Material *mat=nullptr)
 
virtual G4double ComputeIsoCrossSection (G4double kinEnergy, G4double loge, const G4ParticleDefinition *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
 
virtual const G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy, G4double logE)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual void SetVerboseLevel (G4int value)
 
G4double GetMinKinEnergy () const
 
void SetMinKinEnergy (G4double value)
 
G4double GetMaxKinEnergy () const
 
void SetMaxKinEnergy (G4double value)
 
bool ForAllAtomsAndEnergies () const
 
void SetForAllAtomsAndEnergies (G4bool val)
 
const G4StringGetName () const
 
void SetName (const G4String &nam)
 
G4VCrossSectionDataSetoperator= (const G4VCrossSectionDataSet &right)=delete
 
 G4VCrossSectionDataSet (const G4VCrossSectionDataSet &)=delete
 

Additional Inherited Members

- Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel {0}
 
G4String name
 

Detailed Description

Definition at line 62 of file G4BGGNucleonInelasticXS.hh.

Constructor & Destructor Documentation

◆ G4BGGNucleonInelasticXS() [1/2]

G4BGGNucleonInelasticXS::G4BGGNucleonInelasticXS ( const G4ParticleDefinition * p)
explicit

Definition at line 71 of file G4BGGNucleonInelasticXS.cc.

72 : G4VCrossSectionDataSet("BarashenkovGlauberGribov")
73{
74 verboseLevel = 0;
75 fGlauberEnergy = 91.*CLHEP::GeV;
76 fLowEnergy = 14.*CLHEP::MeV;
77
78 fNucleon = new G4NucleonNuclearCrossSection();
79 fGlauber = new G4ComponentGGHadronNucleusXsc();
80 fHadron = new G4HadronNucleonXsc();
81
82 theProton= G4Proton::Proton();
83 isProton = (theProton == p);
85
86 if (0 == theA[0]) { Initialise(); }
87}
static G4Proton * Proton()
Definition G4Proton.cc:90
G4VCrossSectionDataSet(const G4String &nam="")
void SetForAllAtomsAndEnergies(G4bool val)

Referenced by G4BGGNucleonInelasticXS(), G4ParticleHPBGGNucleonInelasticXS::G4ParticleHPBGGNucleonInelasticXS(), and operator=().

◆ ~G4BGGNucleonInelasticXS()

G4BGGNucleonInelasticXS::~G4BGGNucleonInelasticXS ( )
override

Definition at line 91 of file G4BGGNucleonInelasticXS.cc.

92{
93 delete fHadron;
94}

◆ G4BGGNucleonInelasticXS() [2/2]

G4BGGNucleonInelasticXS::G4BGGNucleonInelasticXS ( const G4BGGNucleonInelasticXS & )
delete

Member Function Documentation

◆ BuildPhysicsTable()

void G4BGGNucleonInelasticXS::BuildPhysicsTable ( const G4ParticleDefinition & p)
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 179 of file G4BGGNucleonInelasticXS.cc.

180{
181 if(&p == theProton || &p == G4Neutron::Neutron()) {
182 isProton = (theProton == &p);
183 } else {
185 ed << "This BGG cross section is applicable only to nucleons and not to "
186 << p.GetParticleName() << G4endl;
187 G4Exception("G4BGGNucleonInelasticXS::BuildPhysicsTable", "had001",
188 FatalException, ed);
189 return;
190 }
191}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
#define G4endl
Definition G4ios.hh:67
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
const G4String & GetParticleName() const

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 290 of file G4BGGNucleonInelasticXS.cc.

291{
292 outFile << "The Barashenkov-Glauber-Gribov cross section calculates inelastic\n"
293 << "scattering of protons and neutrons from nuclei using the\n"
294 << "Barashenkov parameterization below 91 GeV and the Glauber-Gribov\n"
295 << "parameterization above 91 GeV. It uses the G4HadronNucleonXsc\n"
296 << "cross section component for hydrogen targets, and the\n"
297 << "G4ComponentGGHadronNucleusXsc component for other targets.\n";
298}

◆ GetElementCrossSection()

G4double G4BGGNucleonInelasticXS::GetElementCrossSection ( const G4DynamicParticle * dp,
G4int Z,
const G4Material * mat )
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 117 of file G4BGGNucleonInelasticXS.cc.

119{
120 G4double cross = 0.0;
121 G4double ekin = dp->GetKineticEnergy();
122 G4int Z = std::min(ZZ, 92);
123 if (1 == Z) {
124 cross = 1.0115*GetIsoCrossSection(dp,1,1);
125 } else {
126 if (ekin <= fLowEnergy) {
127 cross = CoulombFactor(ekin, Z);
128 cross *= (isProton) ? theCoulombFacP[Z] : theCoulombFacN[Z];
129 } else if (ekin > fGlauberEnergy) {
130 cross = fGlauber->GetInelasticGlauberGribov(dp, Z, theA[Z]);
131 cross *= (isProton) ? theGlauberFacP[Z] : theGlauberFacN[Z];
132 } else {
133 cross = fNucleon->GetElementCrossSection(dp, Z);
134 }
135 }
136
137#ifdef G4VERBOSE
138 if (verboseLevel > 1) {
139 G4cout << "G4BGGNucleonInelasticXS::GetCrossSection for "
141 << " Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
142 << " in nucleus Z= " << Z << " A= " << theA[Z]
143 << " XS(b)= " << cross/barn
144 << G4endl;
145 }
146#endif
147 return cross;
148}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr) override
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const

◆ GetIsoCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 153 of file G4BGGNucleonInelasticXS.cc.

158{
159 // this method should be called only for Z = 1
160 fHadron->HadronNucleonXscNS(dp->GetDefinition(), theProton,
161 dp->GetKineticEnergy());
162 G4double cross = A*fHadron->GetInelasticHadronNucleonXsc();
163
164#ifdef G4VERBOSE
165 if(verboseLevel > 1) {
166 G4cout << "G4BGGNucleonInelasticXS::GetIsoCrossSection for "
168 << " Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
169 << " in nucleus Z=1 A=" << A
170 << " XS(b)= " << cross/barn
171 << G4endl;
172 }
173#endif
174 return cross;
175}
const G4double A[17]

Referenced by GetElementCrossSection().

◆ IsElementApplicable()

G4bool G4BGGNucleonInelasticXS::IsElementApplicable ( const G4DynamicParticle * ,
G4int Z,
const G4Material * mat )
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Reimplemented in G4ParticleHPBGGNucleonInelasticXS.

Definition at line 98 of file G4BGGNucleonInelasticXS.cc.

100{
101 return true;
102}

◆ IsIsoApplicable()

G4bool G4BGGNucleonInelasticXS::IsIsoApplicable ( const G4DynamicParticle * ,
G4int Z,
G4int A,
const G4Element * elm,
const G4Material * mat )
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Reimplemented in G4ParticleHPBGGNucleonInelasticXS.

Definition at line 106 of file G4BGGNucleonInelasticXS.cc.

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

◆ operator=()

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

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