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

#include <G4NeutrinoElectronCcXsc.hh>

Inheritance diagram for G4NeutrinoElectronCcXsc:
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
 
 G4NeutrinoElectronCcXsc ()
 
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)
 
 ~G4NeutrinoElectronCcXsc ()
 

Protected Attributes

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

Private Attributes

G4bool isForAllAtomsAndEnergies
 
G4double maxKinEnergy
 
G4double minKinEnergy
 
G4CrossSectionDataSetRegistryregistry
 

Detailed Description

Definition at line 41 of file G4NeutrinoElectronCcXsc.hh.

Constructor & Destructor Documentation

◆ G4NeutrinoElectronCcXsc()

G4NeutrinoElectronCcXsc::G4NeutrinoElectronCcXsc ( )

Definition at line 43 of file G4NeutrinoElectronCcXsc.cc.

44 : G4VCrossSectionDataSet("NuElectronCcXsc")
45{
46 // PDG2016: Gf=1.1663787(6)e-5*(hc)^3/GeV^2
47 // fCofXsc = Gf*Gf*MeC2*2/pi
48
49 fCofXsc = 1.36044e-22;
51 fCofXsc /= halfpi;
52
53 // G4cout<<"fCofXsc = "<<fCofXsc*GeV/cm2<<" cm2/GeV"<<G4endl;
54
55 // G4cout<<"hbarc = "<<hbarc/MeV/fermi<<" MeV*fermi"<<G4endl;
56
57 // PDG2016: sin^2 theta Weinberg
58
59 fSin2tW = 0.23129; // 0.2312;
60
61 fCutEnergy = 0.; // default value
62
63 fBiasingFactor = 1.; // default as physics
64
67}
static constexpr double halfpi
Definition: G4SIunits.hh:57
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:99
G4ParticleDefinition * theMuonMinus
G4ParticleDefinition * theTauMinus
static G4TauMinus * TauMinus()
Definition: G4TauMinus.cc:134
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, source.hepunit::hbarc, G4MuonMinus::MuonMinus(), G4TauMinus::TauMinus(), theMuonMinus, and theTauMinus.

◆ ~G4NeutrinoElectronCcXsc()

G4NeutrinoElectronCcXsc::~G4NeutrinoElectronCcXsc ( )

Definition at line 69 of file G4NeutrinoElectronCcXsc.cc.

70{}

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 G4NeutrinoElectronCcXsc::GetBiasingFactor ( )
inline

Definition at line 60 of file G4NeutrinoElectronCcXsc.hh.

60{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 G4NeutrinoElectronCcXsc::GetCutEnergy ( )
inline

Definition at line 57 of file G4NeutrinoElectronCcXsc.hh.

57{return fCutEnergy;};

References fCutEnergy.

◆ GetElementCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 99 of file G4NeutrinoElectronCcXsc.cc.

102{
103 G4double result = 0., totS, fmass, fmass2, emass=electron_mass_c2, emass2;
104
105 G4double energy = aPart->GetTotalEnergy();
106 G4String pName = aPart->GetDefinition()->GetParticleName();
107
108 emass2 = emass*emass;
109 totS = 2.*energy*emass + emass2;
110
111 if( pName == "anti_nu_e" || pName == "nu_mu")
112 {
113 fmass = theMuonMinus->GetPDGMass();
114 fmass2 = fmass*fmass;
115 result = (1. - fmass2/totS)*(1. - fmass2/totS);
116 }
117 else if( pName == "anti_nu_mu")
118 {
119 fmass = theMuonMinus->GetPDGMass();
120 fmass2 = fmass*fmass;
121
122 result = (1.+ emass2/totS)*(1.+ fmass2/totS);
123 result += (1.- emass2/totS)*(1.- fmass2/totS)/3.;
124 result *= 0.25*(1. - fmass2/totS)*(1. - fmass2/totS);
125 }
126 else if( pName == "nu_tau")
127 {
128 fmass = theTauMinus->GetPDGMass();
129 fmass2 = fmass*fmass;
130 result = (1. - fmass2/totS)*(1. - fmass2/totS);
131 }
132 else if( pName == "anti_nu_tau")
133 {
134 fmass = theTauMinus->GetPDGMass();
135 fmass2 = fmass*fmass;
136
137 result = (1.+ emass2/totS)*(1.+ fmass2/totS);
138 result += (1.- emass2/totS)*(1.- fmass2/totS)/3.;
139 result *= 0.25*(1. - fmass2/totS)*(1. - fmass2/totS);
140 }
141 else
142 {
143 return result;
144 }
145 // if( energy <= electron_mass_c2 ) return result;
146
147 G4double aa = 1.;
148 G4double bb = 1.7;
149 G4double gw = 2.141*GeV;
150 G4double dd = 5000.;
151 G4double mw = 80.385*GeV;
152
153 if( energy > 50.*GeV )
154 {
155 result *= bb;
156 result /= 1.+ aa*totS/mw/mw;
157
158 if( pName == "anti_nu_e")
159 {
160 result *= 1. + dd*gw*gw*totS/( (totS-mw*mw)*(totS-mw*mw)+gw*gw*mw*mw );
161 }
162 }
163 result *= fCofXsc; //*energy;
164 result *= energy + 0.5*emass;
165 result *= ZZ; // incoherent sum over all element electrons
166
167 result *= fBiasingFactor; // biasing up, if set >1
168
169 return result;
170}
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, G4DynamicParticle::GetDefinition(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetPDGMass(), G4DynamicParticle::GetTotalEnergy(), GeV, theMuonMinus, and theTauMinus.

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 G4NeutrinoElectronCcXsc::IsElementApplicable ( const G4DynamicParticle aPart,
G4int  Z,
const G4Material  
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 75 of file G4NeutrinoElectronCcXsc.cc.

76{
77 G4bool result = false;
78 G4String pName = aPart->GetDefinition()->GetParticleName();
79 G4double minEnergy = 0., energy = aPart->GetTotalEnergy();
80 G4double fmass, emass = electron_mass_c2;
81
82 if( pName == "anti_nu_e" || pName == "nu_mu" || pName == "anti_nu_mu" ) fmass = theMuonMinus->GetPDGMass();
83 else if( pName == "nu_tau" || pName == "anti_nu_tau" ) fmass = theTauMinus->GetPDGMass();
84 else fmass = emass;
85
86 minEnergy = (fmass-emass)*(fmass+emass)/emass;
87
88 if( ( pName == "nu_mu" || pName == "anti_nu_mu" ||
89 pName == "nu_tau" || pName == "anti_nu_tau" ) &&
90 energy > minEnergy )
91 {
92 result = true;
93 }
94 return result;
95}
bool G4bool
Definition: G4Types.hh:86

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

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 G4NeutrinoElectronCcXsc::SetBiasingFactor ( G4double  bf)
inline

◆ SetCutEnergy()

void G4NeutrinoElectronCcXsc::SetCutEnergy ( G4double  ec)
inline

Definition at line 56 of file G4NeutrinoElectronCcXsc.hh.

56{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 G4NeutrinoElectronCcXsc::fBiasingFactor
protected

◆ fCofXsc

G4double G4NeutrinoElectronCcXsc::fCofXsc
protected

Definition at line 64 of file G4NeutrinoElectronCcXsc.hh.

Referenced by G4NeutrinoElectronCcXsc(), and GetElementCrossSection().

◆ fCutEnergy

G4double G4NeutrinoElectronCcXsc::fCutEnergy
protected

Definition at line 66 of file G4NeutrinoElectronCcXsc.hh.

Referenced by G4NeutrinoElectronCcXsc(), GetCutEnergy(), and SetCutEnergy().

◆ fSin2tW

G4double G4NeutrinoElectronCcXsc::fSin2tW
protected

Definition at line 65 of file G4NeutrinoElectronCcXsc.hh.

Referenced by G4NeutrinoElectronCcXsc().

◆ 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

◆ theMuonMinus

G4ParticleDefinition* G4NeutrinoElectronCcXsc::theMuonMinus
protected

◆ theTauMinus

G4ParticleDefinition* G4NeutrinoElectronCcXsc::theTauMinus
protected

◆ 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: