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

#include <G4EmElementSelector.hh>

Public Member Functions

 G4EmElementSelector (G4VEmModel *, const G4Material *, G4int bins, G4double emin, G4double emax, G4bool spline=true)
 
 ~G4EmElementSelector ()
 
void Initialise (const G4ParticleDefinition *, G4double cut=0.0)
 
void Dump (const G4ParticleDefinition *p=nullptr)
 
const G4ElementSelectRandomAtom (G4double kineticEnergy) const
 
const G4ElementSelectRandomAtom (const G4double kineticEnergy, const G4double logEKin) const
 
const G4MaterialGetMaterial () const
 

Detailed Description

Definition at line 62 of file G4EmElementSelector.hh.

Constructor & Destructor Documentation

◆ G4EmElementSelector()

G4EmElementSelector::G4EmElementSelector ( G4VEmModel mod,
const G4Material mat,
G4int  bins,
G4double  emin,
G4double  emax,
G4bool  spline = true 
)

Definition at line 54 of file G4EmElementSelector.cc.

59 :
60 model(mod), material(mat), nbins(bins), cutEnergy(-1.0),
61 lowEnergy(emin), highEnergy(emax)
62{
63 G4int n = material->GetNumberOfElements();
64 nElmMinusOne = n - 1;
65 theElementVector = material->GetElementVector();
66 if(nElmMinusOne > 0) {
67 xSections.reserve(n);
68 G4PhysicsLogVector* v0 = new G4PhysicsLogVector(lowEnergy,highEnergy,nbins);
69 xSections.push_back(v0);
70 v0->SetSpline(false);
71 for(G4int i=1; i<n; ++i) {
73 xSections.push_back(v);
74 }
75 }
76 /*
77 G4cout << "G4EmElementSelector for " << mat->GetName() << " n= " << n
78 << " nbins= " << nbins << " Emin= " << lowEnergy
79 << " Emax= " << highEnergy << G4endl;
80 */
81}
int G4int
Definition: G4Types.hh:85
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
void SetSpline(G4bool)

◆ ~G4EmElementSelector()

G4EmElementSelector::~G4EmElementSelector ( )

Definition at line 85 of file G4EmElementSelector.cc.

86{
87 if(nElmMinusOne > 0) {
88 for(G4int i=0; i<=nElmMinusOne; ++i) { delete xSections[i]; }
89 }
90}

Member Function Documentation

◆ Dump()

void G4EmElementSelector::Dump ( const G4ParticleDefinition p = nullptr)

Definition at line 151 of file G4EmElementSelector.cc.

152{
153 G4cout << "======== G4EmElementSelector for the " << model->GetName();
154 if(part) G4cout << " and " << part->GetParticleName();
155 G4cout << " for " << material->GetName() << " ========" << G4endl;
156 if(0 < nElmMinusOne) {
157 for(G4int i=0; i<nElmMinusOne; i++) {
158 G4cout << " " << (*theElementVector)[i]->GetName() << " : " << G4endl;
159 G4cout << *(xSections[i]) << G4endl;
160 }
161 }
162 G4cout << "Last Element in element vector "
163 << (*theElementVector)[nElmMinusOne]->GetName()
164 << G4endl;
165 G4cout << G4endl;
166}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
Definition: G4Material.hh:175
const G4String & GetName() const
Definition: G4VEmModel.hh:827

◆ GetMaterial()

const G4Material * G4EmElementSelector::GetMaterial ( ) const
inline

Definition at line 160 of file G4EmElementSelector.hh.

161{
162 return material;
163}

◆ Initialise()

void G4EmElementSelector::Initialise ( const G4ParticleDefinition part,
G4double  cut = 0.0 
)

Definition at line 94 of file G4EmElementSelector.cc.

96{
97 //G4cout << "G4EmElementSelector initialise for " << material->GetName()
98 // << G4endl;
99 if(0 == nElmMinusOne || cut == cutEnergy) { return; }
100
101 cutEnergy = cut;
102 //G4cout << "cut(keV)= " << cut/keV << G4endl;
103 G4double cross;
104
105 const G4double* theAtomNumDensityVector =
106 material->GetVecNbOfAtomsPerVolume();
107
108 // loop over bins
109 for(G4int j=0; j<=nbins; ++j) {
110 G4double e = (xSections[0])->Energy(j);
111 model->SetupForMaterial(part, material, e);
112 cross = 0.0;
113 //G4cout << "j= " << j << " e(MeV)= " << e/MeV << G4endl;
114 for (G4int i=0; i<=nElmMinusOne; ++i) {
115 cross += theAtomNumDensityVector[i]*
116 model->ComputeCrossSectionPerAtom(part, (*theElementVector)[i], e,
117 cutEnergy, e);
118 xSections[i]->PutValue(j, cross);
119 }
120 }
121
122 // xSections start from null, so use probabilities from the next bin
123 if(0.0 == (*xSections[nElmMinusOne])[0]) {
124 for (G4int i=0; i<=nElmMinusOne; ++i) {
125 xSections[i]->PutValue(0, (*xSections[i])[1]);
126 }
127 }
128 // xSections ends with null, so use probabilities from the previous bin
129 if(0.0 == (*xSections[nElmMinusOne])[nbins]) {
130 for (G4int i=0; i<=nElmMinusOne; ++i) {
131 xSections[i]->PutValue(nbins, (*xSections[i])[nbins-1]);
132 }
133 }
134 // perform normalization
135 for(G4int j=0; j<=nbins; ++j) {
136 cross = (*xSections[nElmMinusOne])[j];
137 // only for positive X-section
138 if(cross > 0.0) {
139 for (G4int i=0; i<nElmMinusOne; ++i) {
140 G4double x = (*xSections[i])[j]/cross;
141 xSections[i]->PutValue(j, x);
142 }
143 }
144 }
145 //G4cout << "======== G4EmElementSelector for the " << model->GetName()
146 // << G4endl;
147}
double G4double
Definition: G4Types.hh:83
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:359
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:449

◆ SelectRandomAtom() [1/2]

const G4Element * G4EmElementSelector::SelectRandomAtom ( const G4double  kineticEnergy,
const G4double  logEKin 
) const
inline

Definition at line 125 of file G4EmElementSelector.hh.

126{
127 const G4Element* element = (*theElementVector)[nElmMinusOne];
128 if (nElmMinusOne > 0) {
129 // 1. Determine energy index (only once)
130 const size_t nBins = (xSections[0])->GetVectorLength();
131 // handle cases below/above the enrgy grid (by ekin, idx that gives a=0/1)
132 // ekin = x[0] if e<=x[0] and idx will be 0 ^ a=0 => so y=y0
133 // ekin = x[N-1] if e>=x[N-1] and idx will be N-2 ^ a=1 => so y=y_{N-1}
134 const G4double ekin = std::max((xSections[0])->Energy(0),
135 std::min((xSections[0])->Energy(nBins-1),e));
136 // compute the lower index of the bin (idx \in [0,N-2] will be guaranted)
137 const size_t idx = (xSections[0])->ComputeLogVectorBin(loge);
138 // 2. Do the linear interp.(robust for corner cases through ekin, idx and a)
139 const G4double x1 = (xSections[0])->Energy(idx);
140 const G4double x2 = (xSections[0])->Energy(idx+1);
141 // note: all corner cases of the previous methods are covered and eventually
142 // gives a=0/1 that results in y=y0\y_{N-1} if e<=x[0]/e>=x[N-1] or
143 // y=y_i/y_{i+1} if e<x[i]/e>=x[i+1] due to small numerical errors
144 const G4double a = std::max(0., std::min(1., (ekin - x1)/(x2 - x1)));
145 const G4double urnd = G4UniformRand();
146 for (G4int i = 0; i < nElmMinusOne; ++i) {
147 const G4double y1 = (*xSections[i])[idx];
148 const G4double y2 = (*xSections[i])[idx+1];
149 if (urnd <= y1 + a*(y2-y1)) {
150 element = (*theElementVector)[i];
151 break;
152 }
153 }
154 }
155 return element;
156}
#define G4UniformRand()
Definition: Randomize.hh:52

◆ SelectRandomAtom() [2/2]

const G4Element * G4EmElementSelector::SelectRandomAtom ( G4double  kineticEnergy) const
inline

Definition at line 107 of file G4EmElementSelector.hh.

108{
109 const G4Element* element = (*theElementVector)[nElmMinusOne];
110 if (nElmMinusOne > 0) {
112 size_t idx(0);
113 for(G4int i=0; i<nElmMinusOne; ++i) {
114 if (x <= (xSections[i])->Value(e, idx)) {
115 element = (*theElementVector)[i];
116 break;
117 }
118 }
119 }
120 return element;
121}

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