Geant4 10.7.0
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 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 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 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 G4NeutronInelasticXS.hh.

Constructor & Destructor Documentation

◆ G4NeutronInelasticXS() [1/2]

G4NeutronInelasticXS::G4NeutronInelasticXS ( )
explicit

Definition at line 65 of file G4NeutronInelasticXS.cc.

67 neutron(G4Neutron::Neutron())
68{
69 verboseLevel = 0;
70 if(verboseLevel > 0){
71 G4cout << "G4NeutronInelasticXS::G4NeutronInelasticXS Initialise for Z < "
72 << MAXZINEL << G4endl;
73 }
75 if(ggXsection == nullptr) ggXsection = new G4ComponentGGHadronNucleusXsc();
77 isMaster = false;
78 temp.resize(13,0.0);
79}
const G4int MAXZINEL
#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 81 of file G4NeutronInelasticXS.cc.

82{
83 if(isMaster) { delete data; data = nullptr; }
84}

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

255{
256 if(verboseLevel > 0) {
257 G4cout << "G4NeutronInelasticXS::BuildPhysicsTable for "
258 << p.GetParticleName() << G4endl;
259 }
260 if(p.GetParticleName() != "neutron") {
262 ed << p.GetParticleName() << " is a wrong particle type -"
263 << " only neutron is allowed";
264 G4Exception("G4NeutronInelasticXS::BuildPhysicsTable(..)","had012",
265 FatalException, ed, "");
266 return;
267 }
268
269 if(nullptr == data) {
270#ifdef G4MULTITHREADED
271 G4MUTEXLOCK(&neutronInelasticXSMutex);
272 if(nullptr == data) {
273#endif
274 isMaster = true;
275 data = new G4ElementData();
276 data->SetName("NeutronInelastic");
277 FindDirectoryPath();
278#ifdef G4MULTITHREADED
279 }
280 G4MUTEXUNLOCK(&neutronInelasticXSMutex);
281#endif
282 }
283
284 // it is possible re-initialisation for the new run
285 if(isMaster) {
286
287 // Access to elements
288 auto theCoupleTable = G4ProductionCutsTable::GetProductionCutsTable();
289 size_t numOfCouples = theCoupleTable->GetTableSize();
290 for(size_t j=0; j<numOfCouples; ++j) {
291 auto mat = theCoupleTable->GetMaterialCutsCouple(j)->GetMaterial();
292 auto elmVec = mat->GetElementVector();
293 size_t numOfElem = mat->GetNumberOfElements();
294 for (size_t ie = 0; ie < numOfElem; ++ie) {
295 G4int Z = std::max(1,std::min(((*elmVec)[ie])->GetZasInt(), MAXZINEL-1));
296 if(nullptr == data->GetElementData(Z)) { Initialise(Z); }
297 }
298 }
299 }
300}
@ 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
G4PhysicsVector * GetElementData(G4int Z)
void SetName(const G4String &nam)
const G4String & GetParticleName() const
static G4ProductionCutsTable * GetProductionCutsTable()

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 86 of file G4NeutronInelasticXS.cc.

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

◆ Default_Name()

static const char * G4NeutronInelasticXS::Default_Name ( )
inlinestatic

Definition at line 66 of file G4NeutronInelasticXS.hh.

66{return "G4NeutronInelasticXS";}

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

◆ GetElementCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 110 of file G4NeutronInelasticXS.cc.

113{
114 G4double xs = 0.0;
115 G4double ekin = aParticle->GetKineticEnergy();
116
117 G4int Z = (ZZ >= MAXZINEL) ? MAXZINEL - 1 : ZZ;
118
119 auto pv = GetPhysicsVector(Z);
120 if(pv == nullptr) { return xs; }
121 // G4cout << "G4NeutronInelasticXS::GetCrossSection e= " << ekin
122 // << " Z= " << Z << G4endl;
123
124 if(ekin <= pv->GetMaxEnergy()) {
125 xs = pv->LogVectorValue(ekin, aParticle->GetLogKineticEnergy());
126 } else {
127 xs = coeff[Z]*ggXsection->GetInelasticElementCrossSection(neutron,
128 ekin, Z, aeff[Z]);
129 }
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 GetLogKineticEnergy() const
G4double GetKineticEnergy() const
G4double GetInelasticElementCrossSection(const G4ParticleDefinition *, G4double kinEnergy, const G4Element *)

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

146{
147 return IsoCrossSection(aParticle->GetKineticEnergy(),
148 aParticle->GetLogKineticEnergy(), Z, A);
149}
G4double IsoCrossSection(G4double ekin, G4double logekin, G4int Z, G4int A)

◆ IsElementApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 96 of file G4NeutronInelasticXS.cc.

98{
99 return true;
100}

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 103 of file G4NeutronInelasticXS.cc.

106{
107 return true;
108}

◆ IsoCrossSection()

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

Definition at line 152 of file G4NeutronInelasticXS.cc.

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

Referenced by 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 202 of file G4NeutronInelasticXS.cc.

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

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