Geant4-11
Public Member Functions | Protected Attributes | Private Member Functions | Private Attributes | Static Private Attributes
G4BGGPionElasticXS Class Referencefinal

#include <G4BGGPionElasticXS.hh>

Inheritance diagram for G4BGGPionElasticXS:
G4VCrossSectionDataSet

Public Member Functions

void BuildPhysicsTable (const G4ParticleDefinition &) final
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
void CrossSectionDescription (std::ostream &) const final
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
bool ForAllAtomsAndEnergies () const
 
 G4BGGPionElasticXS (const G4ParticleDefinition *)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat) final
 
G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr) final
 
G4double GetMaxKinEnergy () const
 
G4double GetMinKinEnergy () const
 
const G4StringGetName () const
 
virtual G4int GetVerboseLevel () const
 
G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *) final
 
G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm, const G4Material *mat) final
 
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)
 
 ~G4BGGPionElasticXS () final
 

Protected Attributes

G4String name
 
G4int verboseLevel
 

Private Member Functions

G4double CoulombFactorPiPlus (G4double kinEnergy, G4int Z)
 
G4double FactorPiMinus (G4double kinEnergy)
 
 G4BGGPionElasticXS (const G4BGGPionElasticXS &)
 
G4BGGPionElasticXSoperator= (const G4BGGPionElasticXS &right)
 

Private Attributes

G4PowfG4pow
 
G4ComponentGGHadronNucleusXscfGlauber
 
G4double fGlauberEnergy
 
G4HadronNucleonXscfHadron
 
G4double fLowEnergy
 
G4double fLowestEnergy
 
G4UPiNuclearCrossSectionfPion
 
G4bool isForAllAtomsAndEnergies
 
G4bool isMaster
 
G4bool isPiplus
 
G4double maxKinEnergy
 
G4double minKinEnergy
 
G4CrossSectionDataSetRegistryregistry
 
const G4ParticleDefinitionthePiPlus
 
const G4ParticleDefinitiontheProton
 

Static Private Attributes

static G4int theA [93] = {0}
 
static G4double theCoulombFacPiMinus [93] = {0.0}
 
static G4double theCoulombFacPiPlus [93] = {0.0}
 
static G4double theGlauberFacPiMinus [93] = {0.0}
 
static G4double theGlauberFacPiPlus [93] = {0.0}
 

Detailed Description

Definition at line 64 of file G4BGGPionElasticXS.hh.

Constructor & Destructor Documentation

◆ G4BGGPionElasticXS() [1/2]

G4BGGPionElasticXS::G4BGGPionElasticXS ( const G4ParticleDefinition p)
explicit

Definition at line 68 of file G4BGGPionElasticXS.cc.

69 : G4VCrossSectionDataSet("BarashenkovGlauberGribov")
70{
71 verboseLevel = 0;
72 fGlauberEnergy = 91.*GeV;
73 fLowEnergy = 20.*MeV;
74 fLowestEnergy = 1.*MeV;
75 SetMinKinEnergy(0.0);
77
78 fPion = nullptr;
79 fGlauber = nullptr;
80 fHadron = nullptr;
81
83
86 isPiplus = (p == thePiPlus);
87 isMaster = false;
89}
static constexpr double GeV
Definition: G4SIunits.hh:203
static constexpr double MeV
Definition: G4SIunits.hh:200
const G4ParticleDefinition * theProton
const G4ParticleDefinition * thePiPlus
G4UPiNuclearCrossSection * fPion
G4HadronNucleonXsc * fHadron
G4ComponentGGHadronNucleusXsc * fGlauber
static G4HadronicParameters * Instance()
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:97
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
static G4Proton * Proton()
Definition: G4Proton.cc:92
G4VCrossSectionDataSet(const G4String &nam="")
void SetMaxKinEnergy(G4double value)
void SetMinKinEnergy(G4double value)
void SetForAllAtomsAndEnergies(G4bool val)

References fG4pow, fGlauber, fGlauberEnergy, fHadron, fLowEnergy, fLowestEnergy, fPion, G4Pow::GetInstance(), GeV, G4HadronicParameters::Instance(), isMaster, isPiplus, MeV, G4PionPlus::PionPlus(), G4Proton::Proton(), G4VCrossSectionDataSet::SetForAllAtomsAndEnergies(), G4VCrossSectionDataSet::SetMaxKinEnergy(), G4VCrossSectionDataSet::SetMinKinEnergy(), thePiPlus, theProton, and G4VCrossSectionDataSet::verboseLevel.

◆ ~G4BGGPionElasticXS()

G4BGGPionElasticXS::~G4BGGPionElasticXS ( )
final

Definition at line 93 of file G4BGGPionElasticXS.cc.

94{
95 delete fHadron;
96}

References fHadron.

◆ G4BGGPionElasticXS() [2/2]

G4BGGPionElasticXS::G4BGGPionElasticXS ( const G4BGGPionElasticXS )
private

Member Function Documentation

◆ BuildPhysicsTable()

void G4BGGPionElasticXS::BuildPhysicsTable ( const G4ParticleDefinition p)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 175 of file G4BGGPionElasticXS.cc.

176{
177 if(fPion) { return; }
178 if(verboseLevel > 1) {
179 G4cout << "G4BGGPionElasticXS::BuildPhysicsTable for "
180 << p.GetParticleName() << G4endl;
181 }
182 if(&p == G4PionPlus::PionPlus() || &p == G4PionMinus::PionMinus()) {
183 isPiplus = (&p == G4PionPlus::PionPlus());
184 } else {
186 ed << "This BGG cross section is applicable only to pions and not to "
187 << p.GetParticleName() << G4endl;
188 G4Exception("G4BGGPionElasticXS::BuildPhysicsTable", "had001",
189 FatalException, ed);
190 return;
191 }
192
196
198
199 if(0 == theA[0]) {
200#ifdef G4MULTITHREADED
201 G4MUTEXLOCK(&pionElasticXSMutex);
202 if(0 == theA[0]) {
203#endif
204 isMaster = true;
205#ifdef G4MULTITHREADED
206 }
207 G4MUTEXUNLOCK(&pionElasticXSMutex);
208#endif
209 } else {
210 return;
211 }
212
213 if(isMaster && 0 == theA[0]) {
214
215 theA[0] = theA[1] = 1;
216 G4ThreeVector mom(0.0,0.0,1.0);
218
220
221 G4double csup, csdn;
222 for(G4int iz=2; iz<93; ++iz) {
223
224 G4int A = G4lrint(nist->GetAtomicMassAmu(iz));
225 theA[iz] = A;
226
227 csup = fGlauber->GetElasticGlauberGribov(&dp, iz, A);
228 csdn = fPion->GetElasticCrossSection(&dp, iz, A);
229 theGlauberFacPiPlus[iz] = csdn/csup;
230 }
231
232 dp.SetDefinition(G4PionMinus::PionMinus());
233 for(G4int iz=2; iz<93; ++iz) {
234 csup = fGlauber->GetElasticGlauberGribov(&dp, iz, theA[iz]);
235 csdn = fPion->GetElasticCrossSection(&dp, iz, theA[iz]);
236 theGlauberFacPiMinus[iz] = csdn/csup;
237
238 if(verboseLevel > 0) {
239 G4cout << "Z= " << iz << " A= " << theA[iz]
240 << " factorPiPlus= " << theGlauberFacPiPlus[iz]
241 << " factorPiMinus= " << theGlauberFacPiMinus[iz]
242 << G4endl;
243 }
244 }
245 theCoulombFacPiPlus[1] = 1.0;
246 theCoulombFacPiMinus[1]= 1.0;
247 dp.SetKineticEnergy(fLowEnergy);
248 dp.SetDefinition(thePiPlus);
249 for(G4int iz=2; iz<93; ++iz) {
252 }
253 dp.SetDefinition(G4PionMinus::PionMinus());
254 for(G4int iz=2; iz<93; ++iz) {
257
258 if(verboseLevel > 0) {
259 G4cout << "Z= " << iz << " A= " << theA[iz]
260 << " CoulombFactorPiPlus= " << theCoulombFacPiPlus[iz]
261 << " CoulombFactorPiMinus= " << theCoulombFacPiMinus[iz]
262 << G4endl;
263 }
264 }
265 }
266}
@ 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
G4double CoulombFactorPiPlus(G4double kinEnergy, G4int Z)
static G4int theA[93]
static G4double theGlauberFacPiPlus[93]
static G4double theCoulombFacPiMinus[93]
G4double FactorPiMinus(G4double kinEnergy)
static G4double theGlauberFacPiMinus[93]
static G4double theCoulombFacPiPlus[93]
G4double GetElasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
static G4NistManager * Instance()
G4double GetAtomicMassAmu(const G4String &symb) const
const G4String & GetParticleName() const
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:97
G4double GetElasticCrossSection(const G4DynamicParticle *aParticle, G4int Z, G4int A) const
void BuildPhysicsTable(const G4ParticleDefinition &) final
int G4lrint(double ad)
Definition: templates.hh:134

References A, G4UPiNuclearCrossSection::BuildPhysicsTable(), CoulombFactorPiPlus(), FactorPiMinus(), FatalException, fGlauber, fGlauberEnergy, fHadron, fLowEnergy, fPion, G4cout, G4endl, G4Exception(), G4lrint(), G4MUTEXLOCK, G4MUTEXUNLOCK, G4NistManager::GetAtomicMassAmu(), G4UPiNuclearCrossSection::GetElasticCrossSection(), G4ComponentGGHadronNucleusXsc::GetElasticGlauberGribov(), G4ParticleDefinition::GetParticleName(), G4NistManager::Instance(), isMaster, isPiplus, G4PionMinus::PionMinus(), G4PionPlus::PionPlus(), G4DynamicParticle::SetDefinition(), G4DynamicParticle::SetKineticEnergy(), theA, theCoulombFacPiMinus, theCoulombFacPiPlus, theGlauberFacPiMinus, theGlauberFacPiPlus, thePiPlus, and G4VCrossSectionDataSet::verboseLevel.

Referenced by G4QMDReaction::G4QMDReaction().

◆ 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().

◆ CoulombFactorPiPlus()

G4double G4BGGPionElasticXS::CoulombFactorPiPlus ( G4double  kinEnergy,
G4int  Z 
)
private

Definition at line 270 of file G4BGGPionElasticXS.cc.

271{
272 return (kinEnergy > 0.0) ?
273 G4NuclearRadii::CoulombFactor(Z, theA[Z], thePiPlus, kinEnergy) : 0.0;
274}
static G4double CoulombFactor(const G4ParticleDefinition *theParticle, const G4ParticleDefinition *nucleon, G4double ekin)

References G4NuclearRadii::CoulombFactor(), theA, thePiPlus, and Z.

Referenced by BuildPhysicsTable(), and GetElementCrossSection().

◆ CrossSectionDescription()

void G4BGGPionElasticXS::CrossSectionDescription ( std::ostream &  outFile) const
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 286 of file G4BGGPionElasticXS.cc.

287{
288 outFile << "The Barashenkov-Glauber-Gribov cross section handles elastic\n"
289 << "scattering of pions from nuclei at all energies. The\n"
290 << "Barashenkov parameterization is used below 91 GeV and the\n"
291 << "Glauber-Gribov parameterization is used above 91 GeV.\n";
292}

◆ DumpPhysicsTable()

void G4VCrossSectionDataSet::DumpPhysicsTable ( const G4ParticleDefinition )
virtualinherited

◆ FactorPiMinus()

G4double G4BGGPionElasticXS::FactorPiMinus ( G4double  kinEnergy)
private

Definition at line 278 of file G4BGGPionElasticXS.cc.

279{
280 return 1.0/std::sqrt(kinEnergy);
281}

Referenced by BuildPhysicsTable(), and GetElementCrossSection().

◆ 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 G4BGGPionElasticXS::GetElementCrossSection ( const G4DynamicParticle dp,
G4int  Z,
const G4Material mat 
)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 119 of file G4BGGPionElasticXS.cc.

121{
122 // this method should be called only for Z > 1
123 G4double cross = 0.0;
125 G4int Z = std::min(ZZ, 92);
126 if(1 == Z) {
127 cross = 1.0115*GetIsoCrossSection(dp,1,1);
128 } else {
129 if(ekin <= fLowEnergy) {
132 } else if(ekin > fGlauberEnergy) {
134 cross *= fGlauber->GetElasticGlauberGribov(dp, Z, theA[Z]);
135 } else {
136 cross = fPion->GetElasticCrossSection(dp, Z, theA[Z]);
137 }
138 }
139 if(verboseLevel > 1) {
140 G4cout << "G4BGGPionElasticXS::GetElementCrossSection for "
142 << " Ekin(GeV)= " << dp->GetKineticEnergy()
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) final
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments

References barn, CoulombFactorPiPlus(), FactorPiMinus(), fGlauber, fGlauberEnergy, fLowEnergy, fLowestEnergy, fPion, G4cout, G4endl, G4DynamicParticle::GetDefinition(), G4UPiNuclearCrossSection::GetElasticCrossSection(), G4ComponentGGHadronNucleusXsc::GetElasticGlauberGribov(), GetIsoCrossSection(), G4DynamicParticle::GetKineticEnergy(), G4ParticleDefinition::GetParticleName(), isPiplus, G4INCL::Math::max(), G4INCL::Math::min(), theA, theCoulombFacPiMinus, theCoulombFacPiPlus, theGlauberFacPiMinus, theGlauberFacPiPlus, G4VCrossSectionDataSet::verboseLevel, and Z.

Referenced by G4QMDReaction::ApplyYourself().

◆ GetIsoCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 151 of file G4BGGPionElasticXS.cc.

156{
157 // this method should be called only for Z = 1
159 dp->GetKineticEnergy());
161
162 if(verboseLevel > 1) {
163 G4cout << "G4BGGPionElasticXS::GetIsoCrossSection for "
165 << " Ekin(GeV)= " << dp->GetKineticEnergy()
166 << " in nucleus Z= " << Z << " A= " << A
167 << " XS(b)= " << cross/barn
168 << G4endl;
169 }
170 return cross;
171}
G4double GetElasticHadronNucleonXsc() const
G4double HadronNucleonXscNS(const G4ParticleDefinition *theParticle, const G4ParticleDefinition *nucleon, G4double ekin)

References A, barn, fHadron, G4cout, G4endl, G4DynamicParticle::GetDefinition(), G4HadronNucleonXsc::GetElasticHadronNucleonXsc(), G4DynamicParticle::GetKineticEnergy(), G4ParticleDefinition::GetParticleName(), G4HadronNucleonXsc::HadronNucleonXscNS(), 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 G4BGGPionElasticXS::IsElementApplicable ( const G4DynamicParticle ,
G4int  Z,
const G4Material  
)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 101 of file G4BGGPionElasticXS.cc.

103{
104 return true;
105}

◆ IsIsoApplicable()

G4bool G4BGGPionElasticXS::IsIsoApplicable ( const G4DynamicParticle ,
G4int  Z,
G4int  A,
const G4Element elm,
const G4Material mat 
)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 109 of file G4BGGPionElasticXS.cc.

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

References Z.

◆ operator=()

G4BGGPionElasticXS & G4BGGPionElasticXS::operator= ( const G4BGGPionElasticXS 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

◆ fG4pow

G4Pow* G4BGGPionElasticXS::fG4pow
private

Definition at line 113 of file G4BGGPionElasticXS.hh.

Referenced by G4BGGPionElasticXS().

◆ fGlauber

G4ComponentGGHadronNucleusXsc* G4BGGPionElasticXS::fGlauber
private

◆ fGlauberEnergy

G4double G4BGGPionElasticXS::fGlauberEnergy
private

◆ fHadron

G4HadronNucleonXsc* G4BGGPionElasticXS::fHadron
private

◆ fLowEnergy

G4double G4BGGPionElasticXS::fLowEnergy
private

◆ fLowestEnergy

G4double G4BGGPionElasticXS::fLowestEnergy
private

Definition at line 102 of file G4BGGPionElasticXS.hh.

Referenced by G4BGGPionElasticXS(), and GetElementCrossSection().

◆ fPion

G4UPiNuclearCrossSection* G4BGGPionElasticXS::fPion
private

◆ isForAllAtomsAndEnergies

G4bool G4VCrossSectionDataSet::isForAllAtomsAndEnergies
privateinherited

◆ isMaster

G4bool G4BGGPionElasticXS::isMaster
private

Definition at line 119 of file G4BGGPionElasticXS.hh.

Referenced by BuildPhysicsTable(), and G4BGGPionElasticXS().

◆ isPiplus

G4bool G4BGGPionElasticXS::isPiplus
private

◆ maxKinEnergy

G4double G4VCrossSectionDataSet::maxKinEnergy
privateinherited

◆ minKinEnergy

G4double G4VCrossSectionDataSet::minKinEnergy
privateinherited

◆ name

G4String G4VCrossSectionDataSet::name
protectedinherited

◆ registry

G4CrossSectionDataSetRegistry* G4VCrossSectionDataSet::registry
privateinherited

◆ theA

G4int G4BGGPionElasticXS::theA = {0}
staticprivate

◆ theCoulombFacPiMinus

G4double G4BGGPionElasticXS::theCoulombFacPiMinus = {0.0}
staticprivate

Definition at line 107 of file G4BGGPionElasticXS.hh.

Referenced by BuildPhysicsTable(), and GetElementCrossSection().

◆ theCoulombFacPiPlus

G4double G4BGGPionElasticXS::theCoulombFacPiPlus = {0.0}
staticprivate

Definition at line 105 of file G4BGGPionElasticXS.hh.

Referenced by BuildPhysicsTable(), and GetElementCrossSection().

◆ theGlauberFacPiMinus

G4double G4BGGPionElasticXS::theGlauberFacPiMinus = {0.0}
staticprivate

Definition at line 106 of file G4BGGPionElasticXS.hh.

Referenced by BuildPhysicsTable(), and GetElementCrossSection().

◆ theGlauberFacPiPlus

G4double G4BGGPionElasticXS::theGlauberFacPiPlus = {0.0}
staticprivate

Definition at line 104 of file G4BGGPionElasticXS.hh.

Referenced by BuildPhysicsTable(), and GetElementCrossSection().

◆ thePiPlus

const G4ParticleDefinition* G4BGGPionElasticXS::thePiPlus
private

◆ theProton

const G4ParticleDefinition* G4BGGPionElasticXS::theProton
private

Definition at line 110 of file G4BGGPionElasticXS.hh.

Referenced by G4BGGPionElasticXS(), and GetIsoCrossSection().

◆ verboseLevel

G4int G4VCrossSectionDataSet::verboseLevel
protectedinherited

Definition at line 168 of file G4VCrossSectionDataSet.hh.

Referenced by G4BGGNucleonElasticXS::BuildPhysicsTable(), BuildPhysicsTable(), G4BGGPionInelasticXS::BuildPhysicsTable(), G4GammaNuclearXS::BuildPhysicsTable(), G4NeutronCaptureXS::BuildPhysicsTable(), G4NeutronElasticXS::BuildPhysicsTable(), G4NeutronInelasticXS::BuildPhysicsTable(), G4ParticleInelasticXS::BuildPhysicsTable(), G4BGGNucleonInelasticXS::BuildPhysicsTable(), G4LENDCrossSection::create_used_target_map(), G4BGGNucleonElasticXS::G4BGGNucleonElasticXS(), G4BGGNucleonInelasticXS::G4BGGNucleonInelasticXS(), G4BGGPionElasticXS(), G4BGGPionInelasticXS::G4BGGPionInelasticXS(), G4GammaNuclearXS::G4GammaNuclearXS(), G4NeutronCaptureXS::G4NeutronCaptureXS(), G4NeutronElasticXS::G4NeutronElasticXS(), G4NeutronInelasticXS::G4NeutronInelasticXS(), G4ParticleInelasticXS::G4ParticleInelasticXS(), G4NeutronCaptureXS::GetElementCrossSection(), G4NeutronElasticXS::GetElementCrossSection(), G4NeutronInelasticXS::GetElementCrossSection(), G4BGGNucleonElasticXS::GetElementCrossSection(), GetElementCrossSection(), G4BGGPionInelasticXS::GetElementCrossSection(), G4BGGNucleonInelasticXS::GetElementCrossSection(), G4GammaNuclearXS::GetElementCrossSection(), G4ParticleInelasticXS::GetElementCrossSection(), G4BGGNucleonElasticXS::GetIsoCrossSection(), GetIsoCrossSection(), G4BGGPionInelasticXS::GetIsoCrossSection(), G4GammaNuclearXS::GetIsoCrossSection(), G4BGGNucleonInelasticXS::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: