Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4VDNAModel.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// This class is used to support PTB models that come from
28// M. Bug et al, Rad. Phys and Chem. 130, 459-479 (2017)
29//
30
31#ifndef G4VDNAModel_HH
32#define G4VDNAModel_HH
33
34#ifdef _MSC_VER
35 #pragma warning(disable : 4503)
36#endif
37
41#include "G4VEmModel.hh"
42
43/*! \class G4VDNAModel
44 * \brief The G4VDNAModel class
45 *
46 * All the models using the DNA material management should inherit from that class.
47 * The goal is to allow the use of the material management system with little code interferences within the model classes.
48 */
50{
51
52public:
53 /*!
54 * \brief G4VDNAModel
55 * Constructeur of the G4VDNAModel class.
56 * \param nam
57 * \param applyToMaterial
58 */
59 G4VDNAModel(const G4String& nam, const G4String& applyToMaterial);
60
61 /*!
62 * \brief ~G4VDNAModel
63 */
64 virtual ~G4VDNAModel();
65
66 /*!
67 * \brief Initialise
68 * Each model must implement an Initialize method.
69 * \param particle
70 * \param cuts
71 */
72 virtual void Initialise(const G4ParticleDefinition* particle,
73 const G4DataVector& cuts,
74 G4ParticleChangeForGamma* fpChangeForGamme=nullptr) =0;
75
76
77 /*!
78 * \brief CrossSectionPerVolume
79 * Every model must implement its own CrossSectionPerVolume method.
80 * It is used by the process to determine the step path and must return a cross section times a number
81 * of molecules per volume unit.
82 * \param material
83 * \param materialName
84 * \param p
85 * \param ekin
86 * \param emin
87 * \param emax
88 * \return crossSection*numberOfMoleculesPerVolumeUnit
89 */
90 virtual G4double CrossSectionPerVolume(const G4Material* material,
91 const G4String& materialName,
92 const G4ParticleDefinition* p,
93 G4double ekin,
94 G4double emin,
95 G4double emax) = 0;
96
97 /*!
98 * \brief SampleSecondaries
99 * Each model must implement SampleSecondaries to decide if a particle will be created after the ModelInterface or
100 * if any charateristic of the incident particle will change.
101 * \param materialName
102 * \param particleChangeForGamma
103 * \param tmin
104 * \param tmax
105 */
106 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
108 const G4String& materialName,
109 const G4DynamicParticle*,
110 G4ParticleChangeForGamma *particleChangeForGamma,
111 G4double tmin = 0,
112 G4double tmax = DBL_MAX) = 0;
113
114 /*!
115 * \brief IsMaterialDefine
116 * Check if the given material is defined in the simulation
117 * \param materialName
118 * \return true if the material is defined in the simulation
119 */
120 G4bool IsMaterialDefine(const G4String &materialName);
121
122 /*!
123 * \brief IsMaterialExistingInModel
124 * Check if the given material is defined in the current model class
125 * \param materialName
126 * \return true if the material is defined in the model
127 */
128 G4bool IsMaterialExistingInModel(const G4String &materialName);
129
130 /*!
131 * \brief IsParticleExistingInModelForMaterial
132 * To check two things:
133 * 1- is the material existing in model ?
134 * 2- if yes, is the particle defined for that material ?
135 * \param particleName
136 * \param materialName
137 * \return true if the particle/material couple is defined in the model
138 */
139 G4bool IsParticleExistingInModelForMaterial(const G4String &particleName, const G4String &materialName);
140
141 /*!
142 * \brief GetName
143 * \return the name of the model
144 */
145 G4String GetName(){return fName;}
146
147 /*!
148 * \brief GetHighEnergyLimit
149 * \param material
150 * \param particle
151 * \return fHighEnergyLimits[material][particle]
152 */
153 G4double GetHighELimit(const G4String& material, const G4String& particle) {return fHighEnergyLimits[material][particle];}
154
155 /*!
156 * \brief GetLowEnergyLimit
157 * \param material
158 * \param particle
159 * \return fLowEnergyLimits[material][particle]
160 */
161 G4double GetLowELimit(const G4String& material, const G4String& particle) {return fLowEnergyLimits[material][particle];}
162
163 /*!
164 * \brief SetHighEnergyLimit
165 * \param material
166 * \param particle
167 * \param lim
168 */
169 void SetHighELimit(const G4String& material, const G4String& particle, G4double lim) {fHighEnergyLimits[material][particle]=lim;}
170
171 /*!
172 * \brief SetLowEnergyLimit
173 * \param material
174 * \param particle
175 * \param lim
176 */
177 void SetLowELimit(const G4String& material, const G4String& particle, G4double lim) {fLowEnergyLimits[material][particle]=lim;}
178
179protected:
180
181 // typedef used to ease the data container reading
182 //
183 typedef std::map<G4String, std::map<G4String,G4DNACrossSectionDataSet*,std::less<G4String> > > TableMapData;
184 typedef std::map<G4String,std::map<G4String, G4double> > RatioMapData;
185 typedef std::map<G4String, G4double>::const_iterator ItCompoMapData;
186
187 // Getters
188 //
189 /*!
190 * \brief GetTableData
191 * \return a pointer to a map with the following structure: [materialName][particleName]=G4DNACrossSectionDataSet*
192 */
193 TableMapData* GetTableData(){return &fTableData;}
194
195 // Setters
196 // ... no setters
197
198 /*!
199 * \brief BuildApplyToMatVect
200 * Build the material name vector which is used to know the materials the user want to include in the model.
201 * \param materials
202 * \return a vector with all the material names
203 */
204 std::vector<G4String> BuildApplyToMatVect(const G4String &materials);
205
206 /*!
207 * \brief ReadAndSaveCSFile
208 * Read and save a "simple" cross section file : use of G4DNACrossSectionDataSet->loadData()
209 * \param materialName
210 * \param particleName
211 * \param file
212 * \param scaleFactor
213 */
214 void ReadAndSaveCSFile(const G4String &materialName, const G4String &particleName, const G4String &file, G4double scaleFactor);
215
216 /*!
217 * \brief RandomSelectShell
218 * Method to randomely select a shell from the data table uploaded.
219 * The size of the table (number of columns) is used to determine the total number of possible shells.
220 * \param k
221 * \param particle
222 * \param materialName
223 * \return the selected shell
224 */
225 G4int RandomSelectShell(G4double k, const G4String &particle, const G4String &materialName);
226
227 /*!
228 * \brief AddCrossSectionData
229 * Method used during the initialization of the model class to add a new material. It adds a material to the model and fills vectors with informations.
230 * \param materialName
231 * \param particleName
232 * \param fileCS
233 * \param fileDiffCS
234 * \param scaleFactor
235 */
236 void AddCrossSectionData(G4String materialName, G4String particleName, G4String fileCS, G4String fileDiffCS, G4double scaleFactor);
237
238 /*!
239 * \brief AddCrossSectionData
240 * Method used during the initialization of the model class to add a new material. It adds a material to the model and fills vectors with informations.
241 * Not every model needs differential cross sections.
242 * \param materialName
243 * \param particleName
244 * \param fileCS
245 * \param scaleFactor
246 */
247 void AddCrossSectionData(G4String materialName, G4String particleName, G4String fileCS, G4double scaleFactor);
248
249 /*!
250 * \brief LoadCrossSectionData
251 * Method to loop on all the registered materials in the model and load the corresponding data.
252 */
253 void LoadCrossSectionData(const G4String &particleName);
254
255 /*!
256 * \brief ReadDiffCSFile
257 * Virtual method that need to be implemented if one wish to use the differential cross sections.
258 * The read method for that kind of information is not standardized yet.
259 * \param materialName
260 * \param particleName
261 * \param path
262 * \param scaleFactor
263 */
264 virtual void ReadDiffCSFile(const G4String& materialName,
265 const G4String& particleName,
266 const G4String& path,
267 const G4double scaleFactor);
268
269 /*!
270 * \brief EnableMaterialAndParticle
271 * \param materialName
272 * \param particleName
273 * Meant to fill fTableData with 0 for the specified material and particle, therefore allowing the ModelInterface class to proceed with the material and particle even if no data
274 * are registered here. The data should obviously be registered somewhere in the child class.
275 * This method is here to allow an easy use of the no-ModelInterface dna models within the ModelInterface system.
276 */
277 void EnableForMaterialAndParticle(const G4String& materialName, const G4String& particleName);
278
279private:
280 /*!
281 * \brief fStringOfMaterials
282 * The user can decide to specify by hand which are the materials the be activated among those implemented in the model.
283 * If the user does then only the specified materials contained in this string variable will be activated.
284 * The string is like: mat1/mat2/mat3/mat4
285 */
286 const G4String fStringOfMaterials;
287
288 /*!
289 * \brief fTableData
290 * It contains the cross section data and can be used like: dataTable=fTableData[material][particle]
291 */
292 TableMapData fTableData;
293
294 std::vector<G4String> fModelMaterials; ///< List the materials that can be activated (and will be by default) within the model.
295 std::vector<G4String> fModelParticles; ///< List the particles that can be activated within the model
296 std::vector<G4String> fModelCSFiles; ///< List the cross section data files
297 std::vector<G4String> fModelDiffCSFiles; ///< List the differential corss section data files
298 std::vector<G4double> fModelScaleFactors; ///< List the model scale factors (they could change with material)
299
300 std::map<G4String, std::map<G4String, G4double> > fLowEnergyLimits; ///< List the low energy limits
301 std::map<G4String, std::map<G4String, G4double> > fHighEnergyLimits; ///< List the high energy limits
302
303 G4String fName; ///< model name
304};
305
306#endif // G4VDNAModel_HH
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
The G4VDNAModel class.
Definition: G4VDNAModel.hh:50
virtual ~G4VDNAModel()
~G4VDNAModel
Definition: G4VDNAModel.cc:41
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4String &materialName, const G4DynamicParticle *, G4ParticleChangeForGamma *particleChangeForGamma, G4double tmin=0, G4double tmax=DBL_MAX)=0
SampleSecondaries Each model must implement SampleSecondaries to decide if a particle will be created...
G4bool IsMaterialExistingInModel(const G4String &materialName)
IsMaterialExistingInModel Check if the given material is defined in the current model class.
Definition: G4VDNAModel.cc:257
G4int RandomSelectShell(G4double k, const G4String &particle, const G4String &materialName)
RandomSelectShell Method to randomely select a shell from the data table uploaded....
Definition: G4VDNAModel.cc:182
G4bool IsParticleExistingInModelForMaterial(const G4String &particleName, const G4String &materialName)
IsParticleExistingInModelForMaterial To check two things: 1- is the material existing in model ?...
Definition: G4VDNAModel.cc:271
virtual void Initialise(const G4ParticleDefinition *particle, const G4DataVector &cuts, G4ParticleChangeForGamma *fpChangeForGamme=nullptr)=0
Initialise Each model must implement an Initialize method.
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4String &materialName, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)=0
CrossSectionPerVolume Every model must implement its own CrossSectionPerVolume method....
TableMapData * GetTableData()
GetTableData.
Definition: G4VDNAModel.hh:193
G4String GetName()
GetName.
Definition: G4VDNAModel.hh:145
void SetHighELimit(const G4String &material, const G4String &particle, G4double lim)
SetHighEnergyLimit.
Definition: G4VDNAModel.hh:169
void ReadAndSaveCSFile(const G4String &materialName, const G4String &particleName, const G4String &file, G4double scaleFactor)
ReadAndSaveCSFile Read and save a "simple" cross section file : use of G4DNACrossSectionDataSet->load...
Definition: G4VDNAModel.cc:174
std::map< G4String, std::map< G4String, G4DNACrossSectionDataSet *, std::less< G4String > > > TableMapData
Definition: G4VDNAModel.hh:183
void AddCrossSectionData(G4String materialName, G4String particleName, G4String fileCS, G4String fileDiffCS, G4double scaleFactor)
AddCrossSectionData Method used during the initialization of the model class to add a new material....
Definition: G4VDNAModel.cc:58
G4bool IsMaterialDefine(const G4String &materialName)
IsMaterialDefine Check if the given material is defined in the simulation.
Definition: G4VDNAModel.cc:237
virtual void ReadDiffCSFile(const G4String &materialName, const G4String &particleName, const G4String &path, const G4double scaleFactor)
ReadDiffCSFile Virtual method that need to be implemented if one wish to use the differential cross s...
Definition: G4VDNAModel.cc:126
void SetLowELimit(const G4String &material, const G4String &particle, G4double lim)
SetLowEnergyLimit.
Definition: G4VDNAModel.hh:177
std::vector< G4String > BuildApplyToMatVect(const G4String &materials)
BuildApplyToMatVect Build the material name vector which is used to know the materials the user want ...
Definition: G4VDNAModel.cc:139
G4double GetHighELimit(const G4String &material, const G4String &particle)
GetHighEnergyLimit.
Definition: G4VDNAModel.hh:153
std::map< G4String, G4double >::const_iterator ItCompoMapData
Definition: G4VDNAModel.hh:185
void LoadCrossSectionData(const G4String &particleName)
LoadCrossSectionData Method to loop on all the registered materials in the model and load the corresp...
Definition: G4VDNAModel.cc:75
std::map< G4String, std::map< G4String, G4double > > RatioMapData
Definition: G4VDNAModel.hh:184
G4double GetLowELimit(const G4String &material, const G4String &particle)
GetLowEnergyLimit.
Definition: G4VDNAModel.hh:161
void EnableForMaterialAndParticle(const G4String &materialName, const G4String &particleName)
EnableMaterialAndParticle.
Definition: G4VDNAModel.cc:134
#define DBL_MAX
Definition: templates.hh:62