Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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"
41 #include "G4ProductionCutsTable.hh"
42 
44 #include "G4Electron.hh"
45 #include "G4Proton.hh"
47 
48 #include "G4LogLogInterpolation.hh"
49 
51 #include "G4VAtomDeexcitation.hh"
52 #include "G4NistManager.hh"
53 
55 {
56 
57 public:
58 
60  const G4String& nam =
61  "DNAEmfietzoglouIonisationModel");
62 
64 
65  virtual void Initialise(const G4ParticleDefinition*,
66  const G4DataVector& = *(new G4DataVector()));
67 
68  virtual G4double CrossSectionPerVolume(const G4Material* material,
69  const G4ParticleDefinition* p,
70  G4double ekin,
71  G4double emin,
72  G4double emax);
73 
74  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
75  const G4MaterialCutsCouple*,
76  const G4DynamicParticle*,
77  G4double tmin,
78  G4double maxEnergy);
79 
80  double DifferentialCrossSection(G4ParticleDefinition * aParticleDefinition,
81  G4double k,
82  G4double energyTransfer,
83  G4int shell);
84 
85  inline void SelectFasterComputation(G4bool input);
86 
87  inline void SelectStationary(G4bool input);
88 
89 protected:
90 
92 
93 private:
94 
95  G4bool fasterCode;
96 
97  G4bool statCode;
98 
99  // Water density table
100  const std::vector<G4double>* fpMolWaterDensity;
101 
102  // Deexcitation manager to produce fluo photons and e-
103  G4VAtomDeexcitation* fAtomDeexcitation;
104 
105  std::map<G4String, G4double, std::less<G4String> > lowEnergyLimit;
106  std::map<G4String, G4double, std::less<G4String> > highEnergyLimit;
107 
108  // TODO :
109 // std::map<const G4ParticleDefinition*,std::pair<G4double,G4double> > fEnergyLimits;
110 
111  G4bool isInitialised;
112  G4int verboseLevel;
113 
114  // Cross section
115 
116  typedef std::map<G4String, G4String, std::less<G4String> > MapFile;
117  MapFile tableFile; // useful ?
118 
119  typedef std::map<G4String, G4DNACrossSectionDataSet*, std::less<G4String> > MapData;
120  MapData tableData;
121 
122  // Final state
123 
125 
126  G4double RandomizeEjectedElectronEnergy(G4ParticleDefinition * aParticleDefinition,
127  G4double incomingParticleEnergy,
128  G4int shell);
129 
130  G4double RandomizeEjectedElectronEnergyFromCumulatedDcs(G4ParticleDefinition * aParticleDefinition,
131  G4double incomingParticleEnergy,
132  G4int shell);
133 
134  G4double RandomTransferedEnergy(G4ParticleDefinition * aParticleDefinition,
135  G4double incomingParticleEnergy,
136  G4int shell);
137 
138  G4double Interpolate(G4double e1,
139  G4double e2,
140  G4double e,
141  G4double xs1,
142  G4double xs2);
143 
144  G4double QuadInterpolator(G4double e11,
145  G4double e12,
146  G4double e21,
147  G4double e22,
148  G4double x11,
149  G4double x12,
150  G4double x21,
151  G4double x22,
152  G4double t1,
153  G4double t2,
154  G4double t,
155  G4double e);
156 
157  typedef std::map<double, std::map<double, double> > TriDimensionMap;
158 
159  TriDimensionMap eDiffCrossSectionData[6];
160  TriDimensionMap eNrjTransfData[6]; // for cumulated dcs
161 
162  TriDimensionMap pDiffCrossSectionData[6];
163 
164  std::vector<double> eTdummyVec;
165 
166  typedef std::map<double, std::vector<double> > VecMap;
167 
168  VecMap eVecm;
169 
170  VecMap eProbaShellMap[6]; // for cumulated dcs
171 
172  // Partial cross section
173 
174  G4int RandomSelect(G4double energy, const G4String& particle);
175 
176  //
177 
180 
181 };
182 
184 {
185  fasterCode = input;
186 }
187 
188 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
189 
191 {
192  statCode = input;
193 }
194 
195 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
196 
197 #endif
G4DNAEmfietzoglouIonisationModel(const G4ParticleDefinition *p=0, const G4String &nam="DNAEmfietzoglouIonisationModel")
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
const char * p
Definition: xmltok.h:285
int G4int
Definition: G4Types.hh:78
double DifferentialCrossSection(G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
bool G4bool
Definition: G4Types.hh:79
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &=*(new G4DataVector()))
static const G4double emax
G4double energy(const ThreeVector &p, const G4double m)
double G4double
Definition: G4Types.hh:76