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

#include <G4NeutronInelasticXS.hh>

+ Inheritance diagram for G4NeutronInelasticXS:

Public Member Functions

 G4NeutronInelasticXS ()
 
 ~G4NeutronInelasticXS () 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 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
 
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
 
G4double ElementCrossSection (G4double kinEnergy, G4double loge, G4int Z)
 
G4double IsoCrossSection (G4double ekin, G4double logekin, G4int Z, G4int A)
 
G4NeutronInelasticXSoperator= (const G4NeutronInelasticXS &right)=delete
 
 G4NeutronInelasticXS (const G4NeutronInelasticXS &)=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 G4NeutronInelasticXS.hh.

Constructor & Destructor Documentation

◆ G4NeutronInelasticXS() [1/2]

G4NeutronInelasticXS::G4NeutronInelasticXS ( )
explicit

Definition at line 62 of file G4NeutronInelasticXS.cc.

64 neutron(G4Neutron::Neutron()),
65 elimit(20*CLHEP::MeV)
66{
67 verboseLevel = 0;
68 if(verboseLevel > 0){
69 G4cout << "G4NeutronInelasticXS::G4NeutronInelasticXS Initialise for Z < "
70 << MAXZINEL << G4endl;
71 }
73 if(ggXsection == nullptr) ggXsection = new G4ComponentGGHadronNucleusXsc();
75}
#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)

◆ ~G4NeutronInelasticXS()

G4NeutronInelasticXS::~G4NeutronInelasticXS ( )
final

Definition at line 77 of file G4NeutronInelasticXS.cc.

78{
79 if(isMaster) { delete data; data = nullptr; }
80}

◆ G4NeutronInelasticXS() [2/2]

G4NeutronInelasticXS::G4NeutronInelasticXS ( const G4NeutronInelasticXS )
delete

Member Function Documentation

◆ BuildPhysicsTable()

void G4NeutronInelasticXS::BuildPhysicsTable ( const G4ParticleDefinition p)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 259 of file G4NeutronInelasticXS.cc.

260{
261 if(verboseLevel > 0) {
262 G4cout << "G4NeutronInelasticXS::BuildPhysicsTable for "
263 << p.GetParticleName() << G4endl;
264 }
265 if(p.GetParticleName() != "neutron") {
267 ed << p.GetParticleName() << " is a wrong particle type -"
268 << " only neutron is allowed";
269 G4Exception("G4NeutronInelasticXS::BuildPhysicsTable(..)","had012",
270 FatalException, ed, "");
271 return;
272 }
273
274 if(nullptr == data) {
275 G4AutoLock l(&nInelasticXSMutex);
276 if(nullptr == data) {
277 isMaster = true;
278 data = new G4ElementData();
279 data->SetName("NeutronInelastic");
280 FindDirectoryPath();
281 }
282 l.unlock();
283 }
284
285 // it is possible re-initialisation for the new run
287 if(isMaster) {
288 // Upload data for elements used in geometry
289 for ( auto & elm : *table ) {
290 G4int Z = std::max( 1, std::min( elm->GetZasInt(), MAXZINEL-1) );
291 if ( nullptr == data->GetElementData(Z) ) { Initialise(Z); }
292 }
293 }
294 // prepare isotope selection
295 std::size_t nIso = temp.size();
296 for ( auto & elm : *table ) {
297 std::size_t n = elm->GetNumberOfIsotopes();
298 if(n > nIso) { nIso = n; }
299 }
300 temp.resize(nIso, 0.0);
301}
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
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
G4PhysicsVector * GetElementData(G4int Z)
void SetName(const G4String &nam)
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
const G4String & GetParticleName() const

◆ ComputeCrossSectionPerElement()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 114 of file G4NeutronInelasticXS.cc.

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

◆ ComputeIsoCrossSection()

G4double G4NeutronInelasticXS::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 142 of file G4NeutronInelasticXS.cc.

147{
148 return IsoCrossSection(ekin, loge, Z, A);
149}
const G4double A[17]
G4double IsoCrossSection(G4double ekin, G4double logekin, G4int Z, G4int A)

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 82 of file G4NeutronInelasticXS.cc.

83{
84 outFile << "G4NeutronInelasticXS calculates the neutron inelastic scattering\n"
85 << "cross section on nuclei using data from the high precision\n"
86 << "neutron database. These data are simplified and smoothed over\n"
87 << "the resonance region in order to reduce CPU time.\n"
88 << "For high energy Glauber-Gribov cross section model is used\n";
89}

◆ Default_Name()

static const char * G4NeutronInelasticXS::Default_Name ( )
inlinestatic

Definition at line 63 of file G4NeutronInelasticXS.hh.

63{return "G4NeutronInelasticXS";}

Referenced by G4INCLXXNeutronBuilder::Build(), and G4NeutronCrossSectionXS::ConstructProcess().

◆ ElementCrossSection()

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

Definition at line 122 of file G4NeutronInelasticXS.cc.

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

Referenced by ComputeCrossSectionPerElement(), and GetElementCrossSection().

◆ GetElementCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 107 of file G4NeutronInelasticXS.cc.

109{
110 return ElementCrossSection(aParticle->GetKineticEnergy(), aParticle->GetLogKineticEnergy(), Z);
111}
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const

◆ GetIsoCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 152 of file G4NeutronInelasticXS.cc.

156{
157 return IsoCrossSection(aParticle->GetKineticEnergy(),
158 aParticle->GetLogKineticEnergy(), Z, A);
159}

◆ IsElementApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 92 of file G4NeutronInelasticXS.cc.

94{
95 return true;
96}

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 99 of file G4NeutronInelasticXS.cc.

102{
103 return true;
104}

◆ IsoCrossSection()

G4double G4NeutronInelasticXS::IsoCrossSection ( G4double  ekin,
G4double  logekin,
G4int  Z,
G4int  A 
)

Definition at line 162 of file G4NeutronInelasticXS.cc.

164{
165 G4double xs = 0.0;
166 G4int Z = (ZZ >= MAXZINEL) ? MAXZINEL - 1 : ZZ;
167
168 /*
169 G4cout << "IsoCrossSection Z= " << Z << " A= " << A
170 << " Amin= " << amin[Z] << " Amax= " << amax[Z]
171 << " E(MeV)= " << ekin << G4endl;
172 */
173 auto pv = GetPhysicsVector(Z);
174
175 // compute isotope cross section if applicable
176 if(ekin <= elimit && amin[Z] < amax[Z] && A >= amin[Z] && A <= amax[Z]) {
177 auto pviso = data->GetComponentDataByIndex(Z, A - amin[Z]);
178 if(nullptr != pviso) {
179 xs = pviso->LogVectorValue(ekin, logekin);
180#ifdef G4VERBOSE
181 if(verboseLevel > 1) {
182 G4cout << "G4NeutronInelasticXS::IsoXS: Ekin(MeV)= "
183 << ekin/CLHEP::MeV
184 << " xs(b)= " << xs/CLHEP::barn
185 << " Z= " << Z << " A= " << A << G4endl;
186 }
187#endif
188 return xs;
189 }
190 }
191
192 // use element x-section
193 xs = (ekin <= pv->GetMaxEnergy()) ? pv->LogVectorValue(ekin, logekin)
194 : coeff[Z]*ggXsection->GetInelasticElementCrossSection(neutron, ekin,
195 Z, aeff[Z]);
196 xs *= A/aeff[Z];
197#ifdef G4VERBOSE
198 if(verboseLevel > 1) {
199 G4cout << "G4NeutronInelasticXS::IsoXS: Z= " << Z << " A= " << A
200 << " Ekin(MeV)= " << ekin/CLHEP::MeV
201 << ", ElmXS(b)= " << xs/CLHEP::barn << G4endl;
202 }
203#endif
204 return xs;
205}
G4PhysicsVector * GetComponentDataByIndex(G4int Z, G4int idx)
G4double LogVectorValue(const G4double energy, const G4double theLogEnergy) const

Referenced by ComputeIsoCrossSection(), GetIsoCrossSection(), and SelectIsotope().

◆ operator=()

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

◆ SelectIsotope()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 207 of file G4NeutronInelasticXS.cc.

209{
210 G4int nIso = (G4int)anElement->GetNumberOfIsotopes();
211 const G4Isotope* iso = anElement->GetIsotope(0);
212
213 //G4cout << "SelectIsotope NIso= " << nIso << G4endl;
214 if(1 == nIso) { return iso; }
215
216 // more than 1 isotope
217 G4int Z = anElement->GetZasInt();
218 //G4cout << "SelectIsotope Z= " << Z << G4endl;
219
220 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
222 G4double sum = 0.0;
223 G4int j;
224
225 // isotope wise cross section not available
226 if(amax[Z] == amin[Z] || Z >= MAXZINEL) {
227 for (j=0; j<nIso; ++j) {
228 sum += abundVector[j];
229 if(q <= sum) {
230 iso = anElement->GetIsotope(j);
231 break;
232 }
233 }
234 return iso;
235 }
236
237 // use isotope cross sections
238 G4int nn = (G4int)temp.size();
239 if(nn < nIso) { temp.resize(nIso, 0.); }
240
241 for (j=0; j<nIso; ++j) {
242 //G4cout << j << "-th isotope " << (*isoVector)[j]->GetN()
243 // << " abund= " << abundVector[j] << G4endl;
244 sum += abundVector[j]*IsoCrossSection(kinEnergy, logE, Z,
245 anElement->GetIsotope(j)->GetN());
246 temp[j] = sum;
247 }
248 sum *= q;
249 for (j = 0; j<nIso; ++j) {
250 if(temp[j] >= sum) {
251 iso = anElement->GetIsotope(j);
252 break;
253 }
254 }
255 return iso;
256}
#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
G4int GetN() const
Definition: G4Isotope.hh:93

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