Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NeutronElasticXS Class Referencefinal

#include <G4NeutronElasticXS.hh>

Inheritance diagram for G4NeutronElasticXS:

Public Member Functions

 G4NeutronElasticXS ()
 ~G4NeutronElasticXS () override=default
G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *) final
G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *, const G4Material *) final
G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *) final
G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso, const G4Element *elm, const G4Material *mat) final
G4double ComputeCrossSectionPerElement (G4double kinEnergy, G4double loge, const G4ParticleDefinition *, const G4Element *, const G4Material *) final
G4double ComputeIsoCrossSection (G4double kinEnergy, G4double loge, const G4ParticleDefinition *, G4int Z, G4int A, const G4Isotope *iso, const G4Element *elm, const G4Material *mat) final
const G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy, G4double logE) final
void BuildPhysicsTable (const G4ParticleDefinition &) final
void CrossSectionDescription (std::ostream &) const final
G4double ElementCrossSection (G4double kinEnergy, G4double loge, G4int Z)
G4NeutronElasticXSoperator= (const G4NeutronElasticXS &right)=delete
 G4NeutronElasticXS (const G4NeutronElasticXS &)=delete
Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
virtual ~G4VCrossSectionDataSet ()
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
virtual void SetVerboseLevel (G4int value)
G4double GetMinKinEnergy () const
void SetMinKinEnergy (G4double value)
G4double GetMaxKinEnergy () const
void SetMaxKinEnergy (G4double value)
bool ForAllAtomsAndEnergies () const
void SetForAllAtomsAndEnergies (G4bool val)
const G4StringGetName () const
void SetName (const G4String &nam)
G4VCrossSectionDataSetoperator= (const G4VCrossSectionDataSet &right)=delete
 G4VCrossSectionDataSet (const G4VCrossSectionDataSet &)=delete

Static Public Member Functions

static const char * Default_Name ()

Additional Inherited Members

Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel {0}
G4String name

Detailed Description

Definition at line 55 of file G4NeutronElasticXS.hh.

Constructor & Destructor Documentation

◆ G4NeutronElasticXS() [1/2]

G4NeutronElasticXS::G4NeutronElasticXS ( )

Definition at line 68 of file G4NeutronElasticXS.cc.

70 neutron(G4Neutron::Neutron())
71{
72 if (verboseLevel > 0) {
73 G4cout << "G4NeutronElasticXS::G4NeutronElasticXS Initialise for Z < "
74 << MAXZEL << G4endl;
75 }
76 ggXsection =
78 if (ggXsection == nullptr)
79 ggXsection = new G4ComponentGGHadronNucleusXsc();
81 if (nullptr == data) {
82 data = new G4ElementData(MAXZEL);
83 data->SetName("nElastic");
84 dataR = new G4ElementData(MAXZEL);
85 dataR->SetName("nRElastic");
86 FindDirectoryPath();
87 }
88}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4VComponentCrossSection * GetComponentCrossSection(const G4String &name)
static G4CrossSectionDataSetRegistry * Instance()
static const char * Default_Name()
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
G4VCrossSectionDataSet(const G4String &nam="")
void SetForAllAtomsAndEnergies(G4bool val)

Referenced by G4NeutronElasticXS(), and operator=().

◆ ~G4NeutronElasticXS()

G4NeutronElasticXS::~G4NeutronElasticXS ( )
overridedefault

◆ G4NeutronElasticXS() [2/2]

G4NeutronElasticXS::G4NeutronElasticXS ( const G4NeutronElasticXS & )
delete

Member Function Documentation

◆ BuildPhysicsTable()

void G4NeutronElasticXS::BuildPhysicsTable ( const G4ParticleDefinition & p)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 210 of file G4NeutronElasticXS.cc.

211{
212 if(verboseLevel > 0){
213 G4cout << "G4NeutronElasticXS::BuildPhysicsTable for "
214 << p.GetParticleName() << G4endl;
215 }
216 if(p.GetParticleName() != "neutron") {
218 ed << p.GetParticleName() << " is a wrong particle type -"
219 << " only neutron is allowed";
220 G4Exception("G4NeutronElasticXS::BuildPhysicsTable(..)","had012",
221 FatalException, ed, "");
222 return;
223 }
224
226
227 // initialise static tables only once
228 std::call_once(applyOnce, [this]() { isInitializer = true; });
229
230 if (isInitializer) {
231 G4AutoLock l(&nElasticXSMutex);
232
233 // Access to elements
235 for ( auto const & elm : *table ) {
236 G4int Z = std::max( 1, std::min( elm->GetZasInt(), MAXZEL-1) );
237 if ( nullptr == data->GetElementData(Z) ) { Initialise(Z); }
238 }
239 l.unlock();
240 }
241}
G4TemplateAutoLock< G4Mutex > G4AutoLock
std::vector< G4Element * > G4ElementTable
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
int G4int
Definition G4Types.hh:85
static const G4ElementTable * GetElementTable()
Definition G4Element.cc:401
static G4HadronicParameters * Instance()
const G4String & GetParticleName() const

◆ ComputeCrossSectionPerElement()

G4double G4NeutronElasticXS::ComputeCrossSectionPerElement ( G4double kinEnergy,
G4double loge,
const G4ParticleDefinition * ,
const G4Element * elm,
const G4Material *  )
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 122 of file G4NeutronElasticXS.cc.

126{
127 return ElementCrossSection(ekin, loge, elm->GetZasInt());
128}
G4int GetZasInt() const
Definition G4Element.hh:120
G4double ElementCrossSection(G4double kinEnergy, G4double loge, G4int Z)

◆ ComputeIsoCrossSection()

G4double G4NeutronElasticXS::ComputeIsoCrossSection ( G4double kinEnergy,
G4double loge,
const G4ParticleDefinition * ,
G4int Z,
G4int A,
const G4Isotope * iso,
const G4Element * elm,
const G4Material * mat )
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 164 of file G4NeutronElasticXS.cc.

169{
170 G4int Z = std::min(ZZ, MAXZEL-1);
171 return ElementCrossSection(ekin, loge, Z)*A/aeff[Z];
172}
const G4double A[17]

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 90 of file G4NeutronElasticXS.cc.

91{
92 outFile << "G4NeutronElasticXS calculates the neutron elastic scattering\n"
93 << "cross section on nuclei using data from the high precision\n"
94 << "neutron database. These data are simplified and smoothed over\n"
95 << "the resonance region in order to reduce CPU time.\n"
96 << "For high energies Glauber-Gribiv cross section is used.\n";
97}

◆ Default_Name()

const char * G4NeutronElasticXS::Default_Name ( )
inlinestatic

Definition at line 63 of file G4NeutronElasticXS.hh.

63{return "G4NeutronElasticXS";}

Referenced by G4NeutronElasticXS().

◆ ElementCrossSection()

G4double G4NeutronElasticXS::ElementCrossSection ( G4double kinEnergy,
G4double loge,
G4int Z )

Definition at line 131 of file G4NeutronElasticXS.cc.

132{
133 G4int Z = std::min(ZZ, MAXZEL-1);
134 G4double xs;
135 G4bool done{false};
136
137 // data from the resonance region
138 if (fRfilesEnabled) {
139 auto pv = GetPhysicsVectorR(Z);
140 if (nullptr != pv && ekin < el_max_r_e[Z]) {
141 xs = pv->LogVectorValue(ekin, loge);
142 done = true;
143 }
144 }
145 // data above the resonance region
146 if (!done) {
147 auto pv = GetPhysicsVector(Z);
148 xs = (ekin <= pv->GetMaxEnergy()) ? pv->LogVectorValue(ekin, loge)
149 : coeff[Z]*ggXsection->GetElasticElementCrossSection(neutron, ekin,
150 Z, aeff[Z]);
151 }
152
153#ifdef G4VERBOSE
154 if (verboseLevel > 1) {
155 G4cout << "Z= " << Z << " Ekin(MeV)= " << ekin/CLHEP::MeV
156 << ", nElmXSel(b)= " << xs/CLHEP::barn
157 << G4endl;
158 }
159#endif
160 return xs;
161}
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86

Referenced by ComputeCrossSectionPerElement(), ComputeIsoCrossSection(), GetElementCrossSection(), and GetIsoCrossSection().

◆ GetElementCrossSection()

G4double G4NeutronElasticXS::GetElementCrossSection ( const G4DynamicParticle * aParticle,
G4int Z,
const G4Material *  )
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 114 of file G4NeutronElasticXS.cc.

116{
117 return ElementCrossSection(aParticle->GetKineticEnergy(),
118 aParticle->GetLogKineticEnergy(), Z);
119}
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const

◆ GetIsoCrossSection()

G4double G4NeutronElasticXS::GetIsoCrossSection ( const G4DynamicParticle * aParticle,
G4int Z,
G4int A,
const G4Isotope * iso,
const G4Element * elm,
const G4Material * mat )
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 175 of file G4NeutronElasticXS.cc.

179{
180 G4int Z = std::min(ZZ, MAXZEL-1);
181 return ElementCrossSection(aParticle->GetKineticEnergy(),
182 aParticle->GetLogKineticEnergy(), Z)*A/aeff[Z];
183
184}

◆ IsElementApplicable()

G4bool G4NeutronElasticXS::IsElementApplicable ( const G4DynamicParticle * ,
G4int Z,
const G4Material *  )
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 100 of file G4NeutronElasticXS.cc.

102{
103 return true;
104}

◆ IsIsoApplicable()

G4bool G4NeutronElasticXS::IsIsoApplicable ( const G4DynamicParticle * ,
G4int Z,
G4int A,
const G4Element * ,
const G4Material *  )
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 106 of file G4NeutronElasticXS.cc.

109{
110 return false;
111}

◆ operator=()

G4NeutronElasticXS & G4NeutronElasticXS::operator= ( const G4NeutronElasticXS & right)
delete

◆ SelectIsotope()

const G4Isotope * G4NeutronElasticXS::SelectIsotope ( const G4Element * anElement,
G4double kinEnergy,
G4double logE )
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 186 of file G4NeutronElasticXS.cc.

188{
189 G4int nIso = (G4int)anElement->GetNumberOfIsotopes();
190 const G4Isotope* iso = anElement->GetIsotope(0);
191
192 if(1 == nIso) { return iso; }
193
194 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
196 G4double sum = 0.0;
197
198 // isotope wise cross section not used
199 for (G4int j=0; j<nIso; ++j) {
200 sum += abundVector[j];
201 if(q <= sum) {
202 iso = anElement->GetIsotope(j);
203 break;
204 }
205 }
206 return iso;
207}
#define G4UniformRand()
Definition Randomize.hh:52
G4double * GetRelativeAbundanceVector() const
Definition G4Element.hh:149
std::size_t GetNumberOfIsotopes() const
Definition G4Element.hh:143
const G4Isotope * GetIsotope(G4int iso) const
Definition G4Element.hh:151

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