Geant4-11
G4DNABornIonisationModel1.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
28#ifndef G4DNABornIonisationModel1_h
29#define G4DNABornIonisationModel1_h 1
30
31#include "G4VEmModel.hh"
34
36#include "G4Electron.hh"
37#include "G4Proton.hh"
39
41
44#include "G4NistManager.hh"
45
46
48{
49
50public:
51
53 const G4String& nam = "DNABornIonisationModel");
54
56
57 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector& = *(new G4DataVector()));
58
60 const G4ParticleDefinition* p,
61 G4double ekin,
62 G4double emin,
64
65 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
67 const G4DynamicParticle*,
68 G4double tmin,
69 G4double maxEnergy);
70
72 G4int /*level*/,
74 G4double /*kineticEnergy*/);
75
76 G4double DifferentialCrossSection(G4ParticleDefinition * aParticleDefinition, G4double k, G4double energyTransfer, G4int shell);
77
79 G4double incomingParticleEnergy, G4int shell, G4double random) ;
80
81 inline void SelectFasterComputation(G4bool input);
82
83 inline void SelectStationary(G4bool input);
84
85 inline void SelectSPScaling(G4bool input);
86
87protected:
88
90
91private:
92
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
106 // TODO :
107// std::map<const G4ParticleDefinition*,std::pair<G4double,G4double> > fEnergyLimits;
108
109
112
113 // Cross section
114
115 typedef std::map<G4String,G4String,std::less<G4String> > MapFile;
116 MapFile tableFile; // useful ?
117
118 typedef std::map<G4String,G4DNACrossSectionDataSet*,std::less<G4String> > MapData;
120
121 // Final state
122
124
125 G4double RandomizeEjectedElectronEnergy(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4int shell) ;
126
127 G4double RandomizeEjectedElectronEnergyFromCumulatedDcs(G4ParticleDefinition * aParticleDefinition, G4double incomingParticleEnergy, G4int shell) ;
128
130
132 G4double e12,
133 G4double e21,
134 G4double e22,
135 G4double x11,
136 G4double x12,
137 G4double x21,
138 G4double x22,
139 G4double t1,
140 G4double t2,
141 G4double t,
142 G4double e);
143
144 typedef std::map<G4double, std::map<G4double, G4double> > TriDimensionMap;
145
147 TriDimensionMap eNrjTransfData[6]; // for cumulated dcs
148
150 TriDimensionMap pNrjTransfData[6]; // for cumulated dcs
151
152 std::vector<G4double> eTdummyVec;
153 std::vector<G4double> pTdummyVec;
154
155 typedef std::map<G4double, std::vector<G4double> > VecMap;
156
159
160 VecMap eProbaShellMap[6]; // for cumulated dcs
161 VecMap pProbaShellMap[6]; // for cumulated dcs
162
163 // Partial cross section
164
165 G4int RandomSelect(G4double energy,const G4String& particle );
166
167 //
168
171
172};
173
175{
176 fasterCode = input;
177}
178
179//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
180
182{
183 statCode = input;
184}
185
186//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
187
189{
190 spScaling = input;
191}
192
193//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
194
195#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, G4double, std::less< G4String > > highEnergyLimit
G4DNABornIonisationModel1(const G4DNABornIonisationModel1 &)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &= *(new G4DataVector()))
G4VAtomDeexcitation * fAtomDeexcitation
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
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)
std::map< G4String, G4double, std::less< G4String > > lowEnergyLimit
G4DNABornIonisationModel1 & operator=(const G4DNABornIonisationModel1 &right)
virtual G4double GetPartialCrossSection(const G4Material *, G4int, const G4ParticleDefinition *, G4double)
const std::vector< G4double > * fpMolWaterDensity
std::map< G4String, G4DNACrossSectionDataSet *, std::less< G4String > > MapData
std::vector< G4double > eTdummyVec
G4double DifferentialCrossSection(G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell)
std::vector< G4double > pTdummyVec
std::map< G4double, std::vector< G4double > > VecMap
std::map< G4double, std::map< G4double, G4double > > TriDimensionMap
G4double RandomizeEjectedElectronEnergy(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell)
G4double Interpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
G4DNAWaterIonisationStructure waterStructure
std::map< G4String, G4String, std::less< G4String > > MapFile
G4double RandomizeEjectedElectronEnergyFromCumulatedDcs(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4DNABornIonisationModel1(const G4ParticleDefinition *p=0, const G4String &nam="DNABornIonisationModel")
G4double TransferedEnergy(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell, G4double random)
G4ParticleChangeForGamma * fParticleChangeForGamma
G4int RandomSelect(G4double energy, const G4String &particle)
G4double energy(const ThreeVector &p, const G4double m)
string material
Definition: eplot.py:19