Geant4-11
Public Member Functions | Protected Attributes | Private Member Functions | Private Attributes | Static Private Attributes
G4BGGNucleonInelasticXS Class Reference

#include <G4BGGNucleonInelasticXS.hh>

Inheritance diagram for G4BGGNucleonInelasticXS:
G4VCrossSectionDataSet G4ParticleHPBGGNucleonInelasticXS

Public Member Functions

void BuildPhysicsTable (const G4ParticleDefinition &) override
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
void CrossSectionDescription (std::ostream &) const override
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
bool ForAllAtomsAndEnergies () const
 
 G4BGGNucleonInelasticXS (const G4ParticleDefinition *)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
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
 
G4double GetMaxKinEnergy () const
 
G4double GetMinKinEnergy () const
 
const G4StringGetName () const
 
virtual G4int GetVerboseLevel () const
 
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
 
virtual const G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy, G4double logE)
 
void SetForAllAtomsAndEnergies (G4bool val)
 
void SetMaxKinEnergy (G4double value)
 
void SetMinKinEnergy (G4double value)
 
void SetName (const G4String &nam)
 
virtual void SetVerboseLevel (G4int value)
 
 ~G4BGGNucleonInelasticXS () override
 

Protected Attributes

G4String name
 
G4int verboseLevel
 

Private Member Functions

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

Private Attributes

G4ComponentGGHadronNucleusXscfGlauber
 
G4double fGlauberEnergy
 
G4HadronNucleonXscfHadron
 
G4double fLowEnergy
 
G4NucleonNuclearCrossSectionfNucleon
 
G4bool isForAllAtomsAndEnergies
 
G4bool isMaster
 
G4bool isProton
 
G4double maxKinEnergy
 
G4double minKinEnergy
 
G4CrossSectionDataSetRegistryregistry
 
const G4ParticleDefinitiontheProton
 

Static Private Attributes

static G4int theA [93] = {0}
 
static G4double theCoulombFacN [93] = {0.0}
 
static G4double theCoulombFacP [93] = {0.0}
 
static G4double theGlauberFacN [93] = {0.0}
 
static G4double theGlauberFacP [93] = {0.0}
 

Detailed Description

Definition at line 63 of file G4BGGNucleonInelasticXS.hh.

Constructor & Destructor Documentation

◆ G4BGGNucleonInelasticXS() [1/2]

G4BGGNucleonInelasticXS::G4BGGNucleonInelasticXS ( const G4ParticleDefinition p)
explicit

Definition at line 74 of file G4BGGNucleonInelasticXS.cc.

75 : G4VCrossSectionDataSet("BarashenkovGlauberGribov")
76{
77 verboseLevel = 0;
78 fGlauberEnergy = 91.*GeV;
79 fLowEnergy = 14.*MeV;
80
81 fNucleon = nullptr;
82 fGlauber = nullptr;
83 fHadron = nullptr;
84
86 isProton = (theProton == p);
87 isMaster = false;
89}
static constexpr double GeV
Definition: G4SIunits.hh:203
static constexpr double MeV
Definition: G4SIunits.hh:200
const G4ParticleDefinition * theProton
G4NucleonNuclearCrossSection * fNucleon
G4ComponentGGHadronNucleusXsc * fGlauber
static G4Proton * Proton()
Definition: G4Proton.cc:92
G4VCrossSectionDataSet(const G4String &nam="")
void SetForAllAtomsAndEnergies(G4bool val)

References fGlauber, fGlauberEnergy, fHadron, fLowEnergy, fNucleon, GeV, isMaster, isProton, MeV, G4Proton::Proton(), G4VCrossSectionDataSet::SetForAllAtomsAndEnergies(), theProton, and G4VCrossSectionDataSet::verboseLevel.

◆ ~G4BGGNucleonInelasticXS()

G4BGGNucleonInelasticXS::~G4BGGNucleonInelasticXS ( )
override

Definition at line 93 of file G4BGGNucleonInelasticXS.cc.

94{
95 delete fHadron;
96}

References fHadron.

◆ G4BGGNucleonInelasticXS() [2/2]

G4BGGNucleonInelasticXS::G4BGGNucleonInelasticXS ( const G4BGGNucleonInelasticXS )
private

Member Function Documentation

◆ BuildPhysicsTable()

void G4BGGNucleonInelasticXS::BuildPhysicsTable ( const G4ParticleDefinition p)
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 177 of file G4BGGNucleonInelasticXS.cc.

178{
179 if(nullptr != fNucleon) { return; }
180 if(&p == theProton || &p == G4Neutron::Neutron()) {
181 isProton = (theProton == &p);
182 } else {
184 ed << "This BGG cross section is applicable only to nucleons and not to "
185 << p.GetParticleName() << G4endl;
186 G4Exception("G4BGGNucleonInelasticXS::BuildPhysicsTable", "had001",
187 FatalException, ed);
188 return;
189 }
190
194
196
197 if(0 == theA[0]) {
198#ifdef G4MULTITHREADED
199 G4MUTEXLOCK(&nucleonInelasticXSMutex);
200 if(0 == theA[0]) {
201#endif
202 isMaster = true;
203#ifdef G4MULTITHREADED
204 }
205 G4MUTEXUNLOCK(&nucleonInelasticXSMutex);
206#endif
207 } else {
208 return;
209 }
210
211 if(isMaster && 0 == theA[0]) {
212
213 theA[0] = theA[1] = 1;
214 G4ThreeVector mom(0.0,0.0,1.0);
216
218 G4double csup, csdn;
219
220 if(verboseLevel > 0) {
221 G4cout << "### G4BGGNucleonInelasticXS::Initialise for "
222 << p.GetParticleName() << G4endl;
223 }
224 for(G4int iz=2; iz<93; ++iz) {
225
226 G4int A = G4lrint(nist->GetAtomicMassAmu(iz));
227 theA[iz] = A;
228
229 csup = fGlauber->GetInelasticGlauberGribov(&dp, iz, A);
230 csdn = fNucleon->GetElementCrossSection(&dp, iz);
231 theGlauberFacP[iz] = csdn/csup;
232 }
233
234 dp.SetDefinition(G4Neutron::Neutron());
235 for(G4int iz=2; iz<93; ++iz) {
236 csup = fGlauber->GetInelasticGlauberGribov(&dp, iz, theA[iz]);
237 csdn = fNucleon->GetElementCrossSection(&dp, iz);
238 theGlauberFacN[iz] = csdn/csup;
239
240 if(verboseLevel > 0) {
241 G4cout << "Z= " << iz << " A= " << theA[iz]
242 << " GFactorP= " << theGlauberFacP[iz]
243 << " GFactorN= " << theGlauberFacN[iz] << G4endl;
244 }
245 }
246 theCoulombFacP[1] = theCoulombFacN[1] = 1.0;
247 dp.SetDefinition(theProton);
248 dp.SetKineticEnergy(fLowEnergy);
249 for(G4int iz=2; iz<93; ++iz) {
252 }
253 dp.SetDefinition(G4Neutron::Neutron());
254 for(G4int iz=2; iz<93; ++iz) {
257
258 if(verboseLevel > 0) {
259 G4cout << "Z= " << iz << " A= " << theA[iz]
260 << " CFactorP= " << theCoulombFacP[iz]
261 << " CFactorN= " << theCoulombFacN[iz] << G4endl;
262 }
263 }
264 }
265}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
#define G4MUTEXLOCK(mutex)
Definition: G4Threading.hh:251
#define G4MUTEXUNLOCK(mutex)
Definition: G4Threading.hh:254
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4double theCoulombFacP[93]
static G4double theCoulombFacN[93]
static G4double theGlauberFacN[93]
static G4double theGlauberFacP[93]
G4double CoulombFactor(G4double kinEnergy, G4int Z)
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
static G4NistManager * Instance()
G4double GetAtomicMassAmu(const G4String &symb) const
void BuildPhysicsTable(const G4ParticleDefinition &) final
G4double GetElementCrossSection(const G4DynamicParticle *aParticle, G4int Z, const G4Material *mat=nullptr) final
const G4String & GetParticleName() const
int G4lrint(double ad)
Definition: templates.hh:134

References A, G4NucleonNuclearCrossSection::BuildPhysicsTable(), CoulombFactor(), FatalException, fGlauber, fGlauberEnergy, fHadron, fLowEnergy, fNucleon, G4cout, G4endl, G4Exception(), G4lrint(), G4MUTEXLOCK, G4MUTEXUNLOCK, G4NistManager::GetAtomicMassAmu(), G4NucleonNuclearCrossSection::GetElementCrossSection(), G4ComponentGGHadronNucleusXsc::GetInelasticGlauberGribov(), G4ParticleDefinition::GetParticleName(), G4NistManager::Instance(), isMaster, isProton, G4Neutron::Neutron(), G4DynamicParticle::SetDefinition(), G4DynamicParticle::SetKineticEnergy(), theA, theCoulombFacN, theCoulombFacP, theGlauberFacN, theGlauberFacP, theProton, and G4VCrossSectionDataSet::verboseLevel.

◆ ComputeCrossSection()

G4double G4VCrossSectionDataSet::ComputeCrossSection ( const G4DynamicParticle part,
const G4Element elm,
const G4Material mat = nullptr 
)
inherited

Definition at line 81 of file G4VCrossSectionDataSet.cc.

84{
85 G4int Z = elm->GetZasInt();
86
87 if (IsElementApplicable(part, Z, mat)) {
88 return GetElementCrossSection(part, Z, mat);
89 }
90
91 // isotope-wise cross section making sum over available
92 // isotope cross sections, which may be incomplete, so
93 // the result is corrected
94 size_t nIso = elm->GetNumberOfIsotopes();
95 G4double fact = 0.0;
96 G4double xsec = 0.0;
97
98 // user-defined isotope abundances
99 const G4IsotopeVector* isoVector = elm->GetIsotopeVector();
100 const G4double* abundVector = elm->GetRelativeAbundanceVector();
101
102 for (size_t j=0; j<nIso; ++j) {
103 const G4Isotope* iso = (*isoVector)[j];
104 G4int A = iso->GetN();
105 if(abundVector[j] > 0.0 && IsIsoApplicable(part, Z, A, elm, mat)) {
106 fact += abundVector[j];
107 xsec += abundVector[j]*GetIsoCrossSection(part, Z, A, iso, elm, mat);
108 }
109 }
110 return (fact > 0.0) ? xsec/fact : 0.0;
111}
std::vector< G4Isotope * > G4IsotopeVector
const G4int Z[17]
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:167
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
G4int GetZasInt() const
Definition: G4Element.hh:132
G4IsotopeVector * GetIsotopeVector() const
Definition: G4Element.hh:163
G4int GetN() const
Definition: G4Isotope.hh:93
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
virtual G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)

References A, G4VCrossSectionDataSet::GetElementCrossSection(), G4VCrossSectionDataSet::GetIsoCrossSection(), G4Element::GetIsotopeVector(), G4Isotope::GetN(), G4Element::GetNumberOfIsotopes(), G4Element::GetRelativeAbundanceVector(), G4Element::GetZasInt(), G4VCrossSectionDataSet::IsElementApplicable(), G4VCrossSectionDataSet::IsIsoApplicable(), and Z.

Referenced by G4VCrossSectionDataSet::GetCrossSection().

◆ CoulombFactor()

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

Definition at line 269 of file G4BGGNucleonInelasticXS.cc.

270{
271 G4double res = 0.0;
272
273 if(kinEnergy <= 0.0) { return res; }
274
275 G4double elog = G4Log(kinEnergy/GeV)/llog10;
276 G4double aa = theA[Z];
277
278 if(isProton) {
279
280 res = G4NuclearRadii::CoulombFactor(Z, theA[Z], theProton, kinEnergy);
281
282 // from G4ProtonInelasticCrossSection
283 if(res > 0.0) {
284 G4double ff1 = 5.6 - 0.016*aa; // slope of the drop at medium energies.
285 G4double ff2 = 1.37 + 1.37/aa; // start of the slope.
286 G4double ff3 = 0.8 + 18./aa - 0.002*aa; // stephight
287 res *= (1.0 + ff3*(1.0 - (1.0/(1+G4Exp(-ff1*(elog + ff2))))));
288 ff1 = 8. - 8./aa - 0.008*aa; // slope of the rise
289 ff2 = 2.34 - 5.4/aa - 0.0028*aa; // start of the rise
290 res /= (1.0 + G4Exp(-ff1*(elog + ff2)));
291 }
292 } else {
293 // from G4NeutronInelasticCrossSection
294 G4double p3 = 0.6 + 13./aa - 0.0005*aa;
295 G4double p4 = 7.2449 - 0.018242*aa;
296 G4double p5 = 1.36 + 1.8/aa + 0.0005*aa;
297 G4double p6 = 1. + 200./aa + 0.02*aa;
298 G4double p7 = 3.0 - (aa-70.)*(aa-200.)/11000.;
299
300 G4double firstexp = G4Exp(-p4*(elog + p5));
301 G4double secondexp = G4Exp(-p6*(elog + p7));
302
303 res = (1.+p3*firstexp/(1. + firstexp))/(1. + secondexp);
304 }
305 return res;
306}
const G4double llog10
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double G4Log(G4double x)
Definition: G4Log.hh:226
static G4double CoulombFactor(const G4ParticleDefinition *theParticle, const G4ParticleDefinition *nucleon, G4double ekin)

References G4NuclearRadii::CoulombFactor(), G4Exp(), G4Log(), GeV, isProton, llog10, theA, theProton, and Z.

Referenced by BuildPhysicsTable(), and GetElementCrossSection().

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 310 of file G4BGGNucleonInelasticXS.cc.

311{
312 outFile << "The Barashenkov-Glauber-Gribov cross section calculates inelastic\n"
313 << "scattering of protons and neutrons from nuclei using the\n"
314 << "Barashenkov parameterization below 91 GeV and the Glauber-Gribov\n"
315 << "parameterization above 91 GeV. It uses the G4HadronNucleonXsc\n"
316 << "cross section component for hydrogen targets, and the\n"
317 << "G4ComponentGGHadronNucleusXsc component for other targets.\n";
318}

◆ DumpPhysicsTable()

void G4VCrossSectionDataSet::DumpPhysicsTable ( const G4ParticleDefinition )
virtualinherited

◆ ForAllAtomsAndEnergies()

bool G4VCrossSectionDataSet::ForAllAtomsAndEnergies ( ) const
inlineinherited

◆ GetCrossSection()

G4double G4VCrossSectionDataSet::GetCrossSection ( const G4DynamicParticle dp,
const G4Element elm,
const G4Material mat = nullptr 
)
inlineinherited

Definition at line 187 of file G4VCrossSectionDataSet.hh.

190{
191 return ComputeCrossSection(dp, elm, mat);
192}
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)

References G4VCrossSectionDataSet::ComputeCrossSection().

◆ GetElementCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 119 of file G4BGGNucleonInelasticXS.cc.

121{
122 G4double cross = 0.0;
123 G4double ekin = dp->GetKineticEnergy();
124 G4int Z = std::min(ZZ, 92);
125 if(1 == Z) {
126 cross = 1.0115*GetIsoCrossSection(dp,1,1);
127 } else {
128 if(ekin <= fLowEnergy) {
129 cross = (isProton) ? theCoulombFacP[Z] : theCoulombFacN[Z];
130 cross *= CoulombFactor(ekin, Z);
131 } else if(ekin > fGlauberEnergy) {
132 cross = (isProton) ? theGlauberFacP[Z] : theGlauberFacN[Z];
133 cross *= fGlauber->GetInelasticGlauberGribov(dp, Z, theA[Z]);
134 } else {
135 cross = fNucleon->GetElementCrossSection(dp, Z);
136 }
137 }
138
139 if(verboseLevel > 1) {
140 G4cout << "G4BGGNucleonInelasticXS::GetCrossSection for "
142 << " Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
143 << " in nucleus Z= " << Z << " A= " << theA[Z]
144 << " XS(b)= " << cross/barn
145 << G4endl;
146 }
147 return cross;
148}
static constexpr double barn
Definition: G4SIunits.hh:85
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
static constexpr double GeV
T min(const T t1, const T t2)
brief Return the smallest of the two arguments

References barn, CoulombFactor(), fGlauber, fGlauberEnergy, fLowEnergy, fNucleon, G4cout, G4endl, G4DynamicParticle::GetDefinition(), G4NucleonNuclearCrossSection::GetElementCrossSection(), G4ComponentGGHadronNucleusXsc::GetInelasticGlauberGribov(), GetIsoCrossSection(), G4DynamicParticle::GetKineticEnergy(), G4ParticleDefinition::GetParticleName(), CLHEP::GeV, isProton, G4INCL::Math::min(), theA, theCoulombFacN, theCoulombFacP, theGlauberFacN, theGlauberFacP, G4VCrossSectionDataSet::verboseLevel, and Z.

◆ 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
161 dp->GetKineticEnergy());
163
164 if(verboseLevel > 1) {
165 G4cout << "G4BGGNucleonInelasticXS::GetIsoCrossSection for "
167 << " Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
168 << " in nucleus Z= " << Z << " A= " << theA[Z]
169 << " XS(b)= " << cross/barn
170 << G4endl;
171 }
172 return cross;
173}
G4double GetInelasticHadronNucleonXsc() const
G4double HadronNucleonXscNS(const G4ParticleDefinition *theParticle, const G4ParticleDefinition *nucleon, G4double ekin)

References A, barn, fHadron, G4cout, G4endl, G4DynamicParticle::GetDefinition(), G4HadronNucleonXsc::GetInelasticHadronNucleonXsc(), G4DynamicParticle::GetKineticEnergy(), G4ParticleDefinition::GetParticleName(), CLHEP::GeV, G4HadronNucleonXsc::HadronNucleonXscNS(), theA, theProton, G4VCrossSectionDataSet::verboseLevel, and Z.

Referenced by GetElementCrossSection().

◆ GetMaxKinEnergy()

G4double G4VCrossSectionDataSet::GetMaxKinEnergy ( ) const
inlineinherited

◆ GetMinKinEnergy()

G4double G4VCrossSectionDataSet::GetMinKinEnergy ( ) const
inlineinherited

◆ GetName()

const G4String & G4VCrossSectionDataSet::GetName ( ) const
inlineinherited

◆ GetVerboseLevel()

G4int G4VCrossSectionDataSet::GetVerboseLevel ( ) const
inlinevirtualinherited

◆ IsElementApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Reimplemented in G4ParticleHPBGGNucleonInelasticXS.

Definition at line 100 of file G4BGGNucleonInelasticXS.cc.

102{
103 return true;
104}

◆ 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 108 of file G4BGGNucleonInelasticXS.cc.

112{
113 return (1 == Z);
114}

References Z.

◆ operator=()

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

◆ SelectIsotope()

const G4Isotope * G4VCrossSectionDataSet::SelectIsotope ( const G4Element anElement,
G4double  kinEnergy,
G4double  logE 
)
virtualinherited

Reimplemented in G4GammaNuclearXS, G4NeutronCaptureXS, G4NeutronElasticXS, G4NeutronInelasticXS, and G4ParticleInelasticXS.

Definition at line 149 of file G4VCrossSectionDataSet.cc.

151{
152 size_t nIso = anElement->GetNumberOfIsotopes();
153 const G4Isotope* iso = anElement->GetIsotope(0);
154
155 // more than 1 isotope
156 if(1 < nIso) {
157 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
158 G4double sum = 0.0;
160 for (size_t j=0; j<nIso; ++j) {
161 sum += abundVector[j];
162 if(q <= sum) {
163 iso = anElement->GetIsotope(j);
164 break;
165 }
166 }
167 }
168 return iso;
169}
#define G4UniformRand()
Definition: Randomize.hh:52
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170

References G4UniformRand, G4Element::GetIsotope(), G4Element::GetNumberOfIsotopes(), and G4Element::GetRelativeAbundanceVector().

◆ SetForAllAtomsAndEnergies()

void G4VCrossSectionDataSet::SetForAllAtomsAndEnergies ( G4bool  val)
inlineinherited

◆ SetMaxKinEnergy()

void G4VCrossSectionDataSet::SetMaxKinEnergy ( G4double  value)
inlineinherited

◆ SetMinKinEnergy()

void G4VCrossSectionDataSet::SetMinKinEnergy ( G4double  value)
inlineinherited

◆ SetName()

void G4VCrossSectionDataSet::SetName ( const G4String nam)
inlineinherited

Definition at line 240 of file G4VCrossSectionDataSet.hh.

241{
242 name = nam;
243}

References G4VCrossSectionDataSet::name.

Referenced by G4ParticleHPInelasticData::G4ParticleHPInelasticData().

◆ SetVerboseLevel()

void G4VCrossSectionDataSet::SetVerboseLevel ( G4int  value)
inlinevirtualinherited

Field Documentation

◆ fGlauber

G4ComponentGGHadronNucleusXsc* G4BGGNucleonInelasticXS::fGlauber
private

◆ fGlauberEnergy

G4double G4BGGNucleonInelasticXS::fGlauberEnergy
private

◆ fHadron

G4HadronNucleonXsc* G4BGGNucleonInelasticXS::fHadron
private

◆ fLowEnergy

G4double G4BGGNucleonInelasticXS::fLowEnergy
private

◆ fNucleon

G4NucleonNuclearCrossSection* G4BGGNucleonInelasticXS::fNucleon
private

◆ isForAllAtomsAndEnergies

G4bool G4VCrossSectionDataSet::isForAllAtomsAndEnergies
privateinherited

◆ isMaster

G4bool G4BGGNucleonInelasticXS::isMaster
private

Definition at line 112 of file G4BGGNucleonInelasticXS.hh.

Referenced by BuildPhysicsTable(), and G4BGGNucleonInelasticXS().

◆ isProton

G4bool G4BGGNucleonInelasticXS::isProton
private

◆ maxKinEnergy

G4double G4VCrossSectionDataSet::maxKinEnergy
privateinherited

◆ minKinEnergy

G4double G4VCrossSectionDataSet::minKinEnergy
privateinherited

◆ name

G4String G4VCrossSectionDataSet::name
protectedinherited

◆ registry

G4CrossSectionDataSetRegistry* G4VCrossSectionDataSet::registry
privateinherited

◆ theA

G4int G4BGGNucleonInelasticXS::theA = {0}
staticprivate

◆ theCoulombFacN

G4double G4BGGNucleonInelasticXS::theCoulombFacN = {0.0}
staticprivate

Definition at line 103 of file G4BGGNucleonInelasticXS.hh.

Referenced by BuildPhysicsTable(), and GetElementCrossSection().

◆ theCoulombFacP

G4double G4BGGNucleonInelasticXS::theCoulombFacP = {0.0}
staticprivate

Definition at line 101 of file G4BGGNucleonInelasticXS.hh.

Referenced by BuildPhysicsTable(), and GetElementCrossSection().

◆ theGlauberFacN

G4double G4BGGNucleonInelasticXS::theGlauberFacN = {0.0}
staticprivate

Definition at line 102 of file G4BGGNucleonInelasticXS.hh.

Referenced by BuildPhysicsTable(), and GetElementCrossSection().

◆ theGlauberFacP

G4double G4BGGNucleonInelasticXS::theGlauberFacP = {0.0}
staticprivate

Definition at line 100 of file G4BGGNucleonInelasticXS.hh.

Referenced by BuildPhysicsTable(), and GetElementCrossSection().

◆ theProton

const G4ParticleDefinition* G4BGGNucleonInelasticXS::theProton
private

◆ verboseLevel

G4int G4VCrossSectionDataSet::verboseLevel
protectedinherited

Definition at line 168 of file G4VCrossSectionDataSet.hh.

Referenced by G4BGGNucleonElasticXS::BuildPhysicsTable(), G4BGGPionElasticXS::BuildPhysicsTable(), G4BGGPionInelasticXS::BuildPhysicsTable(), G4GammaNuclearXS::BuildPhysicsTable(), G4NeutronCaptureXS::BuildPhysicsTable(), G4NeutronElasticXS::BuildPhysicsTable(), G4NeutronInelasticXS::BuildPhysicsTable(), G4ParticleInelasticXS::BuildPhysicsTable(), BuildPhysicsTable(), G4LENDCrossSection::create_used_target_map(), G4BGGNucleonElasticXS::G4BGGNucleonElasticXS(), G4BGGNucleonInelasticXS(), G4BGGPionElasticXS::G4BGGPionElasticXS(), G4BGGPionInelasticXS::G4BGGPionInelasticXS(), G4GammaNuclearXS::G4GammaNuclearXS(), G4NeutronCaptureXS::G4NeutronCaptureXS(), G4NeutronElasticXS::G4NeutronElasticXS(), G4NeutronInelasticXS::G4NeutronInelasticXS(), G4ParticleInelasticXS::G4ParticleInelasticXS(), G4NeutronCaptureXS::GetElementCrossSection(), G4NeutronElasticXS::GetElementCrossSection(), G4NeutronInelasticXS::GetElementCrossSection(), G4BGGNucleonElasticXS::GetElementCrossSection(), G4BGGPionElasticXS::GetElementCrossSection(), G4BGGPionInelasticXS::GetElementCrossSection(), GetElementCrossSection(), G4GammaNuclearXS::GetElementCrossSection(), G4ParticleInelasticXS::GetElementCrossSection(), G4BGGNucleonElasticXS::GetIsoCrossSection(), G4BGGPionElasticXS::GetIsoCrossSection(), G4BGGPionInelasticXS::GetIsoCrossSection(), G4GammaNuclearXS::GetIsoCrossSection(), GetIsoCrossSection(), G4VCrossSectionDataSet::GetVerboseLevel(), G4NeutronElasticXS::Initialise(), G4ParticleInelasticXS::IsoCrossSection(), G4NeutronCaptureXS::IsoCrossSection(), G4NeutronInelasticXS::IsoCrossSection(), G4GammaNuclearXS::RetrieveVector(), G4NeutronCaptureXS::RetrieveVector(), G4NeutronInelasticXS::RetrieveVector(), G4ParticleInelasticXS::RetrieveVector(), and G4VCrossSectionDataSet::SetVerboseLevel().


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