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

#include <G4NeutrinoElectronNcXsc.hh>

Inheritance diagram for G4NeutrinoElectronNcXsc:
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
 
 G4NeutrinoElectronNcXsc ()
 
G4double GetBiasingFactor ()
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
G4double GetCutEnergy ()
 
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
 
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 SetBiasingFactor (G4double bf)
 
void SetCutEnergy (G4double ec)
 
void SetForAllAtomsAndEnergies (G4bool val)
 
void SetMaxKinEnergy (G4double value)
 
void SetMinKinEnergy (G4double value)
 
void SetName (const G4String &nam)
 
virtual void SetVerboseLevel (G4int value)
 
 ~G4NeutrinoElectronNcXsc ()
 

Protected Attributes

G4double fBiasingFactor
 
G4double fCofXsc
 
G4double fCutEnergy
 
G4double fSin2tW
 
G4String name
 
G4int verboseLevel
 

Private Attributes

G4bool isForAllAtomsAndEnergies
 
G4double maxKinEnergy
 
G4double minKinEnergy
 
G4CrossSectionDataSetRegistryregistry
 

Detailed Description

Definition at line 39 of file G4NeutrinoElectronNcXsc.hh.

Constructor & Destructor Documentation

◆ G4NeutrinoElectronNcXsc()

G4NeutrinoElectronNcXsc::G4NeutrinoElectronNcXsc ( )

Definition at line 41 of file G4NeutrinoElectronNcXsc.cc.

42 : G4VCrossSectionDataSet("NuElectronNcXsc")
43{
44 // PDG2016: Gf=1.1663787(6)e-5*(hc)^3/GeV^2
45 // fCofXsc = Gf*Gf*MeC2*2/pi
46
47 fCofXsc = 1.36044e-22;
49 fCofXsc /= halfpi;
50
51 // G4cout<<"hbarc = "<<hbarc/MeV/fermi<<" MeV*fermi"<<G4endl;
52
53 // PDG2016: sin^2 theta Weinberg
54
55 fSin2tW = 0.23129; // 0.2312;
56
57 fCutEnergy = 0.; // default value
58 fBiasingFactor = 1.;
59}
static constexpr double halfpi
Definition: G4SIunits.hh:57
G4VCrossSectionDataSet(const G4String &nam="")
float electron_mass_c2
Definition: hepunit.py:273
float hbarc
Definition: hepunit.py:264

References source.hepunit::electron_mass_c2, fBiasingFactor, fCofXsc, fCutEnergy, fSin2tW, halfpi, and source.hepunit::hbarc.

◆ ~G4NeutrinoElectronNcXsc()

G4NeutrinoElectronNcXsc::~G4NeutrinoElectronNcXsc ( )

Definition at line 61 of file G4NeutrinoElectronNcXsc.cc.

62{}

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

◆ GetBiasingFactor()

G4double G4NeutrinoElectronNcXsc::GetBiasingFactor ( )
inline

Definition at line 58 of file G4NeutrinoElectronNcXsc.hh.

58{return fBiasingFactor;};

References fBiasingFactor.

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

◆ GetCutEnergy()

G4double G4NeutrinoElectronNcXsc::GetCutEnergy ( )
inline

Definition at line 55 of file G4NeutrinoElectronNcXsc.hh.

55{return fCutEnergy;};

References fCutEnergy.

◆ GetElementCrossSection()

G4double G4NeutrinoElectronNcXsc::GetElementCrossSection ( const G4DynamicParticle aPart,
G4int  Z,
const G4Material  
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 89 of file G4NeutrinoElectronNcXsc.cc.

92{
93 G4double result = 0., cofL, cofR, cofL2, cofR2, cofLR;
94
96 G4String pName = aPart->GetDefinition()->GetParticleName();
97
98 if( pName == "nu_e")
99 {
100 cofL = 0.5 + fSin2tW;
101 cofR = fSin2tW;
102 }
103 else if( pName == "anti_nu_e")
104 {
105 cofL = fSin2tW;
106 cofR = 0.5 + fSin2tW;
107 }
108 else if( pName == "nu_mu")
109 {
110 cofL = -0.5 + fSin2tW;
111 cofR = fSin2tW;
112 }
113 else if( pName == "anti_nu_mu")
114 {
115 cofL = fSin2tW;
116 cofR = -0.5 + fSin2tW;
117 }
118 else if( pName == "nu_tau") // vmg: nu_tau as nu_mu ???
119 {
120 cofL = -0.5 + fSin2tW;
121 cofR = fSin2tW;
122 }
123 else if( pName == "anti_nu_tau")
124 {
125 cofL = fSin2tW;
126 cofR = -0.5 + fSin2tW;
127 }
128 else
129 {
130 return result;
131 }
132 // if( energy <= electron_mass_c2 ) return result;
133
134 cofL2 = cofL*cofL;
135 cofR2 = cofR*cofR;
136 cofLR = cofL*cofR;
137
138 if( fCutEnergy > 0. )
139 {
141 G4double tM2 = tM*tM;
142 G4double tM3 = tM*tM2;
143 G4double tC = fCutEnergy;
144 G4double tC2 = tC*tC;
145 G4double tC3 = tC*tC2;
146
147 result = (cofL2+cofR2)*(tM-tC);
148 result -= (cofR2+cofLR*0.5*electron_mass_c2/energy)*(tM2-tC2)/energy;
149 result += cofR2*(tM3-tC3)/energy/energy/3.;
150 }
151 else
152 {
153 G4double rtM = 2.*energy/(electron_mass_c2 + 2.*energy);
154 G4double rtM2 = rtM*rtM;
155 G4double rtM3 = rtM*rtM2;
156
157 result = (cofL2+cofR2)*rtM*energy;
158 result -= (cofR2*energy+cofLR*0.5*electron_mass_c2)*rtM2;
159 result += cofR2*rtM3*energy/3.;
160 }
161 // result = cofL*cofL + cofR*cofR/3.;
162 // G4cout<<"cofL2 + cofR2/3. = "<<result<<G4endl;
163 // result -= 0.5*cofL*cofR*electron_mass_c2/energy;
164
165 G4double aa = 1.;
166 G4double bb = 1.7;
167 G4double gw = 2.141*GeV;
168 G4double dd = 5000.;
169 G4double mw = 80.385*GeV;
170 G4double mz = 91.1876*GeV;
171
173 G4double totS = 2.*energy*emass + emass*emass;
174
175 if( energy > 50.*GeV )
176 {
177 result *= bb;
178 result /= 1.+ aa*totS/mz/mz;
179
180 if( pName == "anti_nu_e")
181 {
182 result *= 1. + dd*gw*gw*totS/( (totS-mw*mw)*(totS-mw*mw)+gw*gw*mw*mw );
183 }
184 }
185
186 result *= fCofXsc; //*energy;
187
188 result *= ZZ; // incoherent sum over all element electrons
189
190 result *= fBiasingFactor;
191
192 return result;
193}
static constexpr double GeV
Definition: G4SIunits.hh:203
G4ParticleDefinition * GetDefinition() const
G4double GetTotalEnergy() const
const G4String & GetParticleName() const
G4double energy(const ThreeVector &p, const G4double m)

References source.hepunit::electron_mass_c2, G4INCL::KinematicsUtils::energy(), fBiasingFactor, fCofXsc, fCutEnergy, fSin2tW, G4DynamicParticle::GetDefinition(), G4ParticleDefinition::GetParticleName(), G4DynamicParticle::GetTotalEnergy(), and GeV.

Referenced by G4NeutrinoElectronTotXsc::GetElementCrossSection().

◆ 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
static constexpr double MeV
Definition: G4SIunits.hh:200
#define G4endl
Definition: G4ios.hh:57
G4double GetKineticEnergy() const
const G4String & GetName() const
Definition: G4Element.hh:127
const G4String & GetName() const
Definition: G4Material.hh:173

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

◆ IsElementApplicable()

G4bool G4NeutrinoElectronNcXsc::IsElementApplicable ( const G4DynamicParticle aPart,
G4int  Z,
const G4Material  
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 67 of file G4NeutrinoElectronNcXsc.cc.

68{
69 G4bool result = false;
70 G4String pName = aPart->GetDefinition()->GetParticleName();
71 G4double minEnergy = 0., energy = aPart->GetTotalEnergy();
72 // Z *= 1;
73 if( fCutEnergy > 0. ) // min detected recoil electron energy
74 {
75 minEnergy = 0.5*(fCutEnergy+sqrt(fCutEnergy*(fCutEnergy+2.*electron_mass_c2)));
76 }
77 if( ( pName == "nu_e" || pName == "anti_nu_e" ||
78 pName == "nu_mu" || pName == "anti_nu_mu" ||
79 pName == "nu_tau" || pName == "anti_nu_tau" ) &&
80 energy > minEnergy )
81 {
82 result = true;
83 }
84 return result;
85}
bool G4bool
Definition: G4Types.hh:86

References source.hepunit::electron_mass_c2, G4INCL::KinematicsUtils::energy(), fCutEnergy, G4DynamicParticle::GetDefinition(), G4ParticleDefinition::GetParticleName(), and G4DynamicParticle::GetTotalEnergy().

Referenced by G4NeutrinoElectronTotXsc::IsElementApplicable().

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

◆ SetBiasingFactor()

void G4NeutrinoElectronNcXsc::SetBiasingFactor ( G4double  bf)
inline

◆ SetCutEnergy()

void G4NeutrinoElectronNcXsc::SetCutEnergy ( G4double  ec)
inline

Definition at line 54 of file G4NeutrinoElectronNcXsc.hh.

54{fCutEnergy=ec;};

References fCutEnergy.

◆ 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

◆ fBiasingFactor

G4double G4NeutrinoElectronNcXsc::fBiasingFactor
protected

◆ fCofXsc

G4double G4NeutrinoElectronNcXsc::fCofXsc
protected

Definition at line 62 of file G4NeutrinoElectronNcXsc.hh.

Referenced by G4NeutrinoElectronNcXsc(), and GetElementCrossSection().

◆ fCutEnergy

G4double G4NeutrinoElectronNcXsc::fCutEnergy
protected

◆ fSin2tW

G4double G4NeutrinoElectronNcXsc::fSin2tW
protected

Definition at line 63 of file G4NeutrinoElectronNcXsc.hh.

Referenced by G4NeutrinoElectronNcXsc(), and GetElementCrossSection().

◆ isForAllAtomsAndEnergies

G4bool G4VCrossSectionDataSet::isForAllAtomsAndEnergies
privateinherited

◆ maxKinEnergy

G4double G4VCrossSectionDataSet::maxKinEnergy
privateinherited

◆ minKinEnergy

G4double G4VCrossSectionDataSet::minKinEnergy
privateinherited

◆ name

G4String G4VCrossSectionDataSet::name
protectedinherited

◆ registry

G4CrossSectionDataSetRegistry* G4VCrossSectionDataSet::registry
privateinherited

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


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