Geant4-11
Public Member Functions | Private Member Functions | Private Attributes
G4HadElementSelector Class Reference

#include <G4HadronXSDataTable.hh>

Public Member Functions

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

Private Member Functions

 G4HadElementSelector (G4HadElementSelector &)=delete
 
G4HadElementSelectoroperator= (const G4HadElementSelector &right)=delete
 

Private Attributes

G4int nElmMinusOne
 
const G4ElementVectortheElementVector
 
std::vector< G4PhysicsVector * > xSections
 

Detailed Description

Definition at line 58 of file G4HadronXSDataTable.hh.

Constructor & Destructor Documentation

◆ G4HadElementSelector() [1/2]

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{
58 nElmMinusOne = n - 1;
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(G4int 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(G4int 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(G4int 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}
static const G4double emax
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
G4double GetCrossSection(const G4DynamicParticle *, const G4Material *)
void SetKineticEnergy(G4double aEnergy)
std::vector< G4PhysicsVector * > xSections
const G4ElementVector * theElementVector
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:186
size_t GetNumberOfElements() const
Definition: G4Material.hh:182
G4double Energy(const std::size_t index) const

References anonymous_namespace{G4HyperonSampler.cc}::bins, emax, G4PhysicsVector::Energy(), G4CrossSectionDataStore::GetCrossSection(), G4Material::GetElementVector(), G4Material::GetNumberOfElements(), CLHEP::detail::n, nElmMinusOne, G4DynamicParticle::SetKineticEnergy(), theElementVector, and xSections.

◆ ~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}

References nElmMinusOne, and xSections.

◆ G4HadElementSelector() [2/2]

G4HadElementSelector::G4HadElementSelector ( G4HadElementSelector )
privatedelete

Member Function Documentation

◆ Dump()

void G4HadElementSelector::Dump ( )

Definition at line 98 of file G4HadronXSDataTable.cc.

99{}

◆ operator=()

G4HadElementSelector & G4HadElementSelector::operator= ( const G4HadElementSelector right)
privatedelete

◆ 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

References G4UniformRand, nElmMinusOne, and xSections.

Field Documentation

◆ nElmMinusOne

G4int G4HadElementSelector::nElmMinusOne
private

◆ theElementVector

const G4ElementVector* G4HadElementSelector::theElementVector
private

Definition at line 91 of file G4HadronXSDataTable.hh.

Referenced by G4HadElementSelector().

◆ xSections

std::vector<G4PhysicsVector*> G4HadElementSelector::xSections
private

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