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

#include <G4HadronXSDataTable.hh>

Public Member Functions

 G4HadElementSelector (G4DynamicParticle *, G4CrossSectionDataStore *, const G4Material *, G4int bins, G4double emin, G4double emax, G4bool spline)
 
 ~G4HadElementSelector ()
 
void Dump ()
 
const G4ElementSelectRandomAtom (G4double e) const
 

Detailed Description

Definition at line 58 of file G4HadronXSDataTable.hh.

Constructor & Destructor Documentation

◆ G4HadElementSelector()

G4HadElementSelector::G4HadElementSelector ( G4DynamicParticle * dp,
G4CrossSectionDataStore * xs,
const G4Material * mat,
G4int bins,
G4double emin,
G4double emax,
G4bool spline )

Definition at line 51 of file G4HadronXSDataTable.cc.

56{
57 std::size_t n = mat->GetNumberOfElements();
58 nElmMinusOne = G4int(n - 1);
59 theElementVector = mat->GetElementVector();
60 if(nElmMinusOne > 0) {
61 G4PhysicsVector* first = nullptr;
62 xSections.resize(n, first);
63 first = new G4PhysicsLogVector(emin,emax,bins,false);
64 xSections[0] = first;
65 for(std::size_t i=1; i<n; ++i) {
66 xSections[i] = new G4PhysicsVector(*first);
67 }
68 std::vector<G4double> temp;
69 temp.resize(n, 0.0);
70 for(G4int j=0; j<=bins; ++j) {
71 G4double cross = 0.0;
72 G4double e = first->Energy(j);
73 dp->SetKineticEnergy(e);
74 for(std::size_t i=0; i<n; ++i) {
75 cross += xs->GetCrossSection(dp, (*theElementVector)[i], mat);
76 temp[i] = cross;
77 }
78 G4double fact = (cross > 0.0) ? 1.0/cross : 0.0;
79 for(std::size_t i=0; i<n; ++i) {
80 G4double y = (i<n-1) ? temp[i]*fact : 1.0;
81 xSections[i]->PutValue(j, y);
82 }
83 }
84 }
85}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
G4double GetCrossSection(const G4DynamicParticle *, const G4Material *)
void SetKineticEnergy(G4double aEnergy)
const G4ElementVector * GetElementVector() const
std::size_t GetNumberOfElements() const
G4double Energy(const std::size_t index) const

◆ ~G4HadElementSelector()

G4HadElementSelector::~G4HadElementSelector ( )

Definition at line 89 of file G4HadronXSDataTable.cc.

90{
91 if(nElmMinusOne > 0) {
92 for(G4int i=0; i<=nElmMinusOne; ++i) { delete xSections[i]; }
93 }
94}

Member Function Documentation

◆ Dump()

void G4HadElementSelector::Dump ( )

Definition at line 98 of file G4HadronXSDataTable.cc.

99{}

◆ SelectRandomAtom()

const G4Element * G4HadElementSelector::SelectRandomAtom ( G4double e) const
inline

Definition at line 70 of file G4HadronXSDataTable.hh.

71 {
72 const G4Element* element = (*theElementVector)[nElmMinusOne];
73 if (nElmMinusOne > 0) {
75 for(G4int i=0; i<nElmMinusOne; ++i) {
76 if (x <= xSections[i]->Value(e)) {
77 element = (*theElementVector)[i];
78 break;
79 }
80 }
81 }
82 return element;
83 }
#define G4UniformRand()
Definition Randomize.hh:52

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