Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAPTBElasticModel.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// Authors: S. Meylan and C. Villagrasa (IRSN, France)
27// Models come from
28// M. Bug et al, Rad. Phys and Chem. 130, 459-479 (2017)
29//
30
31
32#ifndef G4DNAPTBElasticModel_h
33#define G4DNAPTBElasticModel_h 1
34
35#include <map>
37#include "G4VDNAModel.hh"
38#include "G4Electron.hh"
42#include "G4NistManager.hh"
43
44/*!
45 * \brief The G4DNAPTBElasticModel class
46 * This class implements the elastic model for the DNA materials and precursors.
47 */
49{
50
51public:
52
53 /*!
54 * \brief G4DNAPTBElasticModel
55 * Constructor
56 * \param applyToMaterial
57 * \param p
58 * \param nam
59 */
60 G4DNAPTBElasticModel(const G4String &applyToMaterial = "all", const G4ParticleDefinition* p = 0,
61 const G4String& nam = "DNAPTBElasticModel");
62
63 /*!
64 * \brief ~G4DNAPTBElasticModel
65 * Destructor
66 */
67 virtual ~G4DNAPTBElasticModel();
68
69 /*!
70 * \brief Initialise
71 * Mandatory method for every model class. The material/particle for which the model
72 * can be used have to be added here through the AddCrossSectionData method.
73 * Then the LoadCrossSectionData method must be called to trigger the load process.
74 * Scale factors to be applied to the cross section can be defined here.
75 */
76 virtual void Initialise(const G4ParticleDefinition* particle, const G4DataVector&, G4ParticleChangeForGamma* fpChangeForGamme=nullptr);
77
78 /*!
79 * \brief CrossSectionPerVolume
80 * This method is mandatory for any model class. It finds and return the cross section value
81 * for the current material, particle and energy values.
82 * The number of molecule per volume is not used here but in the G4DNAModelInterface class.
83 * \param material
84 * \param materialName
85 * \param p
86 * \param ekin
87 * \param emin
88 * \param emax
89 * \return the cross section value
90 */
91 virtual G4double CrossSectionPerVolume(const G4Material* material,
92 const G4String& materialName,
93 const G4ParticleDefinition* p,
94 G4double ekin,
95 G4double emin,
96 G4double emax);
97
98 /*!
99 * \brief SampleSecondaries
100 * Method called after CrossSectionPerVolume if the process is the one which is selected (according to the sampling on the calculated path length).
101 * Here, the characteristics of the incident and created (if any) particle(s) are set (energy, momentum ...).
102 * \param materialName
103 * \param particleChangeForGamma
104 * \param tmin
105 * \param tmax
106 */
107 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
109 const G4String& materialName,
110 const G4DynamicParticle*,
111 G4ParticleChangeForGamma *particleChangeForGamma,
112 G4double tmin,
113 G4double tmax);
114
115protected:
116
117
118
119private:
120
121 G4int verboseLevel; ///< verbose level
122 std::map<G4String, double > killBelowEnergyTable; ///< map to save the different energy kill limits for the materials
123 G4double fKillBelowEnergy; ///< energy kill limit
124
125 typedef std::map<G4String, std::map<G4String, std::map<double, std::map<double, double> > > > TriDimensionMap;
126 TriDimensionMap diffCrossSectionData; ///< A map: [materialName][particleName]=DiffCrossSectionTable
127
128 typedef std::map<G4String, std::map<G4String, std::map<double, std::vector<double> > > > VecMap;
129 VecMap eValuesVect; /*!< map with vectors containing all the output energy (E) of the differential file */
130 std::map<G4String, std::map<G4String, std::vector<double> > > tValuesVec; ///< map with vectors containing all the incident (T) energy of the differential file
131
132 /*!
133 * \brief ReadDiffCSFile
134 * Method to read the differential cross section files. This method is not standard yet so every model must implement its own.
135 * \param materialName
136 * \param particleName
137 * \param file
138 */
139 void ReadDiffCSFile(const G4String &materialName, const G4String &particleName, const G4String &file, const G4double);
140
141 /*!
142 * \brief Theta
143 * To return an angular theta value from the differential file. This method uses interpolations to calculate
144 * the theta value.
145 * \param fParticleDefinition
146 * \param k
147 * \param integrDiff
148 * \param materialName
149 * \return a theta value
150 */
151 G4double Theta(G4ParticleDefinition * fParticleDefinition, G4double k, G4double integrDiff, const G4String &materialName);
152
153 /*!
154 * \brief LinLinInterpolate
155 * \param e1
156 * \param e2
157 * \param e
158 * \param xs1
159 * \param xs2
160 * \return
161 */
162 G4double LinLinInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
163
164 /*!
165 * \brief LinLogInterpolate
166 * \param e1
167 * \param e2
168 * \param e
169 * \param xs1
170 * \param xs2
171 * \return
172 */
173 G4double LinLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
174
175 /*!
176 * \brief LogLogInterpolate
177 * \param e1
178 * \param e2
179 * \param e
180 * \param xs1
181 * \param xs2
182 * \return
183 */
184 G4double LogLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2);
185
186 /*!
187 * \brief QuadInterpolator
188 * \param e11
189 * \param e12
190 * \param e21
191 * \param e22
192 * \param x11
193 * \param x12
194 * \param x21
195 * \param x22
196 * \param t1
197 * \param t2
198 * \param t
199 * \param e
200 * \return
201 */
202 G4double QuadInterpolator(G4double e11,
203 G4double e12,
204 G4double e21,
205 G4double e22,
206 G4double x11,
207 G4double x12,
208 G4double x21,
209 G4double x22,
210 G4double t1,
211 G4double t2,
212 G4double t,
213 G4double e);
214
215 /*!
216 * \brief RandomizeCosTheta
217 * \param k
218 * \param materialName
219 * \return
220 */
221 G4double RandomizeCosTheta(G4double k, const G4String &materialName);
222
223 // copy constructor and hide assignment operator
224 G4DNAPTBElasticModel(G4DNAPTBElasticModel &); // prevent copy-construction
225 G4DNAPTBElasticModel & operator=(const G4DNAPTBElasticModel &right); // prevent assignement
226};
227
228#endif
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
The G4DNAPTBElasticModel class This class implements the elastic model for the DNA materials and prec...
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4String &materialName, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
CrossSectionPerVolume This method is mandatory for any model class. It finds and return the cross sec...
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4String &materialName, const G4DynamicParticle *, G4ParticleChangeForGamma *particleChangeForGamma, G4double tmin, G4double tmax)
SampleSecondaries Method called after CrossSectionPerVolume if the process is the one which is select...
virtual ~G4DNAPTBElasticModel()
~G4DNAPTBElasticModel Destructor
virtual void Initialise(const G4ParticleDefinition *particle, const G4DataVector &, G4ParticleChangeForGamma *fpChangeForGamme=nullptr)
Initialise Mandatory method for every model class. The material/particle for which the model can be u...
The G4VDNAModel class.
Definition: G4VDNAModel.hh:50