Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PenelopeIonisationCrossSection.hh
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: Luciano Pandola
28//
29// History:
30// -----------
31// 14 Mar 2012 L. Pandola 1st implementation.
32//
33// -------------------------------------------------------------------
34//
35//! Class description:
36//! Access to Penelope ionisation cross sections for e- and e+. To be
37//! registered as a specific model in G4UAtomicDeexcitation.
38//!
39//! NOTICE: working only for e- at the moment (no interface available for
40//! e+)
41//!
42// -------------------------------------------------------------------
43
44#ifndef G4PENELOPEIONISATIONCROSSSECTION_HH
45#define G4PENELOPEIONISATIONCROSSSECTION_HH 1
46
47#include "globals.hh"
49#include <map>
50
56
58{
59public:
60 //! Constructor.
62
63 //! Destructor. Clean all tables.
65
66 //! Purely virtual method from the base interface. Returns the cross
67 //! section for all levels of element Z in material mat at the
68 //! given energy
69 std::vector<G4double> GetCrossSection(G4int Z,
70 G4double incidentEnergy,
71 G4double mass,
72 G4double deltaEnergy,
73 const G4Material* mat);
74
75 //! Purely virtual method from the base interface. Returns the
76 //! cross section for the given shell in the element Z of material
77 //! mat at the specified energy
80 G4double incidentEnergy,
81 G4double mass,
82 const G4Material* mat);
83
84 //! Purely virtual method from the base interface. Returns the
85 //! shell ionisation probabilities for the given Z in the
86 //! material mat at the specified energy.
87 std::vector<G4double> Probabilities(G4int Z,
88 G4double incidentEnergy,
89 G4double mass,
90 G4double deltaEnergy,
91 const G4Material* mat) ;
92 //! Getter/setter for the verbosity level
93 void SetVerbosityLevel(G4int vl){verboseLevel = vl;};
94 G4int GetVerbosityLevel(){return verboseLevel;};
95
96private:
99
100 //Oscillator manager
101 G4PenelopeOscillatorManager* oscManager;
102
103 G4int verboseLevel;
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 look-up table is used to avoid
108 //!recalculation.
109 G4int FindShellIDIndex(const G4Material* mat,G4int Z,G4AtomicShellEnumerator shell);
110 std::map< std::pair<const G4Material*,G4int>, G4DataVector*> *shellIDTable;
111
112 G4int nMaxLevels;
113
114 G4double fLowEnergyLimit;
115 G4double fHighEnergyLimit;
116
117 G4PenelopeIonisationXSHandler* theCrossSectionHandler;
118 const G4AtomicTransitionManager* transitionManager;
119};
120
121#endif
122
G4AtomicShellEnumerator
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
void SetVerbosityLevel(G4int vl)
Getter/setter for the verbosity level.
~G4PenelopeIonisationCrossSection()
Destructor. Clean all tables.
std::vector< G4double > GetCrossSection(G4int Z, G4double incidentEnergy, G4double mass, G4double deltaEnergy, const G4Material *mat)
std::vector< G4double > Probabilities(G4int Z, G4double incidentEnergy, G4double mass, G4double deltaEnergy, const G4Material *mat)
G4double CrossSection(G4int Z, G4AtomicShellEnumerator shell, G4double incidentEnergy, G4double mass, const G4Material *mat)