Geant4 11.2.2
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 (const G4double kineticEnergy, const G4double logEKin) const
 
const G4ElementSelectRandomAtom (G4double kineticEnergy) const
 
const G4MaterialGetMaterial () const
 
G4EmElementSelectoroperator= (const G4EmElementSelector &right)=delete
 
 G4EmElementSelector (const G4EmElementSelector &)=delete
 

Detailed Description

Definition at line 62 of file G4EmElementSelector.hh.

Constructor & Destructor Documentation

◆ G4EmElementSelector() [1/2]

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 = (G4int)material->GetNumberOfElements();
64 nElmMinusOne = n - 1;
65 theElementVector = material->GetElementVector();
66 if(nElmMinusOne > 0) {
67 xSections.reserve(n);
68 auto v0 = new G4PhysicsLogVector(lowEnergy,highEnergy,nbins,false);
69 xSections.push_back(v0);
70 for(G4int i=1; i<n; ++i) {
71 auto v = new G4PhysicsLogVector(*v0);
72 xSections.push_back(v);
73 }
74 }
75 /*
76 G4cout << "G4EmElementSelector for " << mat->GetName() << " n= " << n
77 << " nbins= " << nbins << " Emin= " << lowEnergy
78 << " Emax= " << highEnergy << G4endl;
79 */
80}
int G4int
Definition G4Types.hh:85
const G4ElementVector * GetElementVector() const
std::size_t GetNumberOfElements() const

◆ ~G4EmElementSelector()

G4EmElementSelector::~G4EmElementSelector ( )

Definition at line 84 of file G4EmElementSelector.cc.

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

◆ G4EmElementSelector() [2/2]

G4EmElementSelector::G4EmElementSelector ( const G4EmElementSelector & )
delete

Member Function Documentation

◆ Dump()

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

Definition at line 192 of file G4EmElementSelector.cc.

193{
194 G4cout << "======== G4EmElementSelector for the " << model->GetName();
195 if(part) G4cout << " and " << part->GetParticleName();
196 G4cout << " for " << material->GetName() << " ========" << G4endl;
197 if(0 < nElmMinusOne) {
198 for(G4int i=0; i<nElmMinusOne; i++) {
199 G4cout << " " << (*theElementVector)[i]->GetName() << " : " << G4endl;
200 G4cout << *(xSections[i]) << G4endl;
201 }
202 }
203 G4cout << "Last Element in element vector "
204 << (*theElementVector)[nElmMinusOne]->GetName()
205 << G4endl;
206 G4cout << G4endl;
207}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
const G4String & GetName() const

◆ GetMaterial()

const G4Material * G4EmElementSelector::GetMaterial ( ) const
inline

Definition at line 126 of file G4EmElementSelector.hh.

127{
128 return material;
129}

◆ Initialise()

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

Definition at line 93 of file G4EmElementSelector.cc.

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

◆ operator=()

G4EmElementSelector & G4EmElementSelector::operator= ( const G4EmElementSelector & right)
delete

◆ SelectRandomAtom() [1/2]

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

Definition at line 156 of file G4EmElementSelector.cc.

158{
159 const G4Element* element = (*theElementVector)[nElmMinusOne];
160 if (nElmMinusOne > 0) {
161 // 1. Determine energy index (only once)
162 // handle cases below/above the enrgy grid (by ekin, idx that gives a=0/1)
163 // ekin = x[0] if e<=x[0] and idx will be 0 ^ a=0 => so y=y0
164 // ekin = x[N-1] if e>=x[N-1] and idx will be N-2 ^ a=1 => so y=y_{N-1}
165 G4double ekin = e;
166 std::size_t idx = 0;
167 if(e <= (xSections[0])->Energy(0)) {
168 ekin = (xSections[0])->Energy(0);
169 } else if(e < (xSections[0])->GetMaxEnergy()) {
170 idx = (xSections[0])->ComputeLogVectorBin(loge);
171 } else {
172 ekin = (xSections[0])->GetMaxEnergy();
173 idx = (xSections[0])->GetVectorLength() - 2;
174 }
175 // 2. Do the linear interp.(corner cases are already excluded)
176 const G4double x1 = (xSections[0])->Energy(idx);
177 const G4double a = (ekin - x1)/((xSections[0])->Energy(idx+1) - x1);
178 const G4double urnd = G4UniformRand();
179 for (G4int i = 0; i < nElmMinusOne; ++i) {
180 const G4double y1 = (*xSections[i])[idx];
181 if (urnd <= y1 + a*((*xSections[i])[idx+1] - y1)) {
182 element = (*theElementVector)[i];
183 break;
184 }
185 }
186 }
187 return element;
188}
#define G4UniformRand()
Definition Randomize.hh:52

◆ SelectRandomAtom() [2/2]

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

Definition at line 108 of file G4EmElementSelector.hh.

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

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