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