Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4LivermoreIonisationCrossSection.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// Author: Vladimir Ivanchenko
28//
29// History:
30// --------
31// 31 May 2011 V.Ivanchenko Created
32// 09 Mar 2012 L.Pandola update methods
33//
34//
35
37#include "G4SystemOfUnits.hh"
42
43//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
44
45
47 const G4String& nam) : G4VhShellCrossSection(nam), crossSectionHandler(nullptr)
48{
49 fLowEnergyLimit = 10.0*eV;
50 fHighEnergyLimit = 100.0*GeV;
51
52 transitionManager = G4AtomicTransitionManager::Instance();
53 verboseLevel = 0;
54
55 Initialise();
56}
57
58//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
59
61{
62 delete crossSectionHandler;
63}
64
65//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
66
68{
69 const G4int binForFluo = 20;
70 G4int nbin = G4int(std::log10(fHighEnergyLimit/fLowEnergyLimit) + 0.5);
71 if(nbin <= 0) { nbin = 1; }
72 nbin *= binForFluo;
73
74 // Data on shell ionisation x-sections
75 if (crossSectionHandler) {
76 crossSectionHandler->Clear();
77 delete crossSectionHandler;
78 }
79
81 crossSectionHandler =
82 new G4eCrossSectionHandler(inter,fLowEnergyLimit,fHighEnergyLimit,nbin);
83 crossSectionHandler->LoadShellData("ioni/ion-ss-cs-");
84 //G4cout << "!!! G4LivermoreIonisationCrossSection::Initialise()" << G4endl;
85}
86
89 G4double kinEnergy, G4double,
90 const G4Material*)
91{
92 G4double cross = 0.0;
93 G4int n = G4int(shell);
94 G4int nmax = std::min(9,transitionManager->NumberOfShells(Z));
95 if(Z > 6 && Z < 93 && n < nmax &&
96 kinEnergy >= fLowEnergyLimit && kinEnergy <= fHighEnergyLimit) {
97 cross = crossSectionHandler->FindValue(Z, kinEnergy, n);
98 }
99 return cross;
100}
101
102//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
103
104std::vector<G4double>
106 G4double kinEnergy,
108 const G4Material*)
109{
110 G4int nmax = std::min(9,transitionManager->NumberOfShells(Z));
111 std::vector<G4double> vec(nmax,0.0);
112 for(G4int i=0; i<nmax; ++i) {
113 vec[i] = CrossSection(Z, G4AtomicShellEnumerator(i), kinEnergy);
114 }
115 return vec;
116}
117
118//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
119
120std::vector<G4double>
122 G4double kinEnergy,
123 G4double,
124 G4double,
125 const G4Material*)
126{
127 std::vector<G4double> vec = GetCrossSection(Z, kinEnergy);
128 size_t n = vec.size();
129 size_t i;
130 G4double sum = 0.0;
131 for(i=0; i<n; ++i) { sum += vec[i]; }
132 if(sum > 0.0) {
133 sum = 1.0/sum;
134 for(i=0; i<n; ++i) { vec[i] = vec[i]*sum; }
135 }
136 return vec;
137}
138
139//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4AtomicShellEnumerator
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
static G4AtomicTransitionManager * Instance()
std::vector< G4double > GetCrossSection(G4int Z, G4double incidentEnergy, G4double mass=0.0, G4double deltaEnergy=0.0, const G4Material *mat=0) override
G4double CrossSection(G4int Z, G4AtomicShellEnumerator shell, G4double incidentEnergy, G4double mass=0.0, const G4Material *mat=0) override
G4LivermoreIonisationCrossSection(const G4String &nam="LivermorePIXE")
std::vector< G4double > Probabilities(G4int Z, G4double incidentEnergy, G4double mass=0.0, G4double deltaEnergy=0, const G4Material *mat=0) override
void LoadShellData(const G4String &dataFile)
G4double FindValue(G4int Z, G4double e) const