Geant4-11
G4DNAEmfietzoglouIonisationModel.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// Based on the work described in
27// Rad Res 163, 98-111 (2005)
28// D. Emfietzoglou, H. Nikjoo
29//
30// Authors of the class (2014):
31// I. Kyriakou (kyriak@cc.uoi.gr)
32// D. Emfietzoglou (demfietz@cc.uoi.gr)
33// S. Incerti (incerti@cenbg.in2p3.fr)
34//
35
36#ifndef G4DNAEmfietzoglouIonisationModel_h
37#define G4DNAEmfietzoglouIonisationModel_h 1
38
39#include "G4VEmModel.hh"
43#include "G4NistManager.hh"
44#include "G4Electron.hh"
45#include "G4Proton.hh"
46
51
53{
54
55public:
56
58 const G4String& nam =
59 "DNAEmfietzoglouIonisationModel");
60
62
63 virtual void Initialise(const G4ParticleDefinition*,
64 const G4DataVector& = *(new G4DataVector()));
65
67 const G4ParticleDefinition* p,
68 G4double ekin,
69 G4double emin,
71
72 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
74 const G4DynamicParticle*,
75 G4double tmin,
76 G4double maxEnergy);
77
79 G4double k,
80 G4double energyTransfer,
81 G4int shell);
82
83 inline void SelectFasterComputation(G4bool input);
84
85 inline void SelectStationary(G4bool input);
86
87protected:
88
90
91private:
92
94
96
97 // Water density table
98 const std::vector<G4double>* fpMolWaterDensity;
99
100 // Deexcitation manager to produce fluo photons and e-
102
103 std::map<G4String, G4double, std::less<G4String> > lowEnergyLimit;
104 std::map<G4String, G4double, std::less<G4String> > highEnergyLimit;
105
108
109 // Cross section
110
111 typedef std::map<G4String, G4String, std::less<G4String> > MapFile;
112 MapFile tableFile; // useful ?
113
114 typedef std::map<G4String, G4DNACrossSectionDataSet*, std::less<G4String> > MapData;
116
117 // Final state
118
120
122 G4double incomingParticleEnergy,
123 G4int shell);
124
126 G4double incomingParticleEnergy,
127 G4int shell);
128
130 G4double incomingParticleEnergy,
131 G4int shell);
132
134 G4double e2,
135 G4double e,
136 G4double xs1,
137 G4double xs2);
138
140 G4double e12,
141 G4double e21,
142 G4double e22,
143 G4double x11,
144 G4double x12,
145 G4double x21,
146 G4double x22,
147 G4double t1,
148 G4double t2,
149 G4double t,
150 G4double e);
151
152 typedef std::map<G4double, std::map<G4double, G4double> > TriDimensionMap;
153
155 TriDimensionMap eNrjTransfData[6]; // for cumulated dcs
156
158
159 std::vector<G4double> eTdummyVec;
160
161 typedef std::map<G4double, std::vector<G4double> > VecMap;
162
164
165 VecMap eProbaShellMap[6]; // for cumulated dcs
166
167 // Partial cross section
168
169 G4int RandomSelect(G4double energy, const G4String& particle);
170
171 //
172
175
176};
177
179{
180 fasterCode = input;
181}
182
183//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
184
186{
187 statCode = input;
188}
189
190//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
191
192#endif
static const G4double e1[44]
static const G4double e2[44]
static const G4double emax
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
std::map< G4String, G4DNACrossSectionDataSet *, std::less< G4String > > MapData
std::map< G4double, std::vector< G4double > > VecMap
G4double QuadInterpolator(G4double e11, G4double e12, G4double e21, G4double e22, G4double x11, G4double x12, G4double x21, G4double x22, G4double t1, G4double t2, G4double t, G4double e)
G4DNAEmfietzoglouIonisationModel(const G4DNAEmfietzoglouIonisationModel &)
G4DNAEmfietzoglouIonisationModel(const G4ParticleDefinition *p=0, const G4String &nam="DNAEmfietzoglouIonisationModel")
G4DNAEmfietzoglouIonisationModel & operator=(const G4DNAEmfietzoglouIonisationModel &right)
G4double RandomizeEjectedElectronEnergyFromCumulatedDcs(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &= *(new G4DataVector()))
std::map< G4String, G4String, std::less< G4String > > MapFile
G4double RandomizeEjectedElectronEnergy(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell)
const std::vector< G4double > * fpMolWaterDensity
std::map< G4String, G4double, std::less< G4String > > highEnergyLimit
std::map< G4double, std::map< G4double, G4double > > TriDimensionMap
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
G4double RandomTransferedEnergy(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell)
std::map< G4String, G4double, std::less< G4String > > lowEnergyLimit
G4double DifferentialCrossSection(G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell)
G4DNAEmfietzoglouWaterIonisationStructure waterStructure
G4double Interpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4int RandomSelect(G4double energy, const G4String &particle)
G4double energy(const ThreeVector &p, const G4double m)
string material
Definition: eplot.py:19