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

#include <G4EMDissociationCrossSection.hh>

Inheritance diagram for G4EMDissociationCrossSection:
G4VCrossSectionDataSet

Public Member Functions

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

Protected Attributes

G4String name
 
G4int verboseLevel
 

Private Attributes

G4double epsilon
 
G4bool isForAllAtomsAndEnergies
 
G4double J
 
G4double maxKinEnergy
 
G4double minKinEnergy
 
G4double Qprime
 
G4double r0
 
G4CrossSectionDataSetRegistryregistry
 
G4EMDissociationSpectrumthePhotonSpectrum
 
G4double xd
 

Detailed Description

Definition at line 79 of file G4EMDissociationCrossSection.hh.

Constructor & Destructor Documentation

◆ G4EMDissociationCrossSection()

G4EMDissociationCrossSection::G4EMDissociationCrossSection ( )

Definition at line 79 of file G4EMDissociationCrossSection.cc.

80 : G4VCrossSectionDataSet("Electromagnetic dissociation")
81{
82 // This function makes use of the class which can sample the virtual photon
83 // spectrum, G4EMDissociationSpectrum.
84
86
87 // Define other constants.
88
89 r0 = 1.18 * fermi;
90 J = 36.8 * MeV;
91 Qprime = 17.0 * MeV;
92 epsilon = 0.0768;
93 xd = 0.25;
94}
static constexpr double fermi
Definition: G4SIunits.hh:83
static constexpr double MeV
Definition: G4SIunits.hh:200
G4EMDissociationSpectrum * thePhotonSpectrum
G4VCrossSectionDataSet(const G4String &nam="")

References epsilon, fermi, J, MeV, Qprime, r0, thePhotonSpectrum, and xd.

◆ ~G4EMDissociationCrossSection()

G4EMDissociationCrossSection::~G4EMDissociationCrossSection ( )

Definition at line 98 of file G4EMDissociationCrossSection.cc.

99{
100 delete thePhotonSpectrum;
101}

References thePhotonSpectrum.

Member Function Documentation

◆ BuildPhysicsTable()

void G4VCrossSectionDataSet::BuildPhysicsTable ( const G4ParticleDefinition )
virtualinherited

◆ 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
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[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().

◆ CrossSectionDescription()

void G4VCrossSectionDataSet::CrossSectionDescription ( std::ostream &  outFile) const
virtualinherited

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

◆ GetCrossSectionForProjectile()

G4PhysicsFreeVector * G4EMDissociationCrossSection::GetCrossSectionForProjectile ( G4double  AP,
G4double  ZP,
G4double  ,
G4double  ZT,
G4double  b,
G4double  bmin 
)

Definition at line 169 of file G4EMDissociationCrossSection.cc.

171{
172//
173//
174// Use Wilson et al's approach to calculate the cross-sections due to the E1
175// and E2 moments of the field at the giant dipole and quadrupole resonances
176// respectively, Note that the algorithm is traditionally applied to the
177// EMD break-up of the projectile in the field of the target, as is implemented
178// here.
179//
180// Initialise variables and calculate the energies for the GDR and GQR.
181//
182 G4double AProot3 = G4Pow::GetInstance()->powA(AP,1.0/3.0);
183 G4double u = 3.0 * J / Qprime / AProot3;
184 G4double R0 = r0 * AProot3;
185 G4double E_GDR = hbarc / std::sqrt(0.7*amu_c2*R0*R0/8.0/J*
186 (1.0 + u - (1.0 + epsilon + 3.0*u)/(1.0 + epsilon + u)*epsilon));
187 G4double E_GQR = 63.0 * MeV / AProot3;
188//
189//
190// Determine the virtual photon spectra at these energies.
191//
192 G4double ZTsq = ZT * ZT;
193 G4double nE1 = ZTsq *
195 G4double nE2 = ZTsq *
197//
198//
199// Now calculate the cross-section of the projectile for interaction with the
200// E1 and E2 fields.
201//
202 G4double sE1 = 60.0 * millibarn * MeV * (AP-ZP)*ZP/AP;
203 G4double sE2 = 0.22 * microbarn / MeV * ZP * AProot3 * AProot3;
204 if (AP > 100.0) sE2 *= 0.9;
205 else if (AP > 40.0) sE2 *= 0.6;
206 else sE2 *= 0.3;
207//
208//
209// ... and multiply with the intensity of the virtual photon spectra to get
210// the probability of interaction.
211//
212 G4PhysicsFreeVector *theCrossSectionVector = new G4PhysicsFreeVector(2);
213 theCrossSectionVector->PutValue(0, E_GDR, sE1*nE1);
214 theCrossSectionVector->PutValue(1, E_GQR, sE2*nE2*E_GQR*E_GQR);
215
216 return theCrossSectionVector;
217}
static constexpr double microbarn
Definition: G4SIunits.hh:87
static constexpr double millibarn
Definition: G4SIunits.hh:86
G4double GetGeneralE1Spectrum(G4double, G4double, G4double)
G4double GetGeneralE2Spectrum(G4double, G4double, G4double)
void PutValue(const std::size_t index, const G4double e, const G4double value)
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:230
static const G4double AP[5]
Definition: paraMaker.cc:42
float hbarc
Definition: hepunit.py:264
float amu_c2
Definition: hepunit.py:276

References source.hepunit::amu_c2, anonymous_namespace{paraMaker.cc}::AP, epsilon, G4EMDissociationSpectrum::GetGeneralE1Spectrum(), G4EMDissociationSpectrum::GetGeneralE2Spectrum(), G4Pow::GetInstance(), source.hepunit::hbarc, J, MeV, microbarn, millibarn, G4Pow::powA(), G4PhysicsFreeVector::PutValue(), Qprime, r0, and thePhotonSpectrum.

Referenced by GetCrossSectionForTarget(), and GetElementCrossSection().

◆ GetCrossSectionForTarget()

G4PhysicsFreeVector * G4EMDissociationCrossSection::GetCrossSectionForTarget ( G4double  AP,
G4double  ZP,
G4double  AT,
G4double  ZT,
G4double  b,
G4double  bmin 
)

Definition at line 222 of file G4EMDissociationCrossSection.cc.

224{
225//
226// This is a cheaky little member function to calculate the probability of
227// EMD for the target in the field of the projectile ... just by reversing the
228// A and Z's for the participants.
229//
230 return GetCrossSectionForProjectile (AT, ZT, AP, ZP, b, bmin);
231}
G4PhysicsFreeVector * GetCrossSectionForProjectile(G4double, G4double, G4double, G4double, G4double, G4double)

References anonymous_namespace{paraMaker.cc}::AP, and GetCrossSectionForProjectile().

Referenced by GetElementCrossSection().

◆ GetElementCrossSection()

G4double G4EMDissociationCrossSection::GetElementCrossSection ( const G4DynamicParticle theDynamicParticle,
G4int  Z,
const G4Material  
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 124 of file G4EMDissociationCrossSection.cc.

127{
128 // VI protection for Hydrogen
129 if(1 >= Z) { return 0.0; }
130
131 // Zero cross-section for particles with kinetic energy less than 2 MeV to prevent
132 // possible abort signal from bad arithmetic in GetCrossSectionForProjectile
133 if ( theDynamicParticle->GetKineticEnergy() < 2.0*CLHEP::MeV ) { return 0.0; }
134
135 //
136 // Get relevant information about the projectile and target (A, Z) and
137 // velocity of the projectile.
138 //
139 const G4ParticleDefinition *definitionP = theDynamicParticle->GetDefinition();
140 G4double AP = definitionP->GetBaryonNumber();
141 G4double ZP = definitionP->GetPDGCharge();
142 G4double b = theDynamicParticle->Get4Momentum().beta();
143
145 G4double ZT = (G4double)Z;
146 G4double bmin = thePhotonSpectrum->GetClosestApproach(AP, ZP, AT, ZT, b);
147 //
148 //
149 // Calculate the cross-section for the projectile and then the target. The
150 // information is returned in a G4PhysicsFreeVector, which separates out the
151 // cross-sections for the E1 and E2 moments of the virtual photon field, and
152 // the energies (GDR and GQR).
153 //
154 G4PhysicsFreeVector *theProjectileCrossSections =
155 GetCrossSectionForProjectile (AP, ZP, AT, ZT, b, bmin);
156 G4double crossSection =
157 (*theProjectileCrossSections)[0]+(*theProjectileCrossSections)[1];
158 delete theProjectileCrossSections;
159 G4PhysicsFreeVector *theTargetCrossSections =
160 GetCrossSectionForTarget (AP, ZP, AT, ZT, b, bmin);
161 crossSection +=
162 (*theTargetCrossSections)[0]+(*theTargetCrossSections)[1];
163 delete theTargetCrossSections;
164 return crossSection;
165}
G4ParticleDefinition * GetDefinition() const
G4LorentzVector Get4Momentum() const
G4double GetKineticEnergy() const
G4PhysicsFreeVector * GetCrossSectionForTarget(G4double, G4double, G4double, G4double, G4double, G4double)
G4double GetClosestApproach(const G4double, const G4double, G4double, G4double, G4double)
static G4NistManager * Instance()
G4double GetAtomicMassAmu(const G4String &symb) const
G4double GetPDGCharge() const
static constexpr double MeV

References anonymous_namespace{paraMaker.cc}::AP, CLHEP::HepLorentzVector::beta(), G4DynamicParticle::Get4Momentum(), G4NistManager::GetAtomicMassAmu(), G4ParticleDefinition::GetBaryonNumber(), G4EMDissociationSpectrum::GetClosestApproach(), GetCrossSectionForProjectile(), GetCrossSectionForTarget(), G4DynamicParticle::GetDefinition(), G4DynamicParticle::GetKineticEnergy(), G4ParticleDefinition::GetPDGCharge(), G4NistManager::Instance(), CLHEP::MeV, thePhotonSpectrum, and Z.

◆ GetIsoCrossSection()

G4double G4VCrossSectionDataSet::GetIsoCrossSection ( const G4DynamicParticle dynPart,
G4int  Z,
G4int  A,
const G4Isotope iso = nullptr,
const G4Element elm = nullptr,
const G4Material mat = nullptr 
)
virtualinherited

Reimplemented in G4ChipsAntiBaryonElasticXS, G4ChipsAntiBaryonInelasticXS, G4ChipsHyperonElasticXS, G4ChipsHyperonInelasticXS, G4ChipsKaonMinusElasticXS, G4ChipsKaonMinusInelasticXS, G4ChipsKaonPlusElasticXS, G4ChipsKaonPlusInelasticXS, G4ChipsKaonZeroElasticXS, G4ChipsKaonZeroInelasticXS, G4ChipsNeutronElasticXS, G4ChipsNeutronInelasticXS, G4ChipsPionMinusElasticXS, G4ChipsPionMinusInelasticXS, G4ChipsPionPlusElasticXS, G4ChipsPionPlusInelasticXS, G4ChipsProtonElasticXS, G4ChipsProtonInelasticXS, G4NeutronCaptureXS, G4NeutronElasticXS, G4NeutronInelasticXS, G4IonsShenCrossSection, G4ElNucleusSFcs, G4BGGNucleonElasticXS, G4BGGPionElasticXS, G4BGGPionInelasticXS, G4GammaNuclearXS, G4ParticleInelasticXS, G4BGGNucleonInelasticXS, G4LENDCombinedCrossSection, G4LENDCrossSection, G4LENDGammaCrossSection, G4ParticleHPCaptureData, G4ParticleHPElasticData, G4ParticleHPFissionData, G4ParticleHPInelasticData, G4ParticleHPThermalScatteringData, G4MuNeutrinoNucleusTotXsc, G4ElNeutrinoNucleusTotXsc, and G4PhotoNuclearCrossSection.

Definition at line 130 of file G4VCrossSectionDataSet.cc.

135{
137 ed << "GetIsoCrossSection is not implemented in <" << name << ">\n"
138 << "Particle: " << dynPart->GetDefinition()->GetParticleName()
139 << " Ekin(MeV)= " << dynPart->GetKineticEnergy()/MeV;
140 if(mat) { ed << " material: " << mat->GetName(); }
141 if(elm) { ed << " element: " << elm->GetName(); }
142 ed << " target Z= " << Z << " A= " << A << G4endl;
143 G4Exception("G4VCrossSectionDataSet::GetIsoCrossSection", "had001",
144 FatalException, ed);
145 return 0.0;
146}
@ 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 G4endl
Definition: G4ios.hh:57
const G4String & GetName() const
Definition: G4Element.hh:127
const G4String & GetName() const
Definition: G4Material.hh:173
const G4String & GetParticleName() const

References A, FatalException, G4endl, G4Exception(), G4DynamicParticle::GetDefinition(), G4DynamicParticle::GetKineticEnergy(), G4Element::GetName(), G4Material::GetName(), G4ParticleDefinition::GetParticleName(), MeV, G4VCrossSectionDataSet::name, and Z.

Referenced by G4VCrossSectionDataSet::ComputeCrossSection(), G4GammaNuclearXS::GetIsoCrossSection(), and G4GammaNuclearXS::Initialise().

◆ 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

◆ GetWilsonProbabilityForProtonDissociation()

G4double G4EMDissociationCrossSection::GetWilsonProbabilityForProtonDissociation ( G4double  A,
G4double  Z 
)

Definition at line 236 of file G4EMDissociationCrossSection.cc.

238{
239//
240// This is a simple algorithm to choose whether a proton or neutron is ejected
241// from the nucleus in the EMD interaction.
242//
243 G4double p = 0.0;
244 if (Z < 2.0)
245 p = 0.0; // To avoid to remove one proton from hydrogen isotopes
246 else if (Z < 6.0)
247 p = 0.5;
248 else if (Z < 8.0)
249 p = 0.6;
250 else if (Z < 14.0)
251 p = 0.7;
252 else
253 {
254 G4double p1 = (G4double) Z / (G4double) A;
255 G4double p2 = 1.95*G4Exp(-0.075*Z);
256 if (p1 < p2) p = p1;
257 else p = p2;
258 }
259
260 return p;
261}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179

References A, G4Exp(), and Z.

◆ IsElementApplicable()

G4bool G4EMDissociationCrossSection::IsElementApplicable ( const G4DynamicParticle part,
G4int  Z,
const G4Material  
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 105 of file G4EMDissociationCrossSection.cc.

107{
108//
109// The condition for the applicability of this class is that the projectile
110// must be an ion and the target must have more than one nucleon. In reality
111// the value of A for either the projectile or target could be much higher,
112// since for cases where both he projectile and target are medium to small
113// Z, the probability of the EMD process is, I think, VERY small.
114//
115 if (G4ParticleTable::GetParticleTable()->GetIonTable()->IsIon(part->GetDefinition())) {
116 return true;
117 } else {
118 return false;
119 }
120}
static G4ParticleTable * GetParticleTable()

References G4DynamicParticle::GetDefinition(), and G4ParticleTable::GetParticleTable().

◆ IsIsoApplicable()

G4bool G4VCrossSectionDataSet::IsIsoApplicable ( const G4DynamicParticle ,
G4int  Z,
G4int  A,
const G4Element elm = nullptr,
const G4Material mat = nullptr 
)
virtualinherited

◆ 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

◆ epsilon

G4double G4EMDissociationCrossSection::epsilon
private

◆ isForAllAtomsAndEnergies

G4bool G4VCrossSectionDataSet::isForAllAtomsAndEnergies
privateinherited

◆ J

G4double G4EMDissociationCrossSection::J
private

◆ maxKinEnergy

G4double G4VCrossSectionDataSet::maxKinEnergy
privateinherited

◆ minKinEnergy

G4double G4VCrossSectionDataSet::minKinEnergy
privateinherited

◆ name

G4String G4VCrossSectionDataSet::name
protectedinherited

◆ Qprime

G4double G4EMDissociationCrossSection::Qprime
private

◆ r0

G4double G4EMDissociationCrossSection::r0
private

◆ registry

G4CrossSectionDataSetRegistry* G4VCrossSectionDataSet::registry
privateinherited

◆ thePhotonSpectrum

G4EMDissociationSpectrum* G4EMDissociationCrossSection::thePhotonSpectrum
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(), G4BGGNucleonInelasticXS::BuildPhysicsTable(), G4LENDCrossSection::create_used_target_map(), G4BGGNucleonElasticXS::G4BGGNucleonElasticXS(), G4BGGNucleonInelasticXS::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(), G4BGGNucleonInelasticXS::GetElementCrossSection(), G4GammaNuclearXS::GetElementCrossSection(), G4ParticleInelasticXS::GetElementCrossSection(), G4BGGNucleonElasticXS::GetIsoCrossSection(), G4BGGPionElasticXS::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().

◆ xd

G4double G4EMDissociationCrossSection::xd
private

Definition at line 104 of file G4EMDissociationCrossSection.hh.

Referenced by G4EMDissociationCrossSection().


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