Geant4 11.2.2
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
 
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
G4double GetMaxEnergyDoppler() const
static G4ParticleHPManager * GetInstance()
G4VCrossSectionDataSet(const G4String &nam="")

◆ ~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 265 of file G4CrossSectionHP.cc.

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

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

◆ 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 160 of file G4CrossSectionHP.hh.

161{
162 return emax;
163}

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 214 of file G4CrossSectionHP.cc.

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

◆ SetBinSearch()

void G4CrossSectionHP::SetBinSearch ( G4int n)
inlineprotected

Definition at line 165 of file G4CrossSectionHP.hh.

166{
167 if (n > 0) { binSearch = n; }
168}

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