Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNARelativisticIonisationModel.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: G4DNARelativisticIonisationModel.hh 90057 2015-05-11 22:25:50Z matkara $
27//
28
29#ifndef G4DNARelativisticIonisationModel_h
30#define G4DNARelativisticIonisationModel_h 1
31
32#include "G4VEmModel.hh"
37
39#include "G4Electron.hh"
40#include "G4Proton.hh"
41#include "G4NistManager.hh"
42
44//#include "G4DNAWaterExcitationStructure.hh"
45
47{
48public:
50 const G4String& nam = "DNARelativisticIonisationModel");
51
53
55 =(const G4DNARelativisticIonisationModel &right) = delete;
57
59 const G4DataVector& = *(new G4DataVector())) override;
60
62 const G4ParticleDefinition* p,
63 G4double ekin,
64 G4double emin,
65 G4double emax) override;
66
67 virtual G4double GetTotalCrossSection (const G4Material* material,
69 G4double kineticEnergy);
71 G4int level,
73 G4double kineticEnergy) override;
74 virtual G4double GetDifferentialCrossSection(const G4Material* material,
75 const G4ParticleDefinition* particle,
76 G4double kineticEnergy,
77 G4double secondaryEnergy,
78 G4int level);
79
80 void SampleSecondaries(std::vector<G4DynamicParticle*>*,
82 const G4DynamicParticle*,
83 G4double tmin,
84 G4double maxEnergy) override;
85 virtual void LoadAtomicStates(G4int z, const char *path);
86 inline void SelectStationary (G4bool input){statCode = input;};
87 inline void SelectFasterComputation(G4bool input){fasterCode = input;};
88
89
90protected:
91
93
94private:
95
96 std::vector <G4int > iState [99];
97 std::vector <G4int > iShell [99];
98 std::vector <G4int > iSubShell [99];
99 std::vector <G4double> Nelectrons[99];
100 std::vector <G4double> Ebinding [99];
101 std::vector <G4double> Ekinetic [99];
102
103 std::map <G4int, std::vector<G4double> > eVecEZ;
104
105 using DeauxDimensionVecMapZ = std::map<G4int, std::map<G4double, std::vector<G4double>>>;
106 DeauxDimensionVecMapZ eVecEjeEZ;
107
108 using TriDimensionVecMapZ = std::map<G4int, std::map<G4int, std::map<G4double, std::vector<G4double>>>>;
109 TriDimensionVecMapZ eProbaShellMapZ;
110
111 using QuadDimensionMapZ = std::map<G4int, std::map<G4int, std::map<G4double, std::map<G4double, G4double>>>>;
112 QuadDimensionMapZ eDiffCrossSectionDataZ;
113 QuadDimensionMapZ eEjectedEnergyDataZ;
114
115
116 G4double fLowEnergyLimit=0.;
117 G4double fHighEnergyLimit=0.;
118
119 G4bool isInitialised=false;
120 G4bool statCode=false;
121 G4bool fasterCode=false;
122 G4int verboseLevel=0;
123
124 const std::vector<G4double>* fMaterialDensity=nullptr;
125 const G4ParticleDefinition* fParticleDefinition=nullptr;
126 G4VAtomDeexcitation* fAtomDeexcitation=nullptr;
127
128 G4int RandomSelect(const G4Material* material,
130 G4double kineticEnergy);
131
132 G4double GetEjectedElectronEnergy (
133 const G4Material* material,
134 const G4ParticleDefinition* ,
135 G4double energy,
136 G4int shell );
137 G4ThreeVector GetEjectedElectronDirection(
138 const G4ParticleDefinition* ,
139 G4double energy,G4double secondaryenergy);
140
141 G4double Interpolate (G4double e1 ,G4double e2 ,G4double e ,
142 G4double xs1, G4double xs2);
143 G4double QuadInterpolator(G4double e11,G4double e12,G4double e21,G4double e22,
144 G4double x11,G4double x12,G4double x21,G4double x22,
145 G4double t1 ,G4double t2 ,G4double t ,G4double e);
146
147};
148
149//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
150
151#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4DNARelativisticIonisationModel(const G4DNARelativisticIonisationModel &)=delete
virtual void LoadAtomicStates(G4int z, const char *path)
virtual G4double GetTotalCrossSection(const G4Material *material, const G4ParticleDefinition *, G4double kineticEnergy)
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
void Initialise(const G4ParticleDefinition *, const G4DataVector &= *(new G4DataVector())) override
G4double GetPartialCrossSection(const G4Material *material, G4int level, const G4ParticleDefinition *, G4double kineticEnergy) override
G4DNARelativisticIonisationModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="DNARelativisticIonisationModel")
virtual G4double GetDifferentialCrossSection(const G4Material *material, const G4ParticleDefinition *particle, G4double kineticEnergy, G4double secondaryEnergy, G4int level)
G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override