Geant4-11
G4DNAOneStepThermalizationModel.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: Mathieu Karamitros
28
29// The code is developed in the framework of the ESA AO7146
30//
31// We would be very happy hearing from you, send us your feedback! :)
32//
33// In order for Geant4-DNA to be maintained and still open-source,
34// article citations are crucial.
35// If you use Geant4-DNA chemistry and you publish papers about your software,
36// in addition to the general paper on Geant4-DNA:
37//
38// Int. J. Model. Simul. Sci. Comput. 1 (2010) 157–178
39//
40// we would be very happy if you could please also cite the following
41// reference papers on chemistry:
42//
43// J. Comput. Phys. 274 (2014) 841-882
44// Prog. Nucl. Sci. Tec. 2 (2011) 503-508
45
46#ifndef G4DNAOneStepThermalizationModel_hh
47#define G4DNAOneStepThermalizationModel_hh
48
49#include <memory>
50#include "G4VEmModel.hh"
51
52class G4ITNavigator;
53class G4Navigator;
54
55namespace DNA{
56 namespace Penetration{
57 //-----------------------
58 /*
59 * Article: Jintana Meesungnoen, Jean-Paul Jay-Gerin,
60 * Abdelali Filali-Mouhim, and Samlee Mankhetkorn (2002)
61 * Low-Energy Electron Penetration Range in Liquid Water.
62 * Radiation Research: November 2002, Vol. 158, No. 5, pp.657-660.
63 */
65 static void GetPenetration(G4double energy,
66 G4ThreeVector& displacement);
67 static double GetRmean(double energy);
68 //-----
69 // Polynomial fit of Meesungnoen, 2002
70 static const double gCoeff[13];
71 };
72
74 static void GetPenetration(G4double energy,
75 G4ThreeVector& displacement);
76 static double GetRmean(double energy);
77 //-----
78 // Polynomial fit of Meesungnoen, 2002
79 static const double gCoeff[7];
80 };
81
82 //-----------------------
83 /*
84 * Article: Kreipl M S, Friedland W, Paretzke H G (2009) Time- and
85 * space-resolved Monte Carlo study of water radiolysis
86 * for photon, electron and ion irradiation.
87 * Radiat Environ Biophys 48:11-20
88 */
89
90 struct Kreipl2009{
91 static void GetPenetration(G4double energy,
92 G4ThreeVector& displacement);
93 };
94
95 //-----------------------
96 /*
97 * Article: Terrissol M, Beaudre A (1990) Simulation of space and time
98 * evolution of radiolytic species induced by electrons in water.
99 * Radiat Prot Dosimetry 31:171–175
100 */
102 static void GetPenetration(G4double energy,
103 G4ThreeVector& displacement);
104 static double GetRmean(double energy);
105 static double Get3DStdDeviation(double energy);
106 //-----
107 // Terrisol, 1990
108 static const double gEnergies_T1990[11];
109 static const double gStdDev_T1990[11];
110 };
111
112 //-----------------------
113 /*
114 * Article: Ritchie RH, Hamm RN, Turner JE, Bolch WE (1994) Interaction of
115 * low-energy electrons with condensed matter: relevance for track
116 * structure.
117 * Computational approaches in molecular radiation biology, Plenum,
118 * New York, Vol. 63, pp. 155–166
119 * Note: also used in Ballarini et al., 2000
120 */
122 static void GetPenetration(G4double energy,
123 G4ThreeVector& displacement);
124 static double GetRmean(double energy);
125 };
126 }
127}
128
136template<typename MODEL=DNA::Penetration::Meesungnoen2002>
138{
139public:
140 typedef MODEL Model;
142 const G4String& nam =
143 "DNAOneStepThermalizationModel");
145
146 virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
147
149 const G4ParticleDefinition* p,
150 G4double ekin,
151 G4double emin,
152 G4double emax);
153
154 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
156 const G4DynamicParticle*,
157 G4double tmin,
158 G4double maxEnergy);
159
160 inline void SetVerbose(int flag){
161 fVerboseLevel = flag;
162 }
163
165 G4ThreeVector& displacement);
166
167 double GetRmean(double energy);
168
169protected:
170 const std::vector<G4double>* fpWaterDensity;
171
175 std::unique_ptr<G4Navigator> fpNavigator;
176
177private:
181};
182
184
186
187// typedef G4TDNAOneStepThermalizationModel<DNA::Penetration::Terrisol1990> G4DNAOneStepThermalizationModel;
188// Note: if you use the above distribution, it would be
189// better to follow the electrons down to 6 eV and only then apply
190// the one step thermalization
191
193{
194public:
197 static G4VEmModel* Create(const G4String& penetrationModel);
198
206};
207
208#endif
static const G4double emax
G4TDNAOneStepThermalizationModel< DNA::Penetration::Meesungnoen2002 > G4DNAOneStepThermalizationModel
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
static G4VEmModel * GetMacroDefinedModel()
One step thermalization model can be chosen via macro using /process/dna/e-SolvationSubType Ritchie19...
static G4VEmModel * Create(const G4String &penetrationModel)
const std::vector< G4double > * fpWaterDensity
void GetPenetration(G4double energy, G4ThreeVector &displacement)
G4TDNAOneStepThermalizationModel & operator=(const G4TDNAOneStepThermalizationModel &right)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
G4TDNAOneStepThermalizationModel(const G4TDNAOneStepThermalizationModel &)
G4TDNAOneStepThermalizationModel(const G4ParticleDefinition *p=0, const G4String &nam="DNAOneStepThermalizationModel")
G4double energy(const ThreeVector &p, const G4double m)
string material
Definition: eplot.py:19
static void GetPenetration(G4double energy, G4ThreeVector &displacement)
static void GetPenetration(G4double energy, G4ThreeVector &displacement)
static void GetPenetration(G4double energy, G4ThreeVector &displacement)
static void GetPenetration(G4double energy, G4ThreeVector &displacement)
static double GetRmean(double energy)
static void GetPenetration(G4double energy, G4ThreeVector &displacement)
static double Get3DStdDeviation(double energy)