Geant4 10.7.0
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{
52
53public:
54
55 G4LivermorePhotoElectricModel(const G4String& nam = "LivermorePhElectric");
56
58
59 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
60
63 G4double energy,
64 G4double cutEnergy = 0.0,
65 G4double maxEnergy = DBL_MAX);
66
69 G4double energy,
70 G4double Z,
71 G4double A=0,
72 G4double cut=0,
73 G4double emax=DBL_MAX);
74
75 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
77 const G4DynamicParticle*,
78 G4double tmin,
79 G4double maxEnergy);
80
81
82 virtual void InitialiseForElement(const G4ParticleDefinition*, G4int Z);
83
84 inline void SetLimitNumberOfShells(G4int);
86
88 (const G4LivermorePhotoElectricModel &right) = delete;
90
91protected:
92
94
95private:
96
97 void ReadData(G4int Z);
98
99 const G4String& FindDirectoryPath();
100
101 const G4ParticleDefinition* theGamma;
102 const G4ParticleDefinition* theElectron;
103
104 G4int verboseLevel;
105 G4int maxZ;
106 G4int nShellLimit;
107 G4bool fDeexcitationActive;
108 G4bool isInitialised;
109
110 static G4LPhysicsFreeVector* fCrossSection[99];
111 static G4LPhysicsFreeVector* fCrossSectionLE[99];
112 static std::vector<G4double>* fParamHigh[99];
113 static std::vector<G4double>* fParamLow[99];
114 static G4int fNShells[99];
115 static G4int fNShellsUsed[99];
116 static G4ElementData* fShellCrossSection;
117 static G4Material* fWater;
118 static G4double fWaterEnergyLimit;
119 static G4String fDataDirectory;
120
121#ifdef G4MULTITHREADED
122 static G4Mutex livPhotoeffMutex;
123#endif
124
125 G4VAtomDeexcitation* fAtomDeexcitation;
126
127 G4double fCurrSection;
128 std::vector<G4double> fSandiaCof;
129};
130
131inline
133{
134 nShellLimit = n;
135}
136
137//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
138
139#endif
double A(double temperature)
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
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
virtual G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double energy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
G4ParticleChangeForGamma * fParticleChange
G4LivermorePhotoElectricModel(const G4LivermorePhotoElectricModel &)=delete
G4double GetBindingEnergy(G4int Z, G4int shell)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double energy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
#define DBL_MAX
Definition: templates.hh:62