Geant4-11
Public Member Functions | Static Public Member Functions | Protected Attributes | Private Member Functions | Private Attributes | Static Private Attributes
G4NeutronInelasticXS Class Referencefinal

#include <G4NeutronInelasticXS.hh>

Inheritance diagram for G4NeutronInelasticXS:
G4VCrossSectionDataSet

Public Member Functions

void BuildPhysicsTable (const G4ParticleDefinition &) final
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
void CrossSectionDescription (std::ostream &) const final
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
bool ForAllAtomsAndEnergies () const
 
 G4NeutronInelasticXS ()
 
 G4NeutronInelasticXS (const G4NeutronInelasticXS &)=delete
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
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 GetMaxKinEnergy () const
 
G4double GetMinKinEnergy () const
 
const G4StringGetName () const
 
virtual G4int GetVerboseLevel () const
 
G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *) final
 
G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *, const G4Material *) final
 
G4double IsoCrossSection (G4double ekin, G4double logekin, G4int Z, G4int A)
 
G4NeutronInelasticXSoperator= (const G4NeutronInelasticXS &right)=delete
 
const G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy, G4double logE) final
 
void SetForAllAtomsAndEnergies (G4bool val)
 
void SetMaxKinEnergy (G4double value)
 
void SetMinKinEnergy (G4double value)
 
void SetName (const G4String &nam)
 
virtual void SetVerboseLevel (G4int value)
 
 ~G4NeutronInelasticXS () final
 

Static Public Member Functions

static const char * Default_Name ()
 

Protected Attributes

G4String name
 
G4int verboseLevel
 

Private Member Functions

const G4StringFindDirectoryPath ()
 
const G4PhysicsVectorGetPhysicsVector (G4int Z)
 
void Initialise (G4int Z)
 
void InitialiseOnFly (G4int Z)
 
G4PhysicsVectorRetrieveVector (std::ostringstream &in, G4bool warn)
 

Private Attributes

G4VComponentCrossSectionggXsection = nullptr
 
G4bool isForAllAtomsAndEnergies
 
G4bool isMaster = false
 
G4double maxKinEnergy
 
G4double minKinEnergy
 
const G4ParticleDefinitionneutron
 
G4CrossSectionDataSetRegistryregistry
 
std::vector< G4doubletemp
 

Static Private Attributes

static G4double coeff [MAXZINEL] = {1.0}
 
static G4ElementDatadata = nullptr
 
static G4String gDataDirectory = ""
 
static const G4int MAXZINEL = 93
 

Detailed Description

Definition at line 56 of file G4NeutronInelasticXS.hh.

Constructor & Destructor Documentation

◆ G4NeutronInelasticXS() [1/2]

G4NeutronInelasticXS::G4NeutronInelasticXS ( )
explicit

Definition at line 60 of file G4NeutronInelasticXS.cc.

63{
64 verboseLevel = 0;
65 if(verboseLevel > 0){
66 G4cout << "G4NeutronInelasticXS::G4NeutronInelasticXS Initialise for Z < "
67 << MAXZINEL << G4endl;
68 }
72}
#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 const G4int MAXZINEL
const G4ParticleDefinition * neutron
G4VComponentCrossSection * ggXsection
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
G4VCrossSectionDataSet(const G4String &nam="")
void SetForAllAtomsAndEnergies(G4bool val)

References G4cout, G4endl, G4CrossSectionDataSetRegistry::GetComponentCrossSection(), ggXsection, G4CrossSectionDataSetRegistry::Instance(), MAXZINEL, G4VCrossSectionDataSet::SetForAllAtomsAndEnergies(), and G4VCrossSectionDataSet::verboseLevel.

◆ ~G4NeutronInelasticXS()

G4NeutronInelasticXS::~G4NeutronInelasticXS ( )
final

Definition at line 74 of file G4NeutronInelasticXS.cc.

75{
76 if(isMaster) { delete data; data = nullptr; }
77}
static G4ElementData * data

References data, and isMaster.

◆ 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 247 of file G4NeutronInelasticXS.cc.

248{
249 if(verboseLevel > 0) {
250 G4cout << "G4NeutronInelasticXS::BuildPhysicsTable for "
251 << p.GetParticleName() << G4endl;
252 }
253 if(p.GetParticleName() != "neutron") {
255 ed << p.GetParticleName() << " is a wrong particle type -"
256 << " only neutron is allowed";
257 G4Exception("G4NeutronInelasticXS::BuildPhysicsTable(..)","had012",
258 FatalException, ed, "");
259 return;
260 }
261
262 if(nullptr == data) {
263#ifdef G4MULTITHREADED
264 G4MUTEXLOCK(&neutronInelasticXSMutex);
265 if(nullptr == data) {
266#endif
267 isMaster = true;
268 data = new G4ElementData();
269 data->SetName("NeutronInelastic");
271#ifdef G4MULTITHREADED
272 }
273 G4MUTEXUNLOCK(&neutronInelasticXSMutex);
274#endif
275 }
276
277 // it is possible re-initialisation for the new run
279 if(isMaster) {
280 // Upload data for elements used in geometry
281 for ( auto & elm : *table ) {
282 G4int Z = std::max( 1, std::min( elm->GetZasInt(), MAXZINEL-1) );
283 if ( nullptr == data->GetElementData(Z) ) { Initialise(Z); }
284 }
285 }
286 // prepare isotope selection
287 size_t nIso = temp.size();
288 for ( auto & elm : *table ) {
289 size_t n = elm->GetNumberOfIsotopes();
290 if(n > nIso) { nIso = n; }
291 }
292 temp.resize(nIso, 0.0);
293}
std::vector< G4Element * > G4ElementTable
@ 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
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:397
const G4String & FindDirectoryPath()
std::vector< G4double > temp
const G4String & GetParticleName() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments

References data, FatalException, FindDirectoryPath(), G4cout, G4endl, G4Exception(), G4MUTEXLOCK, G4MUTEXUNLOCK, G4ElementData::GetElementData(), G4Element::GetElementTable(), G4ParticleDefinition::GetParticleName(), Initialise(), isMaster, G4INCL::Math::max(), MAXZINEL, G4INCL::Math::min(), CLHEP::detail::n, G4ElementData::SetName(), temp, G4VCrossSectionDataSet::verboseLevel, and Z.

◆ 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
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 G4NeutronInelasticXS::CrossSectionDescription ( std::ostream &  outFile) const
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 79 of file G4NeutronInelasticXS.cc.

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

◆ Default_Name()

static const char * G4NeutronInelasticXS::Default_Name ( )
inlinestatic

Definition at line 64 of file G4NeutronInelasticXS.hh.

64{return "G4NeutronInelasticXS";}

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

◆ DumpPhysicsTable()

void G4VCrossSectionDataSet::DumpPhysicsTable ( const G4ParticleDefinition )
virtualinherited

◆ FindDirectoryPath()

const G4String & G4NeutronInelasticXS::FindDirectoryPath ( )
private

Definition at line 295 of file G4NeutronInelasticXS.cc.

296{
297 // check environment variable
298 // build the complete string identifying the file with the data set
299 if(gDataDirectory.empty()) {
300 char* path = std::getenv("G4PARTICLEXSDATA");
301 if (nullptr != path) {
302 std::ostringstream ost;
303 ost << path << "/neutron/inel";
304 gDataDirectory = ost.str();
305 } else {
306 G4Exception("G4NeutronInelasticXS::Initialise(..)","had013",
308 "Environment variable G4PARTICLEXSDATA is not defined");
309 }
310 }
311 return gDataDirectory;
312}
static G4String gDataDirectory

References FatalException, G4Exception(), and gDataDirectory.

Referenced by BuildPhysicsTable(), and Initialise().

◆ ForAllAtomsAndEnergies()

bool G4VCrossSectionDataSet::ForAllAtomsAndEnergies ( ) const
inlineinherited

◆ 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().

◆ GetElementCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 103 of file G4NeutronInelasticXS.cc.

106{
107 G4double xs = 0.0;
108 G4double ekin = aParticle->GetKineticEnergy();
109
110 G4int Z = (ZZ >= MAXZINEL) ? MAXZINEL - 1 : ZZ;
111
112 auto pv = GetPhysicsVector(Z);
113 if(pv == nullptr) { return xs; }
114 // G4cout << "G4NeutronInelasticXS::GetCrossSection e= " << ekin
115 // << " Z= " << Z << G4endl;
116
117 if(ekin <= pv->GetMaxEnergy()) {
118 xs = pv->LogVectorValue(ekin, aParticle->GetLogKineticEnergy());
119 } else {
121 ekin, Z, aeff[Z]);
122 }
123
124#ifdef G4VERBOSE
125 if(verboseLevel > 1) {
126 G4cout << "Z= " << Z << " Ekin(MeV)= " << ekin/CLHEP::MeV
127 << ", ElmXSinel(b)= " << xs/CLHEP::barn
128 << G4endl;
129 }
130#endif
131 return xs;
132}
static const G4double aeff[]
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const
const G4PhysicsVector * GetPhysicsVector(G4int Z)
static G4double coeff[MAXZINEL]
G4double GetInelasticElementCrossSection(const G4ParticleDefinition *, G4double kinEnergy, const G4Element *)
static constexpr double barn
Definition: SystemOfUnits.h:86
static constexpr double MeV

References aeff, CLHEP::barn, coeff, G4cout, G4endl, G4VComponentCrossSection::GetInelasticElementCrossSection(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetLogKineticEnergy(), GetPhysicsVector(), ggXsection, MAXZINEL, CLHEP::MeV, neutron, G4VCrossSectionDataSet::verboseLevel, and Z.

◆ 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 134 of file G4NeutronInelasticXS.cc.

139{
140 return IsoCrossSection(aParticle->GetKineticEnergy(),
141 aParticle->GetLogKineticEnergy(), Z, A);
142}
G4double IsoCrossSection(G4double ekin, G4double logekin, G4int Z, G4int A)

References A, G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetLogKineticEnergy(), IsoCrossSection(), and Z.

◆ GetMaxKinEnergy()

G4double G4VCrossSectionDataSet::GetMaxKinEnergy ( ) const
inlineinherited

◆ GetMinKinEnergy()

G4double G4VCrossSectionDataSet::GetMinKinEnergy ( ) const
inlineinherited

◆ GetName()

const G4String & G4VCrossSectionDataSet::GetName ( ) const
inlineinherited

◆ GetPhysicsVector()

const G4PhysicsVector * G4NeutronInelasticXS::GetPhysicsVector ( G4int  Z)
inlineprivate

Definition at line 123 of file G4NeutronInelasticXS.hh.

124{
125 const G4PhysicsVector* pv = data->GetElementData(Z);
126 if(pv == nullptr) {
128 pv = data->GetElementData(Z);
129 }
130 return pv;
131}

References data, G4ElementData::GetElementData(), InitialiseOnFly(), and Z.

Referenced by GetElementCrossSection(), and IsoCrossSection().

◆ GetVerboseLevel()

G4int G4VCrossSectionDataSet::GetVerboseLevel ( ) const
inlinevirtualinherited

◆ Initialise()

void G4NeutronInelasticXS::Initialise ( G4int  Z)
private

Definition at line 327 of file G4NeutronInelasticXS.cc.

328{
329 if(nullptr != data->GetElementData(Z)) { return; }
330
331 // upload element data
332 std::ostringstream ost;
333 ost << FindDirectoryPath() << Z;
334 G4PhysicsVector* v = RetrieveVector(ost, true);
336 /*
337 G4cout << "G4NeutronInelasticXS::Initialise for Z= " << Z
338 << " A= " << Amean << " Amin= " << amin[Z]
339 << " Amax= " << amax[Z] << G4endl;
340 */
341 // upload isotope data
342 if(amin[Z] < amax[Z]) {
343 G4int nmax = amax[Z] - amin[Z] + 1;
345
346 for(G4int A=amin[Z]; A<=amax[Z]; ++A) {
347 std::ostringstream ost1;
348 ost1 << gDataDirectory << Z << "_" << A;
349 G4PhysicsVector* v1 = RetrieveVector(ost1, false);
350 data->AddComponent(Z, A, v1);
351 }
352 }
353
354 // smooth transition
355 G4double sig1 = (*v)[v->GetVectorLength()-1];
356 G4double ehigh= v->GetMaxEnergy();
358 ehigh, Z, aeff[Z]);
359 coeff[Z] = (sig2 > 0.) ? sig1/sig2 : 1.0;
360}
static const G4int amax[]
static const G4int amin[]
void InitialiseForComponent(G4int Z, G4int nComponents=0)
void InitialiseForElement(G4int Z, G4PhysicsVector *v)
void AddComponent(G4int Z, G4int id, G4PhysicsVector *v)
G4PhysicsVector * RetrieveVector(std::ostringstream &in, G4bool warn)
G4double GetMaxEnergy() const
std::size_t GetVectorLength() const

References A, G4ElementData::AddComponent(), aeff, amax, amin, coeff, data, FindDirectoryPath(), gDataDirectory, G4ElementData::GetElementData(), G4VComponentCrossSection::GetInelasticElementCrossSection(), G4PhysicsVector::GetMaxEnergy(), G4PhysicsVector::GetVectorLength(), ggXsection, G4ElementData::InitialiseForComponent(), G4ElementData::InitialiseForElement(), neutron, RetrieveVector(), and Z.

Referenced by BuildPhysicsTable(), and InitialiseOnFly().

◆ InitialiseOnFly()

void G4NeutronInelasticXS::InitialiseOnFly ( G4int  Z)
private

Definition at line 314 of file G4NeutronInelasticXS.cc.

315{
316#ifdef G4MULTITHREADED
317 G4MUTEXLOCK(&neutronInelasticXSMutex);
318 if(nullptr == data->GetElementData(Z)) {
319#endif
320 Initialise(Z);
321#ifdef G4MULTITHREADED
322 }
323 G4MUTEXUNLOCK(&neutronInelasticXSMutex);
324#endif
325}

References data, G4MUTEXLOCK, G4MUTEXUNLOCK, G4ElementData::GetElementData(), Initialise(), and Z.

Referenced by GetPhysicsVector().

◆ IsElementApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 89 of file G4NeutronInelasticXS.cc.

91{
92 return true;
93}

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 96 of file G4NeutronInelasticXS.cc.

99{
100 return true;
101}

◆ IsoCrossSection()

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

Definition at line 145 of file G4NeutronInelasticXS.cc.

147{
148 G4double xs = 0.0;
149 G4int Z = (ZZ >= MAXZINEL) ? MAXZINEL - 1 : ZZ;
150
151 /*
152 G4cout << "IsoCrossSection Z= " << Z << " A= " << A
153 << " Amin= " << amin[Z] << " Amax= " << amax[Z]
154 << " E(MeV)= " << ekin << G4endl;
155 */
156 auto pv = GetPhysicsVector(Z);
157 if(pv == nullptr) { return xs; }
158
159 // compute isotope cross section if applicable
160 const G4double emax = pv->GetMaxEnergy();
161 if(ekin <= emax && amin[Z] < amax[Z] && A >= amin[Z] && A <= amax[Z]) {
162 auto pviso = data->GetComponentDataByIndex(Z, A - amin[Z]);
163 if(pviso) {
164 xs = pviso->LogVectorValue(ekin, logekin);
165#ifdef G4VERBOSE
166 if(verboseLevel > 1) {
167 G4cout << "G4NeutronInelasticXS::IsoXS: Ekin(MeV)= "
168 << ekin/CLHEP::MeV
169 << " xs(b)= " << xs/CLHEP::barn
170 << " Z= " << Z << " A= " << A << G4endl;
171 }
172#endif
173 return xs;
174 }
175 }
176
177 // use element x-section
178 if(ekin <= emax) {
179 xs = pv->LogVectorValue(ekin, logekin);
180 } else {
182 ekin, Z, aeff[Z]);
183 }
184 xs *= A/aeff[Z];
185#ifdef G4VERBOSE
186 if(verboseLevel > 1) {
187 G4cout << "G4NeutronInelasticXS::IsoXS: Z= " << Z << " A= " << A
188 << " Ekin(MeV)= " << ekin/CLHEP::MeV
189 << ", ElmXS(b)= " << xs/CLHEP::barn << G4endl;
190 }
191#endif
192 return xs;
193}
static const G4double emax
G4PhysicsVector * GetComponentDataByIndex(G4int Z, G4int idx)
G4double LogVectorValue(const G4double energy, const G4double theLogEnergy) const

References A, aeff, amax, amin, CLHEP::barn, coeff, data, emax, G4cout, G4endl, G4ElementData::GetComponentDataByIndex(), G4VComponentCrossSection::GetInelasticElementCrossSection(), GetPhysicsVector(), ggXsection, G4PhysicsVector::LogVectorValue(), MAXZINEL, CLHEP::MeV, neutron, G4VCrossSectionDataSet::verboseLevel, and Z.

Referenced by GetIsoCrossSection(), and SelectIsotope().

◆ operator=()

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

◆ RetrieveVector()

G4PhysicsVector * G4NeutronInelasticXS::RetrieveVector ( std::ostringstream &  in,
G4bool  warn 
)
private

Definition at line 363 of file G4NeutronInelasticXS.cc.

364{
365 G4PhysicsLogVector* v = nullptr;
366 std::ifstream filein(ost.str().c_str());
367 if (!filein.is_open()) {
368 if(warn) {
370 ed << "Data file <" << ost.str().c_str()
371 << "> is not opened!";
372 G4Exception("G4NeutronInelasticXS::RetrieveVector(..)","had014",
373 FatalException, ed, "Check G4PARTICLEXSDATA");
374 }
375 } else {
376 if(verboseLevel > 1) {
377 G4cout << "File " << ost.str()
378 << " is opened by G4NeutronInelasticXS" << G4endl;
379 }
380 // retrieve data from DB
381 v = new G4PhysicsLogVector();
382 if(!v->Retrieve(filein, true)) {
384 ed << "Data file <" << ost.str().c_str()
385 << "> is not retrieved!";
386 G4Exception("G4NeutronInelasticXS::RetrieveVector(..)","had015",
387 FatalException, ed, "Check G4PARTICLEXSDATA");
388 }
389 }
390 return v;
391}
G4bool Retrieve(std::ifstream &fIn, G4bool ascii=false)

References FatalException, G4cout, G4endl, G4Exception(), G4PhysicsVector::Retrieve(), and G4VCrossSectionDataSet::verboseLevel.

Referenced by Initialise().

◆ SelectIsotope()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 195 of file G4NeutronInelasticXS.cc.

197{
198 size_t nIso = anElement->GetNumberOfIsotopes();
199 const G4Isotope* iso = anElement->GetIsotope(0);
200
201 //G4cout << "SelectIsotope NIso= " << nIso << G4endl;
202 if(1 == nIso) { return iso; }
203
204 // more than 1 isotope
205 G4int Z = anElement->GetZasInt();
206 //G4cout << "SelectIsotope Z= " << Z << G4endl;
207
208 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
210 G4double sum = 0.0;
211 size_t j;
212
213 // isotope wise cross section not available
214 if(amax[Z] == amin[Z] || Z >= MAXZINEL) {
215 for (j=0; j<nIso; ++j) {
216 sum += abundVector[j];
217 if(q <= sum) {
218 iso = anElement->GetIsotope(j);
219 break;
220 }
221 }
222 return iso;
223 }
224
225 // use isotope cross sections
226 size_t nn = temp.size();
227 if(nn < nIso) { temp.resize(nIso, 0.); }
228
229 for (j=0; j<nIso; ++j) {
230 //G4cout << j << "-th isotope " << (*isoVector)[j]->GetN()
231 // << " abund= " << abundVector[j] << G4endl;
232 sum += abundVector[j]*IsoCrossSection(kinEnergy, logE, Z,
233 anElement->GetIsotope(j)->GetN());
234 temp[j] = sum;
235 }
236 sum *= q;
237 for (j = 0; j<nIso; ++j) {
238 if(temp[j] >= sum) {
239 iso = anElement->GetIsotope(j);
240 break;
241 }
242 }
243 return iso;
244}
#define G4UniformRand()
Definition: Randomize.hh:52
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170

References amax, amin, G4UniformRand, G4Element::GetIsotope(), G4Isotope::GetN(), G4Element::GetNumberOfIsotopes(), G4Element::GetRelativeAbundanceVector(), G4Element::GetZasInt(), IsoCrossSection(), MAXZINEL, G4InuclParticleNames::nn, temp, and Z.

◆ 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

◆ coeff

G4double G4NeutronInelasticXS::coeff = {1.0}
staticprivate

Definition at line 114 of file G4NeutronInelasticXS.hh.

Referenced by GetElementCrossSection(), Initialise(), and IsoCrossSection().

◆ data

G4ElementData * G4NeutronInelasticXS::data = nullptr
staticprivate

◆ gDataDirectory

G4String G4NeutronInelasticXS::gDataDirectory = ""
staticprivate

Definition at line 115 of file G4NeutronInelasticXS.hh.

Referenced by FindDirectoryPath(), and Initialise().

◆ ggXsection

G4VComponentCrossSection* G4NeutronInelasticXS::ggXsection = nullptr
private

◆ isForAllAtomsAndEnergies

G4bool G4VCrossSectionDataSet::isForAllAtomsAndEnergies
privateinherited

◆ isMaster

G4bool G4NeutronInelasticXS::isMaster = false
private

Definition at line 110 of file G4NeutronInelasticXS.hh.

Referenced by BuildPhysicsTable(), and ~G4NeutronInelasticXS().

◆ maxKinEnergy

G4double G4VCrossSectionDataSet::maxKinEnergy
privateinherited

◆ MAXZINEL

const G4int G4NeutronInelasticXS::MAXZINEL = 93
staticprivate

◆ minKinEnergy

G4double G4VCrossSectionDataSet::minKinEnergy
privateinherited

◆ name

G4String G4VCrossSectionDataSet::name
protectedinherited

◆ neutron

const G4ParticleDefinition* G4NeutronInelasticXS::neutron
private

Definition at line 106 of file G4NeutronInelasticXS.hh.

Referenced by GetElementCrossSection(), Initialise(), and IsoCrossSection().

◆ registry

G4CrossSectionDataSetRegistry* G4VCrossSectionDataSet::registry
privateinherited

◆ temp

std::vector<G4double> G4NeutronInelasticXS::temp
private

Definition at line 108 of file G4NeutronInelasticXS.hh.

Referenced by BuildPhysicsTable(), and SelectIsotope().

◆ 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(), 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(), G4ParticleInelasticXS::G4ParticleInelasticXS(), G4NeutronCaptureXS::GetElementCrossSection(), G4NeutronElasticXS::GetElementCrossSection(), 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(), IsoCrossSection(), G4GammaNuclearXS::RetrieveVector(), G4NeutronCaptureXS::RetrieveVector(), RetrieveVector(), G4ParticleInelasticXS::RetrieveVector(), and G4VCrossSectionDataSet::SetVerboseLevel().


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