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