Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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(0)
48{
49 fLowEnergyLimit = 10.0*eV;
50 fHighEnergyLimit = 100.0*GeV;
51
52 transitionManager = G4AtomicTransitionManager::Instance();
53
54 verboseLevel = 0;
55
56 Initialise();
57}
58
59//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
60
62{
63 delete crossSectionHandler;
64}
65
66//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
67
69{
70 const G4int binForFluo = 20;
71 G4int nbin = G4int(std::log10(fHighEnergyLimit/fLowEnergyLimit) + 0.5);
72 if(nbin <= 0) { nbin = 1; }
73 nbin *= binForFluo;
74
75 // Data on shell ionisation x-sections
76 if (crossSectionHandler) {
77 crossSectionHandler->Clear();
78 delete crossSectionHandler;
79 }
80
82 crossSectionHandler =
83 new G4eCrossSectionHandler(inter,fLowEnergyLimit,fHighEnergyLimit,nbin);
84 crossSectionHandler->LoadShellData("ioni/ion-ss-cs-");
85 //G4cout << "!!! G4LivermoreIonisationCrossSection::Initialise()" << G4endl;
86}
87
90 G4double kinEnergy, G4double,
91 const G4Material*)
92{
93 G4double cross = 0.0;
94 G4int n = G4int(shell);
95 G4int nmax = std::min(9,transitionManager->NumberOfShells(Z));
96 if(Z > 6 && Z < 93 && n < nmax &&
97 kinEnergy >= fLowEnergyLimit && kinEnergy <= fHighEnergyLimit) {
98 //G4cout << "Z= " << Z << " n= " << n << " E(MeV)= " << kinEnergy/MeV << G4endl;
99 cross = crossSectionHandler->FindValue(Z, kinEnergy, n);
100 }
101 return cross;
102}
103
104//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
105
106std::vector<G4double>
108 G4double kinEnergy,
110 const G4Material*)
111{
112 G4int nmax = std::min(9,transitionManager->NumberOfShells(Z));
113 std::vector<G4double> vec(nmax,0.0);
114 for(G4int i=0; i<nmax; ++i) {
115 vec[i] = CrossSection(Z, G4AtomicShellEnumerator(i), kinEnergy);
116 }
117 return vec;
118}
119
120//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
121
122std::vector<G4double>
124 G4double kinEnergy,
125 G4double,
126 G4double,
127 const G4Material*)
128{
129 std::vector<G4double> vec = GetCrossSection(Z, kinEnergy);
130 size_t n = vec.size();
131 size_t i;
132 G4double sum = 0.0;
133 for(i=0; i<n; ++i) { sum += vec[i]; }
134 if(sum > 0.0) {
135 sum = 1.0/sum;
136 for(i=0; i<n; ++i) { vec[i] = vec[i]*sum; }
137 }
138 return vec;
139}
140
141//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4AtomicShellEnumerator
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
static G4AtomicTransitionManager * Instance()
G4double CrossSection(G4int Z, G4AtomicShellEnumerator shell, G4double incidentEnergy, G4double mass=0.0, const G4Material *mat=0)
G4LivermoreIonisationCrossSection(const G4String &nam="LivermorePIXE")
std::vector< G4double > Probabilities(G4int Z, G4double incidentEnergy, G4double mass=0.0, G4double deltaEnergy=0, const G4Material *mat=0)
std::vector< G4double > GetCrossSection(G4int Z, G4double incidentEnergy, G4double mass=0.0, G4double deltaEnergy=0.0, const G4Material *mat=0)
void LoadShellData(const G4String &dataFile)
G4double FindValue(G4int Z, G4double e) const