Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmElementSelector.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class file
30//
31//
32// File name: G4EmElementSelector
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 29.05.2008
37//
38// Modifications:
39//
40// Class Description:
41//
42// Generic helper class for the random selection of an element
43
44// -------------------------------------------------------------------
45//
46
48#include "G4VEmModel.hh"
49#include "G4SystemOfUnits.hh"
50
51//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
52//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
53
55 const G4Material* mat,
56 G4int bins,
57 G4double emin,
58 G4double emax,
59 G4bool):
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}
81
82//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
83
85{
86 if(nElmMinusOne > 0) {
87 for(G4int i=0; i<=nElmMinusOne; ++i) { delete xSections[i]; }
88 }
89}
90
91//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
92
94 G4double cut)
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}
152
153//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
154
155const G4Element*
157 const G4double loge) const
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}
189
190//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
191
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}
208
209//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
210
211
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
void Initialise(const G4ParticleDefinition *, G4double cut=0.0)
void Dump(const G4ParticleDefinition *p=nullptr)
G4EmElementSelector(G4VEmModel *, const G4Material *, G4int bins, G4double emin, G4double emax, G4bool spline=true)
const G4Element * SelectRandomAtom(const G4double kineticEnergy, const G4double logEKin) const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:185
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:201
const G4String & GetName() const
Definition: G4Material.hh:172
const G4String & GetParticleName() const
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:284
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:383
const G4String & GetName() const
Definition: G4VEmModel.hh:816