Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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// $Id$
27//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class file
31//
32//
33// File name: G4EmElementSelector
34//
35// Author: Vladimir Ivanchenko
36//
37// Creation date: 29.05.2008
38//
39// Modifications:
40//
41// Class Description:
42//
43// Generic helper class for the random selection of an element
44
45// -------------------------------------------------------------------
46//
47
49#include "G4VEmModel.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 spline):
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 element = (*theElementVector)[0];
67 if(nElmMinusOne > 0) {
68 xSections.reserve(n);
69 G4PhysicsLogVector* v0 = new G4PhysicsLogVector(lowEnergy,highEnergy,nbins);
70 xSections.push_back(v0);
71 v0->SetSpline(spline);
72 for(G4int i=1; i<n; ++i) {
74 xSections.push_back(v);
75 }
76 }
77 //G4cout << "G4EmElementSelector for " << mat->GetName() << " n= " << n << G4endl;
78}
79
80//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
81
83{
84 if(nElmMinusOne > 0) {
85 for(G4int i=0; i<=nElmMinusOne; ++i) {
86 delete xSections[i];
87 }
88 }
89}
90
91//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
92
94 G4double cut)
95{
96 //G4cout << "G4EmElementSelector initialise for " << material->GetName() << G4endl;
97 if(0 == nElmMinusOne || cut == cutEnergy) { return; }
98
99 cutEnergy = cut;
100 //G4cout << "cut(keV)= " << cut/keV << G4endl;
101 G4double cross;
102
103 const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
104
105 // loop over bins
106 for(G4int j=0; j<=nbins; ++j) {
107 G4double e = (xSections[0])->Energy(j);
108 model->SetupForMaterial(part, material, e);
109 cross = 0.0;
110 //G4cout << "j= " << j << " e(MeV)= " << e/MeV << G4endl;
111 for (G4int i=0; i<=nElmMinusOne; ++i) {
112 cross += theAtomNumDensityVector[i]*
113 model->ComputeCrossSectionPerAtom(part, (*theElementVector)[i], e,
114 cutEnergy, e);
115 xSections[i]->PutValue(j, cross);
116 }
117 }
118
119 // xSections start from null, so use probabilities from the next bin
120 if(0.0 == (*xSections[nElmMinusOne])[0]) {
121 for (G4int i=0; i<=nElmMinusOne; ++i) {
122 xSections[i]->PutValue(0, (*xSections[i])[1]);
123 }
124 }
125 // xSections ends with null, so use probabilities from the previous bin
126 if(0.0 == (*xSections[nElmMinusOne])[nbins]) {
127 for (G4int i=0; i<=nElmMinusOne; ++i) {
128 xSections[i]->PutValue(nbins, (*xSections[i])[nbins-1]);
129 }
130 }
131 // perform normalization
132 for(G4int j=0; j<=nbins; ++j) {
133 cross = (*xSections[nElmMinusOne])[j];
134 // only for positive X-section
135 if(cross > 0.0) {
136 for (G4int i=0; i<nElmMinusOne; ++i) {
137 G4double x = (*xSections[i])[j]/cross;
138 xSections[i]->PutValue(j, x);
139 }
140 }
141 }
142 //G4cout << "======== G4EmElementSelector for the " << model->GetName();
143}
144
145//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
146
148{
149 G4cout << "======== G4EmElementSelector for the " << model->GetName();
150 if(part) G4cout << " and " << part->GetParticleName();
151 G4cout << " for " << material->GetName() << " ========" << G4endl;
152 if(0 < nElmMinusOne) {
153 for(G4int i=0; i<nElmMinusOne; i++) {
154 G4cout << " " << (*theElementVector)[i]->GetName() << " : " << G4endl;
155 G4cout << *(xSections[i]) << G4endl;
156 }
157 }
158 G4cout << "Last Element in element vector "
159 << (*theElementVector)[nElmMinusOne]->GetName()
160 << G4endl;
161 G4cout << G4endl;
162}
163
164//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
165
166
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
void Initialise(const G4ParticleDefinition *, G4double cut=0.0)
void Dump(const G4ParticleDefinition *p=0)
G4EmElementSelector(G4VEmModel *, const G4Material *, G4int bins, G4double emin, G4double emax, G4bool spline=true)
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:189
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:205
const G4String & GetName() const
Definition: G4Material.hh:177
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:240
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:311
const G4String & GetName() const
Definition: G4VEmModel.hh:655