Geant4 10.7.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 () 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
 
const G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy, G4double logE) final
 
void BuildPhysicsTable (const G4ParticleDefinition &) final
 
void CrossSectionDescription (std::ostream &) const final
 
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 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 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 G4int GetVerboseLevel () 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
 

Static Public Member Functions

static const char * Default_Name ()
 

Additional Inherited Members

- Protected Member Functions inherited from G4VCrossSectionDataSet
void SetName (const G4String &)
 
- Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel
 

Detailed Description

Definition at line 58 of file G4NeutronElasticXS.hh.

Constructor & Destructor Documentation

◆ G4NeutronElasticXS() [1/2]

G4NeutronElasticXS::G4NeutronElasticXS ( )
explicit

Definition at line 67 of file G4NeutronElasticXS.cc.

69 ggXsection(nullptr),
70 neutron(G4Neutron::Neutron()),
71 isMaster(false)
72{
73 // verboseLevel = 0;
74 if(verboseLevel > 0){
75 G4cout << "G4NeutronElasticXS::G4NeutronElasticXS Initialise for Z < "
76 << MAXZEL << G4endl;
77 }
79 if(ggXsection == nullptr) ggXsection = new G4ComponentGGHadronNucleusXsc();
81}
const G4int MAXZEL
#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 83 of file G4NeutronElasticXS.cc.

84{
85 if(isMaster) {
86 for(G4int i=0; i<MAXZEL; ++i) {
87 delete data[i];
88 data[i] = nullptr;
89 }
90 }
91}
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 184 of file G4NeutronElasticXS.cc.

185{
186 if(verboseLevel > 0){
187 G4cout << "G4NeutronElasticXS::BuildPhysicsTable for "
188 << p.GetParticleName() << G4endl;
189 }
190 if(p.GetParticleName() != "neutron") {
192 ed << p.GetParticleName() << " is a wrong particle type -"
193 << " only neutron is allowed";
194 G4Exception("G4NeutronElasticXS::BuildPhysicsTable(..)","had012",
195 FatalException, ed, "");
196 return;
197 }
198 if(0. == coeff[0]) {
199#ifdef G4MULTITHREADED
200 G4MUTEXLOCK(&neutronElasticXSMutex);
201 if(0. == coeff[0]) {
202#endif
203 coeff[0] = 1.0;
204 isMaster = true;
205 FindDirectoryPath();
206#ifdef G4MULTITHREADED
207 }
208 G4MUTEXUNLOCK(&neutronElasticXSMutex);
209#endif
210 }
211
212 // it is possible re-initialisation for the second run
213 if(isMaster) {
214
215 // Access to elements
216 auto theCoupleTable = G4ProductionCutsTable::GetProductionCutsTable();
217 size_t numOfCouples = theCoupleTable->GetTableSize();
218 for(size_t j=0; j<numOfCouples; ++j) {
219 auto mat = theCoupleTable->GetMaterialCutsCouple(j)->GetMaterial();
220 auto elmVec = mat->GetElementVector();
221 size_t numOfElem = mat->GetNumberOfElements();
222 for (size_t ie = 0; ie < numOfElem; ++ie) {
223 G4int Z = std::max(1,std::min(((*elmVec)[ie])->GetZasInt(), MAXZEL-1));
224 if(data[Z] == nullptr) { Initialise(Z); }
225 }
226 }
227 }
228}
@ 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 G4MUTEXLOCK(mutex)
Definition: G4Threading.hh:251
#define G4MUTEXUNLOCK(mutex)
Definition: G4Threading.hh:254
const G4String & GetParticleName() const
static G4ProductionCutsTable * GetProductionCutsTable()

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 93 of file G4NeutronElasticXS.cc.

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

◆ Default_Name()

static const char * G4NeutronElasticXS::Default_Name ( )
inlinestatic

Definition at line 66 of file G4NeutronElasticXS.hh.

66{return "G4NeutronElasticXS";}

◆ GetElementCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 117 of file G4NeutronElasticXS.cc.

119{
120 G4double xs = 0.0;
121 G4double ekin = aParticle->GetKineticEnergy();
122
123 G4int Z = (ZZ >= MAXZEL) ? MAXZEL - 1 : ZZ;
124
125 auto pv = GetPhysicsVector(Z);
126 if(pv == nullptr) { return xs; }
127 // G4cout << "G4NeutronElasticXS::GetCrossSection e= " << ekin
128 // << " Z= " << Z << G4endl;
129
130 if(ekin <= pv->Energy(1)) {
131 xs = (*pv)[1];
132 } else if(ekin <= pv->GetMaxEnergy()) {
133 xs = pv->LogVectorValue(ekin, aParticle->GetLogKineticEnergy());
134 } else {
135 xs = coeff[Z]*ggXsection->GetElasticElementCrossSection(neutron,
136 ekin, Z, aeff[Z]);
137 }
138
139#ifdef G4VERBOSE
140 if(verboseLevel > 1) {
141 G4cout << "Z= " << Z << " Ekin(MeV)= " << ekin/CLHEP::MeV
142 << ", nElmXSel(b)= " << xs/CLHEP::barn
143 << G4endl;
144 }
145#endif
146 return xs;
147}
double G4double
Definition: G4Types.hh:83
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const
G4double GetElasticElementCrossSection(const G4ParticleDefinition *, G4double kinEnergy, const G4Element *)

Referenced by GetIsoCrossSection().

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

154{
155 return GetElementCrossSection(aParticle, Z, mat) * A/aeff[Z];
156}
double A(double temperature)
G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *) final

◆ IsElementApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 103 of file G4NeutronElasticXS.cc.

105{
106 return true;
107}

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 109 of file G4NeutronElasticXS.cc.

112{
113 return true;
114}

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

160{
161 size_t nIso = anElement->GetNumberOfIsotopes();
162 const G4Isotope* iso = anElement->GetIsotope(0);
163
164 //G4cout << "SelectIsotope NIso= " << nIso << G4endl;
165 if(1 == nIso) { return iso; }
166
167 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
169 G4double sum = 0.0;
170 size_t j;
171
172 // isotope wise cross section not used
173 for (j=0; j<nIso; ++j) {
174 sum += abundVector[j];
175 if(q <= sum) {
176 iso = anElement->GetIsotope(j);
177 break;
178 }
179 }
180 return iso;
181}
#define G4UniformRand()
Definition: Randomize.hh:52
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:166
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:169
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158

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