Geant4 10.7.0
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 = 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}
82
83//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
84
86{
87 if(nElmMinusOne > 0) {
88 for(G4int i=0; i<=nElmMinusOne; ++i) { delete xSections[i]; }
89 }
90}
91
92//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
93
95 G4double cut)
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}
148
149//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
150
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}
167
168//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
169
170
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
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 G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
const G4String & GetName() const
Definition: G4Material.hh:175
const G4String & GetParticleName() const
void SetSpline(G4bool)
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
const G4String & GetName() const
Definition: G4VEmModel.hh:827