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

#include <G4HadronXSDataTable.hh>

Public Member Functions

 G4HadronXSDataTable ()
 
 ~G4HadronXSDataTable ()
 
void Initialise (G4DynamicParticle *, G4CrossSectionDataStore *, G4int bins, G4double emin, G4double emax, G4bool spline)
 
const G4PhysicsVectorHasData (size_t idx) const
 
G4double GetCrossSection (G4double e, size_t idx) const
 
const G4ElementSelectRandomAtom (G4double e, size_t idx) const
 
void Dump ()
 

Detailed Description

Definition at line 95 of file G4HadronXSDataTable.hh.

Constructor & Destructor Documentation

◆ G4HadronXSDataTable()

G4HadronXSDataTable::G4HadronXSDataTable ( )
explicit

Definition at line 104 of file G4HadronXSDataTable.cc.

104 : nMaterials(0)
105{}

◆ ~G4HadronXSDataTable()

G4HadronXSDataTable::~G4HadronXSDataTable ( )

Definition at line 149 of file G4HadronXSDataTable.cc.

150{
151 for(size_t i=0; i<nMaterials; ++i) {
152 delete xsData[i];
153 delete elmSelectors[i];
154 }
155}

Member Function Documentation

◆ Dump()

void G4HadronXSDataTable::Dump ( )

Definition at line 159 of file G4HadronXSDataTable.cc.

160{}

◆ GetCrossSection()

G4double G4HadronXSDataTable::GetCrossSection ( G4double  e,
size_t  idx 
) const
inline

Definition at line 113 of file G4HadronXSDataTable.hh.

114 {
115 return xsData[idx]->Value(e);
116 };

◆ HasData()

const G4PhysicsVector * G4HadronXSDataTable::HasData ( size_t  idx) const
inline

Definition at line 108 of file G4HadronXSDataTable.hh.

109 {
110 return xsData[idx];
111 };

◆ Initialise()

void G4HadronXSDataTable::Initialise ( G4DynamicParticle dp,
G4CrossSectionDataStore xs,
G4int  bins,
G4double  emin,
G4double  emax,
G4bool  spline 
)

Definition at line 109 of file G4HadronXSDataTable.cc.

113{
115 if(nn > nMaterials) {
116 G4PhysicsLogVector* first = nullptr;
117 G4int sbins = std::max(10, bins/5);
119 for(size_t i=nMaterials; i<nn; ++i) {
120 const G4Material* mat = (*mtable)[i];
121 G4PhysicsVector* v = nullptr;
122 G4HadElementSelector* es = nullptr;
123 // create real vector only for complex materials
124 if(mat->GetNumberOfElements() > 1) {
125 if(!first) {
126 first = new G4PhysicsLogVector(emin, emax, bins);
127 first->SetSpline(spline);
128 v = first;
129 } else {
130 v = new G4PhysicsVector(*first);
131 }
132 for(G4int j=0; j<=bins; ++j) {
133 G4double e = first->Energy(j);
134 dp->SetKineticEnergy(e);
135 G4double cros = xs->ComputeCrossSection(dp, mat);
136 v->PutValue(j, cros);
137 }
138 elmSelectors[i] = new G4HadElementSelector(dp, xs, mat, sbins, emin, emax, spline);
139 }
140 xsData.push_back(v);
141 elmSelectors.push_back(es);
142 }
143 nMaterials = nn;
144 }
145}
std::vector< G4Material * > G4MaterialTable
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Material *)
void SetKineticEnergy(G4double aEnergy)
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:644
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:637
G4double Energy(std::size_t index) const
void PutValue(std::size_t index, G4double theValue)
void SetSpline(G4bool)

◆ SelectRandomAtom()

const G4Element * G4HadronXSDataTable::SelectRandomAtom ( G4double  e,
size_t  idx 
) const
inline

Definition at line 118 of file G4HadronXSDataTable.hh.

119 {
120 return elmSelectors[idx]->SelectRandomAtom(e);
121 };

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