Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4PenelopeIonisationXSHandler.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// $Id$
27//
28// Author: Luciano Pandola
29//
30// History:
31// -----------
32// 09 Mar 2012 L. Pandola 1st implementation.
33//
34// -------------------------------------------------------------------
35//
36//! Class description:
37//! This class is meant to calculate, store and provide the shell-per-shell
38//! and the total ionisation cross section calculated by the Penelope model.
39//! The information is provided to other physics models, notably: Penelope
40//! ionisation model and Penelope "PIXE" model. Cross sections are calculated
41//! per material-cut couple and stored as G4PenelopeCrossSection objects.
42//!
43// -------------------------------------------------------------------
44
45#ifndef G4PENELOPEIONISATIONXSHANDLER_HH
46#define G4PENELOPEIONISATIONXSHANDLER_HH 1
47
48#include "globals.hh"
49#include "G4DataVector.hh"
50#include "G4Material.hh"
51#include <map>
52
59
61{
62
63public:
64 //! Constructor. nBins is the number of intervals in the
65 //! energy grid. By default the energy grid goes from 100 eV
66 //! to 100 GeV.
67 G4PenelopeIonisationXSHandler(size_t nBins=200);
68
69 //! Destructor. Clean all tables.
71
72 //!Returns the density coeection for the material at the given energy
74
75 //! Returns the table of cross sections for the given particle, given
76 //! material and given cut as a G4PenelopeCrossSection* pointer.
78 const G4Material*,G4double cut);
79 //!Setter for the verbosity level
80 void SetVerboseLevel(G4int vl){verboseLevel = vl;};
81
82private:
85
86
87 void BuildXSTable(const G4Material*,G4double cut,
89
90
91 void BuildDeltaTable(const G4Material*);
92
93 G4DataVector* ComputeShellCrossSectionsElectron(G4PenelopeOscillator* ,
94 G4double energy,G4double cut,
95 G4double delta);
96
97 G4DataVector* ComputeShellCrossSectionsPositron(G4PenelopeOscillator* ,
98 G4double energy,G4double cut,
99 G4double delta);
100
101 //Oscillator manager
102 G4PenelopeOscillatorManager* oscManager;
103
104 //G4PenelopeCrossSection takes care of the logs
105 std::map< std::pair<const G4Material*,G4double>, G4PenelopeCrossSection*> *XSTableElectron;
106 std::map< std::pair<const G4Material*,G4double>, G4PenelopeCrossSection*> *XSTablePositron;
107
108 //delta vs. log(energy)
109 std::map<const G4Material*,G4PhysicsFreeVector*> *theDeltaTable;
110
111 //energy grid
112 G4PhysicsLogVector* energyGrid;
113 size_t nBins;
114
115 G4int verboseLevel;
116
117};
118
119#endif
120
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
G4double GetDensityCorrection(const G4Material *, G4double energy)
Returns the density coeection for the material at the given energy.
virtual ~G4PenelopeIonisationXSHandler()
Destructor. Clean all tables.
G4PenelopeCrossSection * GetCrossSectionTableForCouple(const G4ParticleDefinition *, const G4Material *, G4double cut)
void SetVerboseLevel(G4int vl)
Setter for the verbosity level.