Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LivermorePhotoElectricModel.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: Sebastien Incerti
28// 30 October 2008
29// on base of G4LowEnergyPhotoElectric developed by A.Forti and M.G.Pia
30//
31// 15 Mar 2010 L. Pandola, removed methods to set explicitly fluorescence cuts.
32// Main cuts from G4ProductionCutsTable are always used
33// 30 May 2011 A Mantero & V Ivanchenko Migration to model design for deexcitation
34// 22 Oct 2012 A & V Ivanchenko Migration data structure to G4PhysicsVector
35// 1 June 2017 M Bandieramonte
36
37
38#ifndef G4LivermorePhotoElectricModel_h
39#define G4LivermorePhotoElectricModel_h 1
40
41#include "G4VEmModel.hh"
42#include "G4ElementData.hh"
43#include "G4Threading.hh"
44#include <vector>
45
49
51{
52public:
53 explicit G4LivermorePhotoElectricModel(const G4String& nam = "LivermorePhElectric");
54
56
57 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
58
61 G4double energy,
62 G4double cutEnergy = 0.0,
63 G4double maxEnergy = DBL_MAX) override;
64
67 G4double energy,
68 G4double Z,
69 G4double A=0,
70 G4double cut=0,
71 G4double emax=DBL_MAX) override;
72
73 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
75 const G4DynamicParticle*,
76 G4double tmin,
77 G4double maxEnergy) override;
78
79
80 void InitialiseForElement(const G4ParticleDefinition*, G4int Z) override;
81
82 inline void SetLimitNumberOfShells(G4int);
84
86 (const G4LivermorePhotoElectricModel &right) = delete;
88
89protected:
91
92private:
93 void ReadData(G4int Z);
94 const G4String& FindDirectoryPath();
95
96 const G4ParticleDefinition* theGamma;
97 const G4ParticleDefinition* theElectron;
98
99 static G4PhysicsFreeVector* fCrossSection[101]; // 101 because Z range is 1-100
100 static G4PhysicsFreeVector* fCrossSectionLE[101];
101 static std::vector<G4double>* fParamHigh[101];
102 static std::vector<G4double>* fParamLow[101];
103 static G4int fNShells[101];
104 static G4int fNShellsUsed[101];
105 static G4ElementData* fShellCrossSection;
106 static G4Material* fWater;
107 static G4double fWaterEnergyLimit;
108 static G4String fDataDirectory;
109
110#ifdef G4MULTITHREADED
111 static G4Mutex livPhotoeffMutex;
112#endif
113
114 G4VAtomDeexcitation* fAtomDeexcitation;
115 std::vector<G4double> fSandiaCof;
116
117 G4double fCurrSection;
118 G4int verboseLevel;
119 G4int maxZ;
120 G4int nShellLimit;
121 G4bool fDeexcitationActive;
122 G4bool isInitialised;
123
124};
125
126inline
128{
129 nShellLimit = n;
130}
131
132//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
133
134#endif
std::mutex G4Mutex
Definition: G4Threading.hh:81
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double energy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX) override
G4ParticleChangeForGamma * fParticleChange
G4LivermorePhotoElectricModel(const G4LivermorePhotoElectricModel &)=delete
G4double GetBindingEnergy(G4int Z, G4int shell)
void InitialiseForElement(const G4ParticleDefinition *, G4int Z) override
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double energy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) override
#define DBL_MAX
Definition: templates.hh:62