Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4DNAOneStepThermalizationModel.hh
이 파일의 문서화 페이지로 가기
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: G4DNAOneStepThermalizationModel.hh 110873 2018-06-22 13:11:22Z gcosmo $
27 //
28 // Author: Mathieu Karamitros
29 
30 // The code is developed in the framework of the ESA AO7146
31 //
32 // We would be very happy hearing from you, send us your feedback! :)
33 //
34 // In order for Geant4-DNA to be maintained and still open-source,
35 // article citations are crucial.
36 // If you use Geant4-DNA chemistry and you publish papers about your software,
37 // in addition to the general paper on Geant4-DNA:
38 //
39 // Int. J. Model. Simul. Sci. Comput. 1 (2010) 157–178
40 //
41 // we would be very happy if you could please also cite the following
42 // reference papers on chemistry:
43 //
44 // J. Comput. Phys. 274 (2014) 841-882
45 // Prog. Nucl. Sci. Tec. 2 (2011) 503-508
46 
47 #ifndef G4DNAOneStepThermalizationModel_hh
48 #define G4DNAOneStepThermalizationModel_hh
49 
50 #include <memory>
51 #include "G4VEmModel.hh"
52 
53 class G4ITNavigator;
54 class G4Navigator;
55 
56 namespace DNA{
57  namespace Penetration{
58  //-----------------------
59  /*
60  * Article: Jintana Meesungnoen, Jean-Paul Jay-Gerin,
61  * Abdelali Filali-Mouhim, and Samlee Mankhetkorn (2002)
62  * Low-Energy Electron Penetration Range in Liquid Water.
63  * Radiation Research: November 2002, Vol. 158, No. 5, pp.657-660.
64  */
66  static void GetPenetration(G4double energy,
67  G4ThreeVector& displacement);
68  static double GetRmean(double energy);
69  //-----
70  // Polynomial fit of Meesungnoen, 2002
71  static const double gCoeff[13];
72  };
73 
74  //-----------------------
75  /*
76  * Article: Terrissol M, Beaudre A (1990) Simulation of space and time
77  * evolution of radiolytic species induced by electrons in water.
78  * Radiat Prot Dosimetry 31:171–175
79  */
80  struct Terrisol1990{
81  static void GetPenetration(G4double energy,
82  G4ThreeVector& displacement);
83  static double GetRmean(double energy);
84  static double Get3DStdDeviation(double energy);
85  //-----
86  // Terrisol, 1990
87  static const double gEnergies_T1990[11];
88  static const double gStdDev_T1990[11];
89  };
90 
91  //-----------------------
92  /*
93  * Article: Ritchie RH, Hamm RN, Turner JE, Bolch WE (1994) Interaction of
94  * low-energy electrons with condensed matter: relevance for track
95  * structure.
96  * Computational approaches in molecular radiation biology, Plenum,
97  * New York, Vol. 63, pp. 155–166
98  * Note: also used in Ballarini et al., 2000
99  */
100  struct Ritchie1994{
101  static void GetPenetration(G4double energy,
102  G4ThreeVector& displacement);
103  static double GetRmean(double energy);
104  };
105  }
106 }
107 
115 template<typename MODEL=DNA::Penetration::Meesungnoen2002>
117 {
118 public:
119  typedef MODEL Model;
121  const G4String& nam =
122  "DNAOneStepThermalizationModel");
124 
125  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
126 
127  virtual G4double CrossSectionPerVolume(const G4Material* material,
128  const G4ParticleDefinition* p,
129  G4double ekin,
130  G4double emin,
131  G4double emax);
132 
133  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
134  const G4MaterialCutsCouple*,
135  const G4DynamicParticle*,
136  G4double tmin,
137  G4double maxEnergy);
138 
139  inline void SetVerbose(int flag){
140  fVerboseLevel = flag;
141  }
142 
144  G4ThreeVector& displacement);
145 
146  double GetRmean(double energy);
147 
148 protected:
149  const std::vector<G4double>* fpWaterDensity;
150 
154  std::unique_ptr<G4Navigator> fpNavigator;
155 
156 private:
160 };
161 
162 #include "G4DNAOneStepThermalizationModel.hpp"
163 
165 
166 // typedef G4TDNAOneStepThermalizationModel<DNA::Penetration::Terrisol1990> G4DNAOneStepThermalizationModel;
167 // Note: if you use the above distribution, it would be
168 // better to follow the electrons down to 6 eV and only then apply
169 // the one step thermalization
170 
172 {
173 public:
176  static G4VEmModel* Create(const G4String& penetrationModel);
177 
185 };
186 
187 #endif
static double Get3DStdDeviation(double energy)
G4TDNAOneStepThermalizationModel< DNA::Penetration::Meesungnoen2002 > G4DNAOneStepThermalizationModel
const char * p
Definition: xmltok.h:285
static const G4double emax
double GetRmean(double energy)
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
const std::vector< G4double > * fpWaterDensity
double energy
Definition: plottest35.C:25
void GetPenetration(G4double energy, G4ThreeVector &displacement)
static void GetPenetration(G4double energy, G4ThreeVector &displacement)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4TDNAOneStepThermalizationModel(const G4ParticleDefinition *p=0, const G4String &nam="DNAOneStepThermalizationModel")
static G4VEmModel * Create(const G4String &penetrationModel)
int G4int
Definition: G4Types.hh:78
static void GetPenetration(G4double energy, G4ThreeVector &displacement)
static double GetRmean(double energy)
static G4VEmModel * GetMacroDefinedModel()
One step thermalization model can be chosen via macro using /process/dna/e-SolvationSubType Ritchie19...
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4TDNAOneStepThermalizationModel & operator=(const G4TDNAOneStepThermalizationModel &right)
static void GetPenetration(G4double energy, G4ThreeVector &displacement)