Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4CrossSectionHP Class Reference

#include <G4CrossSectionHP.hh>

+ Inheritance diagram for G4CrossSectionHP:

Public Member Functions

 G4CrossSectionHP (const G4ParticleDefinition *, const G4String &nameData, const G4String &nameDir, G4double emaxHP, G4int zmin, G4int zmax)
 
 ~G4CrossSectionHP () override=default
 
G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *, const G4Material *) override
 
G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *, const G4Element *, const G4Material *) override
 
G4double ComputeIsoCrossSection (G4double kinEnergy, G4double loge, const G4ParticleDefinition *, G4int Z, G4int A, const G4Isotope *iso, const G4Element *elm, const G4Material *mat) override
 
const G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy, G4double logE) override
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
void DumpPhysicsTable (const G4ParticleDefinition &) override
 
G4CrossSectionHPoperator= (const G4CrossSectionHP &right)=delete
 
 G4CrossSectionHP (const G4CrossSectionHP &)=delete
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, 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 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
 

Protected Member Functions

G4double GetMaxHPEnergy () const
 
void SetBinSearch (G4int n)
 

Additional Inherited Members

- Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel {0}
 
G4String name
 

Detailed Description

Definition at line 50 of file G4CrossSectionHP.hh.

Constructor & Destructor Documentation

◆ G4CrossSectionHP() [1/2]

G4CrossSectionHP::G4CrossSectionHP ( const G4ParticleDefinition * p,
const G4String & nameData,
const G4String & nameDir,
G4double emaxHP,
G4int zmin,
G4int zmax )
explicit

Definition at line 61 of file G4CrossSectionHP.cc.

65 : G4VCrossSectionDataSet(nameData),
66 fParticle(p),
67 fNeutron(G4Neutron::Neutron()),
69 emax(emaxHP),
70 emaxT(fManagerHP->GetMaxEnergyDoppler()),
71 elimit(1.0e-11*CLHEP::eV),
72 logElimit(G4Log(elimit)),
73 minZ(zmin),
74 maxZ(zmax),
75 fDataName(nameData),
76 fDataDirectory(nameDir)
77{
78 // verboseLevel = 1;
79 if (verboseLevel > 1) {
80 G4cout << "G4CrossSectionHP::G4CrossSectionHP: Initialise for "
81 << fDataName << " " << minZ << " < Z < " << maxZ
82 << " EmaxT(MeV)=" << emaxT << G4endl;
83 G4cout << "Data directory: " << fDataDirectory << G4endl;
84 }
86 auto data = ptr->GetElementDataByName(fDataName);
87 if (nullptr == data) {
88 data = new G4ElementData(maxZ - minZ + 1);
89 data->SetName(fDataName);
90 }
91 fData = data;
92}
G4double G4Log(G4double x)
Definition G4Log.hh:227
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4ElementDataRegistry * Instance()
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
static G4ParticleHPManager * GetInstance()
G4VCrossSectionDataSet(const G4String &nam="")

Referenced by G4CrossSectionHP(), G4NeutronHPCaptureXS::G4NeutronHPCaptureXS(), G4NeutronHPElasticXS::G4NeutronHPElasticXS(), G4NeutronHPFissionXS::G4NeutronHPFissionXS(), G4NeutronHPInelasticXS::G4NeutronHPInelasticXS(), G4ParticleHPInelasticXS::G4ParticleHPInelasticXS(), and operator=().

◆ ~G4CrossSectionHP()

G4CrossSectionHP::~G4CrossSectionHP ( )
overridedefault

◆ G4CrossSectionHP() [2/2]

G4CrossSectionHP::G4CrossSectionHP ( const G4CrossSectionHP & )
delete

Member Function Documentation

◆ BuildPhysicsTable()

void G4CrossSectionHP::BuildPhysicsTable ( const G4ParticleDefinition & p)
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 270 of file G4CrossSectionHP.cc.

271{
272 if (verboseLevel > 1) {
273 G4cout << "G4CrossSectionHP::BuildPhysicsTable for "
274 << p.GetParticleName() << " and " << fDataName << G4endl;
275 }
276
277 // it is possible re-initialisation for the second run
279
280 // Access to elements
281 for ( auto const & elm : *table ) {
282 G4int Z = elm->GetZasInt();
283 if (Z >= minZ && Z <= maxZ && nullptr == fData->GetElementData(Z - minZ)) {
284 Initialise(Z);
285 }
286 }
287
288 // prepare isotope selection
289 std::size_t nmax = 0;
290 std::size_t imax = 0;
291 for ( auto const & mat : *G4Material::GetMaterialTable() ) {
292 std::size_t n = 0;
293 for ( auto const & elm : *mat->GetElementVector() ) {
294 std::size_t niso = elm->GetNumberOfIsotopes();
295 n += niso;
296 if(niso > imax) { imax = niso; }
297 }
298 if (n > nmax) { nmax = n; }
299 }
300 fTemp.resize(imax, 0.0);
301 fZA.clear();
302 fZA.reserve(nmax);
303 fIsoXS.resize(nmax, 0.0);
304}
std::vector< G4Element * > G4ElementTable
int G4int
Definition G4Types.hh:85
static const G4ElementTable * GetElementTable()
Definition G4Element.cc:401
static G4MaterialTable * GetMaterialTable()
const G4String & GetParticleName() const

◆ ComputeIsoCrossSection()

G4double G4CrossSectionHP::ComputeIsoCrossSection ( G4double kinEnergy,
G4double loge,
const G4ParticleDefinition * ,
G4int Z,
G4int A,
const G4Isotope * iso,
const G4Element * elm,
const G4Material * mat )
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 120 of file G4CrossSectionHP.cc.

126{
127 G4double ekin = e;
128 G4double loge = le;
129 if (ekin < elimit) {
130 ekin = elimit;
131 loge = logElimit;
132 }
133 if (mat != fCurrentMat) { PrepareCache(mat); }
134
135 return IsoCrossSection(ekin, loge, Z, A, mat->GetTemperature());
136}
double G4double
Definition G4Types.hh:83
const G4double A[17]
G4double GetTemperature() const

◆ DumpPhysicsTable()

void G4CrossSectionHP::DumpPhysicsTable ( const G4ParticleDefinition & )
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 306 of file G4CrossSectionHP.cc.

307{
308 if (fManagerHP->GetVerboseLevel() == 0 || fPrinted)
309 return;
310
311 //
312 // Dump element based cross section
313 // range 10e-5 eV to 20 MeV
314 // 10 point per decade
315 // in barn
316 //
317 fPrinted = true;
318 G4cout << G4endl;
319 G4cout << "HP Cross Section " << fDataName << " for "
320 << fParticle->GetParticleName() << G4endl;
321 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
322 G4cout << G4endl;
323 G4cout << "Name of Element" << G4endl;
324 G4cout << "Energy[eV] XS[barn]" << G4endl;
325 G4cout << G4endl;
326
328 for ( auto const & elm : *table ) {
329 G4int Z = elm->GetZasInt();
330 if (Z >= minZ && Z <= maxZ && nullptr != fData->GetElementData(Z - minZ)) {
331 G4cout << "---------------------------------------------------" << G4endl;
332 G4cout << elm->GetName() << G4endl;
333 std::size_t n = fData->GetNumberOfComponents(Z);
334 for (size_t i=0; i<n; ++i) {
335 G4cout << "## Z=" << Z << " A=" << fData->GetComponentID(Z - minZ, i);
336 G4cout << *(fData->GetComponentDataByIndex(Z - minZ, i)) << G4endl;
337 }
338 }
339 }
340}

◆ GetIsoCrossSection()

G4double G4CrossSectionHP::GetIsoCrossSection ( const G4DynamicParticle * dp,
G4int Z,
G4int A,
const G4Isotope * ,
const G4Element * ,
const G4Material * mat )
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 102 of file G4CrossSectionHP.cc.

107{
108 G4double ekin = dp->GetKineticEnergy();
109 G4double loge = dp->GetLogKineticEnergy();
110 if (ekin < elimit) {
111 ekin = elimit;
112 loge = logElimit;
113 }
114 if (mat != fCurrentMat) { PrepareCache(mat); }
115
116 return IsoCrossSection(ekin, loge, Z, A, mat->GetTemperature());
117}
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const

◆ GetMaxHPEnergy()

G4double G4CrossSectionHP::GetMaxHPEnergy ( ) const
inlineprotected

Definition at line 158 of file G4CrossSectionHP.hh.

159{
160 return emax;
161}

Referenced by G4NeutronHPCaptureXS::IsElementApplicable(), and G4NeutronHPFissionXS::IsElementApplicable().

◆ IsIsoApplicable()

G4bool G4CrossSectionHP::IsIsoApplicable ( const G4DynamicParticle * dp,
G4int Z,
G4int A,
const G4Element * ,
const G4Material *  )
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 94 of file G4CrossSectionHP.cc.

98{
99 return (dp->GetKineticEnergy() <= emax);
100}

◆ operator=()

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

◆ SelectIsotope()

const G4Isotope * G4CrossSectionHP::SelectIsotope ( const G4Element * elm,
G4double kinEnergy,
G4double logE )
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 219 of file G4CrossSectionHP.cc.

221{
222 std::size_t nIso = elm->GetNumberOfIsotopes();
223 const G4Isotope* iso = elm->GetIsotope(0);
224
225 //G4cout << "SelectIsotope NIso= " << nIso << G4endl;
226 if(1 == nIso) { return iso; }
227
228 // more than 1 isotope
229 G4int Z = elm->GetZasInt();
230 if (Z >= minZ && Z <= maxZ && nullptr == fData->GetElementData(Z - minZ)) {
231 Initialise(Z);
232 }
233
234 const G4double* abundVector = elm->GetRelativeAbundanceVector();
236 G4double sum = 0.0;
237
238 // is there cached isotope wise cross section?
239 std::size_t j;
240 if (Z < minZ || Z > maxZ || !CheckCache(Z) ||
241 0 == fData->GetNumberOfComponents(Z - minZ)) {
242 for (j = 0; j<nIso; ++j) {
243 sum += abundVector[j];
244 if(q <= sum) {
245 iso = elm->GetIsotope((G4int)j);
246 break;
247 }
248 }
249 return iso;
250 }
251 std::size_t nn = fTemp.size();
252 if (nn < nIso) { fTemp.resize(nIso, 0.); }
253
254 // reuse cache
255 for (j=0; j<nIso; ++j) {
256 sum += abundVector[j]*
257 GetCrossSection(Z - minZ, elm->GetIsotope((G4int)j)->GetN());
258 fTemp[j] = sum;
259 }
260 sum *= q;
261 for (j = 0; j<nIso; ++j) {
262 if (fTemp[j] >= sum) {
263 iso = elm->GetIsotope((G4int)j);
264 break;
265 }
266 }
267 return iso;
268}
#define G4UniformRand()
Definition Randomize.hh:52
G4double * GetRelativeAbundanceVector() const
Definition G4Element.hh:149
std::size_t GetNumberOfIsotopes() const
Definition G4Element.hh:143
const G4Isotope * GetIsotope(G4int iso) const
Definition G4Element.hh:151
G4int GetZasInt() const
Definition G4Element.hh:120
G4int GetN() const
Definition G4Isotope.hh:83

◆ SetBinSearch()

void G4CrossSectionHP::SetBinSearch ( G4int n)
inlineprotected

Definition at line 163 of file G4CrossSectionHP.hh.

164{
165 if (n > 0) { binSearch = n; }
166}

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