Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4PenelopeIonisationCrossSection.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 // Author: Luciano Pandola
29 //
30 // History:
31 // --------
32 // 14 Mar 2012 L Pandola First complete implementation for e-
33 //
34 //
37 //
39 #include "G4SystemOfUnits.hh"
43 #include "G4Electron.hh"
45 
47  G4VhShellCrossSection("Penelope"),shellIDTable(0),
48  theCrossSectionHandler(0)
49 {
51  nMaxLevels = 9;
52 
53  // Verbosity scale:
54  // 0 = nothing
55  // 1 = calculation of cross sections, file openings, sampling of atoms
56  // 2 = entering in methods
57  verboseLevel = 0;
58 
59  fLowEnergyLimit = 10.0*eV;
60  fHighEnergyLimit = 100.0*GeV;
61 
62  transitionManager = G4AtomicTransitionManager::Instance();
63 }
64 
65 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
66 
68 {
69  if (theCrossSectionHandler)
70  delete theCrossSectionHandler;
71 }
72 
73 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
74 
77  G4double incidentEnergy,
78  G4double ,
79  const G4Material* material)
80 {
81  if (verboseLevel > 1)
82  G4cout << "Entering in method G4PenelopeIonisationCrossSection::CrossSection()" << G4endl;
83 
84  G4double cross = 0.;
85 
86  //Material pointer is not available
87  if (!material)
88  {
89  //CRASH!
91  ed << "The method has been called with a NULL G4Material pointer" << G4endl;
92  G4Exception("G4PenelopeIonisationCrossSection::CrossSection()","em2042",
93  FatalException,ed);
94  return cross;
95  }
96 
97  if (!theCrossSectionHandler)
98  theCrossSectionHandler = new G4PenelopeIonisationXSHandler();
99 
100  G4int nmax = std::min(nMaxLevels,transitionManager->NumberOfShells(Z));
101 
102  if(G4int(shell) < nmax &&
103  incidentEnergy >= fLowEnergyLimit && incidentEnergy <= fHighEnergyLimit)
104  {
105  //The shells in Penelope are organized per *material*, rather than per
106  //element, so given a material one has to find the proper index for the
107  //given Z and shellID. An appropriate lookup table is used to avoid
108  //recalculation.
109  G4int index = FindShellIDIndex(material,Z,shell);
110 
111  //Index is not available!
112  if (index < 0)
113  return cross;
114 
115  G4PenelopeCrossSection* theXS =
116  theCrossSectionHandler->GetCrossSectionTableForCouple(G4Electron::Electron(),
117  material,
118  0.);
119 
120  //Cross check that everything is fine:
121  G4PenelopeOscillator* theOsc = oscManager->GetOscillatorIonisation(material,index);
122  if (theOsc->GetParentZ() != Z || theOsc->GetShellFlag()-1 != G4int(shell))
123  {
124  //something went wrong!
126  ed << "There is something wrong here: it looks like the index is wrong" << G4endl;
127  ed << "Requested: shell " << G4int(shell) << " and Z = " << Z << G4endl;
128  ed << "Retrieved: " << theOsc->GetShellFlag()-1 << " and Z = " << theOsc->GetParentZ() << G4endl;
129  G4Exception("G4PenelopeIonisationCrossSection::CrossSection()","em2043",
130  JustWarning,ed);
131  return cross;
132  }
133 
134 
135  G4double crossPerMolecule = (theXS) ? theXS->GetShellCrossSection(index,incidentEnergy) : 0.;
136 
137  //Now it must be converted in cross section per atom. I need the number of
138  //atoms of the given Z per molecule.
139  G4double atomsPerMolec = oscManager->GetNumberOfZAtomsPerMolecule(material,Z);
140  if (atomsPerMolec)
141  cross = crossPerMolecule/atomsPerMolec;
142 
143 
144  if (verboseLevel > 0)
145  {
146  G4cout << "Cross section of shell " << G4int(shell) << " and Z= " << Z;
147  G4cout << " of material: " << material->GetName() << " and energy = " << incidentEnergy/keV << " keV" << G4endl;
148  G4cout << "--> " << cross/barn << " barn" << G4endl;
149  G4cout << "Shell binding energy: " << theOsc->GetIonisationEnergy()/eV << " eV;" ;
150  G4cout << " resonance energy: " << theOsc->GetResonanceEnergy()/eV << "eV" << G4endl;
151  if (verboseLevel > 2)
152  {
153  G4cout << "Cross section per molecule: " << crossPerMolecule/barn << " barn" << G4endl;
154  G4cout << "Atoms " << Z << " per molecule: " << atomsPerMolec << G4endl;
155  }
156  }
157  }
158 
159  return cross;
160 }
161 
162 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
163 std::vector<G4double>
165  G4double kinEnergy,
166  G4double, G4double,
167  const G4Material* mat)
168 {
169  G4int nmax = std::min(nMaxLevels,transitionManager->NumberOfShells(Z));
170  std::vector<G4double> vec(nmax,0.0);
171  for(G4int i=0; i<nmax; ++i) {
172  vec[i] = CrossSection(Z, G4AtomicShellEnumerator(i), kinEnergy,0.,mat);
173  }
174  return vec;
175 }
176 
177 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
178 
179 std::vector<G4double>
181  G4double kinEnergy,
182  G4double,
183  G4double,
184  const G4Material* mat)
185 {
186  std::vector<G4double> vec = GetCrossSection(Z, kinEnergy,0,0,mat);
187  size_t n = vec.size();
188  size_t i=0;
189  G4double sum = 0.0;
190  for(i=0; i<n; ++i) { sum += vec[i]; }
191  if(sum > 0.0) {
192  sum = 1.0/sum;
193  for(i=0; i<n; ++i) { vec[i] = vec[i]*sum; }
194  }
195  return vec;
196 }
197 
198 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
199 
200 G4int G4PenelopeIonisationCrossSection::FindShellIDIndex(const G4Material* mat,
201  G4int Z,
203 {
204  if (verboseLevel > 1)
205  G4cout << "Entering in method G4PenelopeIonisationCrossSection::FindShellIDIndex()" << G4endl;
206 
207  if (!shellIDTable)
208  shellIDTable = new std::map< std::pair<const G4Material*,G4int>, G4DataVector*>;
209 
210  std::pair<const G4Material*,G4int> theKey = std::make_pair(mat,Z);
211  G4int result = -1;
212  G4int ishell = G4int(shell);
213 
214  if (shellIDTable->count(theKey)) //table already built, and containing the element
215  {
216  if (verboseLevel > 2)
217  G4cout << "FindShellIDIndex: Table already built for " << mat->GetName() << G4endl;
218  G4DataVector* theVec = shellIDTable->find(theKey)->second;
219 
220  if (ishell>=0 && ishell < (G4int) theVec->size()) //check we are not off-boundary
221  result = (G4int) (*theVec)[ishell];
222  else
223  {
225  ed << "Shell ID: " << ishell << " not available for material " << mat->GetName() << " and Z = " <<
226  Z << G4endl;
227  G4Exception("G4PenelopeIonisationCrossSection::FindShellIDIndex()","em2041",JustWarning,
228  ed);
229  return -1;
230  }
231  }
232  else
233  {
234  if (verboseLevel > 2)
235  G4cout << "FindShellIDIndex: Table to be built for " << mat->GetName() << G4endl;
236  //Not contained: look for it
237  G4PenelopeOscillatorTable* theTable = oscManager->GetOscillatorTableIonisation(mat);
238  size_t numberOfOscillators = theTable->size();
239 
240  //oscillator loop
241  //initialize everything at -1
242  G4DataVector* dat = new G4DataVector(nMaxLevels,-1);
243  for (size_t iosc=0;iosc<numberOfOscillators;iosc++)
244  {
245  G4PenelopeOscillator* theOsc = (*theTable)[iosc];
246  //level is found!
247  if (theOsc->GetParentZ() == Z)
248  {
249  //individual shells relative to the given material
250  G4int shFlag = theOsc->GetShellFlag();
251  //Notice: GetShellFlag() starts from 1, the G4AtomicShellEnumerator from 0
252  if (shFlag < 30)
253  (*dat)[shFlag-1] = (G4double) iosc; //index of the given shell
254  if ((shFlag-1) == ishell) // this is what we were looking for
255  result = (G4int) iosc;
256  }
257  }
258  shellIDTable->insert(std::make_pair(theKey,dat));
259  }
260 
261  if (verboseLevel > 1)
262  G4cout << "Leaving method G4PenelopeIonisationCrossSection::FindShellIDIndex() with index = " << result << G4endl;
263 
264  return result;
265 
266 }