Geant4 11.1.1
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 () final
 
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 ()
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
 
virtual G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
virtual G4double ComputeCrossSectionPerElement (G4double kinEnergy, G4double loge, const G4ParticleDefinition *, const G4Element *, const G4Material *mat=nullptr)
 
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 G4double ComputeIsoCrossSection (G4double kinEnergy, G4double loge, const G4ParticleDefinition *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
 
virtual const G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy, G4double logE)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual void CrossSectionDescription (std::ostream &) const
 
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
 
G4String name
 

Detailed Description

Definition at line 55 of file G4NeutronElasticXS.hh.

Constructor & Destructor Documentation

◆ G4NeutronElasticXS() [1/2]

G4NeutronElasticXS::G4NeutronElasticXS ( )
explicit

Definition at line 64 of file G4NeutronElasticXS.cc.

66 neutron(G4Neutron::Neutron())
67{
68 // verboseLevel = 0;
69 if(verboseLevel > 0){
70 G4cout << "G4NeutronElasticXS::G4NeutronElasticXS Initialise for Z < "
71 << MAXZEL << G4endl;
72 }
74 if(ggXsection == nullptr) ggXsection = new G4ComponentGGHadronNucleusXsc();
76}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4VComponentCrossSection * GetComponentCrossSection(const G4String &name)
static G4CrossSectionDataSetRegistry * Instance()
static const char * Default_Name()
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
void SetForAllAtomsAndEnergies(G4bool val)

◆ ~G4NeutronElasticXS()

G4NeutronElasticXS::~G4NeutronElasticXS ( )
final

Definition at line 78 of file G4NeutronElasticXS.cc.

79{
80 if(isMaster) {
81 for(G4int i=0; i<MAXZEL; ++i) {
82 delete data[i];
83 data[i] = nullptr;
84 }
85 }
86}
int G4int
Definition: G4Types.hh:85

◆ 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 192 of file G4NeutronElasticXS.cc.

193{
194 if(verboseLevel > 0){
195 G4cout << "G4NeutronElasticXS::BuildPhysicsTable for "
196 << p.GetParticleName() << G4endl;
197 }
198 if(p.GetParticleName() != "neutron") {
200 ed << p.GetParticleName() << " is a wrong particle type -"
201 << " only neutron is allowed";
202 G4Exception("G4NeutronElasticXS::BuildPhysicsTable(..)","had012",
203 FatalException, ed, "");
204 return;
205 }
206 if(0. == coeff[0]) {
207 G4AutoLock l(&nElasticXSMutex);
208 if(0. == coeff[0]) {
209 coeff[0] = 1.0;
210 isMaster = true;
211 FindDirectoryPath();
212 }
213 l.unlock();
214 }
215
216 // it is possible re-initialisation for the second run
217 if(isMaster) {
218
219 // Access to elements
221 for ( auto & elm : *table ) {
222 G4int Z = std::max( 1, std::min( elm->GetZasInt(), MAXZEL-1) );
223 if ( nullptr == data[Z] ) { Initialise(Z); }
224 }
225 }
226}
std::vector< G4Element * > G4ElementTable
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
const G4int Z[17]
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
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 119 of file G4NeutronElasticXS.cc.

123{
124 return ElementCrossSection(ekin, loge, elm->GetZasInt());
125}
G4int GetZasInt() const
Definition: G4Element.hh:132
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 147 of file G4NeutronElasticXS.cc.

152{
153 return ElementCrossSection(ekin, loge, Z)*A/aeff[Z];
154}
const G4double A[17]

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 88 of file G4NeutronElasticXS.cc.

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

◆ Default_Name()

static const char * G4NeutronElasticXS::Default_Name ( )
inlinestatic

Definition at line 63 of file G4NeutronElasticXS.hh.

63{return "G4NeutronElasticXS";}

◆ ElementCrossSection()

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

Definition at line 127 of file G4NeutronElasticXS.cc.

128{
129 G4int Z = (ZZ >= MAXZEL) ? MAXZEL - 1 : ZZ;
130 auto pv = GetPhysicsVector(Z);
131
132 G4double xs = (ekin <= pv->GetMaxEnergy()) ? pv->LogVectorValue(ekin, loge)
133 : coeff[Z]*ggXsection->GetElasticElementCrossSection(neutron, ekin,
134 Z, aeff[Z]);
135
136#ifdef G4VERBOSE
137 if(verboseLevel > 1) {
138 G4cout << "Z= " << Z << " Ekin(MeV)= " << ekin/CLHEP::MeV
139 << ", nElmXSel(b)= " << xs/CLHEP::barn
140 << G4endl;
141 }
142#endif
143 return xs;
144}
double G4double
Definition: G4Types.hh:83
G4double GetElasticElementCrossSection(const G4ParticleDefinition *, G4double kinEnergy, const G4Element *)

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 112 of file G4NeutronElasticXS.cc.

114{
115 return ElementCrossSection(aParticle->GetKineticEnergy(), aParticle->GetLogKineticEnergy(), Z);
116}
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 157 of file G4NeutronElasticXS.cc.

161{
162 return ElementCrossSection(aParticle->GetKineticEnergy(),
163 aParticle->GetLogKineticEnergy(), Z)*A/aeff[Z];
164
165}

◆ IsElementApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 98 of file G4NeutronElasticXS.cc.

100{
101 return true;
102}

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 104 of file G4NeutronElasticXS.cc.

107{
108 return false;
109}

◆ 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 167 of file G4NeutronElasticXS.cc.

169{
170 G4int nIso = (G4int)anElement->GetNumberOfIsotopes();
171 const G4Isotope* iso = anElement->GetIsotope(0);
172
173 //G4cout << "SelectIsotope NIso= " << nIso << G4endl;
174 if(1 == nIso) { return iso; }
175
176 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
178 G4double sum = 0.0;
179
180 // isotope wise cross section not used
181 for (G4int j=0; j<nIso; ++j) {
182 sum += abundVector[j];
183 if(q <= sum) {
184 iso = anElement->GetIsotope(j);
185 break;
186 }
187 }
188 return iso;
189}
#define G4UniformRand()
Definition: Randomize.hh:52
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:167
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159

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