Geant4 11.2.2
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#ifndef G4LivermorePhotoElectricModel_h
38#define G4LivermorePhotoElectricModel_h 1
39
40#include "G4ElementData.hh"
41#include "G4VEmModel.hh"
42
43#include <vector>
44
48
50{
51public:
52 explicit G4LivermorePhotoElectricModel(const G4String& nam = "LivermorePhElectric");
53
55
56 void Initialise(const G4ParticleDefinition*, const G4DataVector&) override;
57
59 G4double energy, G4double cutEnergy = 0.0,
60 G4double maxEnergy = DBL_MAX) override;
61
63 G4double Z, G4double A = 0, G4double cut = 0,
64 G4double emax = DBL_MAX) override;
65
66 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
68 G4double tmin, G4double maxEnergy) override;
69
70 void InitialiseForElement(const G4ParticleDefinition*, G4int Z) override;
71
72 [[maybe_unused]] inline void SetLimitNumberOfShells(G4int n) { nShellLimit = n; };
74
76 (const G4LivermorePhotoElectricModel& right) = delete;
78
79protected:
81
82private:
83 void ReadData(const G4int Z);
84
85 const G4String& FindDirectoryPath();
86
87 void InitialiseOnFly(G4int Z);
88
89 const G4ParticleDefinition* theGamma;
90 const G4ParticleDefinition* theElectron;
91 G4VAtomDeexcitation* fAtomDeexcitation{nullptr};
92
93 static constexpr G4int ZMAXPE{101}; // 101 because Z range is 1-100
94 static G4ElementData* fCrossSection;
95 static G4ElementData* fCrossSectionLE;
96 static std::vector<G4double>* fParamHigh[ZMAXPE];
97 static std::vector<G4double>* fParamLow[ZMAXPE];
98 static G4int fNShells[ZMAXPE];
99 static G4int fNShellsUsed[ZMAXPE];
100 static G4Material* fWater;
101 static G4double fWaterEnergyLimit;
102 static G4String fDataDirectory;
103
104 std::vector<G4double> fSandiaCof;
105
106 G4double fCurrSection{0.0};
107 G4int verboseLevel;
108 G4int nShellLimit{100};
109 G4bool fDeexcitationActive{false};
110 G4bool isInitializer{false};
111};
112
113//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
114
115#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
G4double ComputeCrossSectionPerAtom(const G4ParticleDefinition *, G4double energy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX) override
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
G4LivermorePhotoElectricModel(const G4String &nam="LivermorePhElectric")
G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double energy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX) override
#define DBL_MAX
Definition templates.hh:62