Geant4-11
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
50{
51
52public:
59 G4VDNAModel(const G4String& nam, const G4String& applyToMaterial);
60
64 virtual ~G4VDNAModel();
65
72 virtual void Initialise(const G4ParticleDefinition* particle,
73 const G4DataVector& cuts,
74 G4ParticleChangeForGamma* fpChangeForGamme=nullptr) =0;
75
76
91 const G4String& materialName,
92 const G4ParticleDefinition* p,
93 G4double ekin,
94 G4double emin,
95 G4double emax) = 0;
96
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
120 G4bool IsMaterialDefine(const G4String &materialName);
121
128 G4bool IsMaterialExistingInModel(const G4String &materialName);
129
139 G4bool IsParticleExistingInModelForMaterial(const G4String &particleName, const G4String &materialName);
140
146
153 G4double GetHighELimit(const G4String& material, const G4String& particle) {return fHighEnergyLimits[material][particle];}
154
161 G4double GetLowELimit(const G4String& material, const G4String& particle) {return fLowEnergyLimits[material][particle];}
162
169 void SetHighELimit(const G4String& material, const G4String& particle, G4double lim) {fHighEnergyLimits[material][particle]=lim;}
170
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 //
194
195 // Setters
196 // ... no setters
197
204 std::vector<G4String> BuildApplyToMatVect(const G4String &materials);
205
214 void ReadAndSaveCSFile(const G4String &materialName, const G4String &particleName, const G4String &file, G4double scaleFactor);
215
225 G4int RandomSelectShell(G4double k, const G4String &particle, const G4String &materialName);
226
236 void AddCrossSectionData(G4String materialName, G4String particleName, G4String fileCS, G4String fileDiffCS, G4double scaleFactor);
237
247 void AddCrossSectionData(G4String materialName, G4String particleName, G4String fileCS, G4double scaleFactor);
248
253 void LoadCrossSectionData(const G4String &particleName);
254
264 virtual void ReadDiffCSFile(const G4String& materialName,
265 const G4String& particleName,
266 const G4String& path,
267 const G4double scaleFactor);
268
277 void EnableForMaterialAndParticle(const G4String& materialName, const G4String& particleName);
278
279private:
287
293
294 std::vector<G4String> fModelMaterials;
295 std::vector<G4String> fModelParticles;
296 std::vector<G4String> fModelCSFiles;
297 std::vector<G4String> fModelDiffCSFiles;
298 std::vector<G4double> fModelScaleFactors;
299
300 std::map<G4String, std::map<G4String, G4double> > fLowEnergyLimits;
301 std::map<G4String, std::map<G4String, G4double> > fHighEnergyLimits;
302
304};
305
306#endif // G4VDNAModel_HH
static const G4double emax
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
std::vector< G4String > fModelDiffCSFiles
List the differential corss section data files.
Definition: G4VDNAModel.hh:297
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
G4String fName
model name
Definition: G4VDNAModel.hh:303
G4VDNAModel(const G4String &nam, const G4String &applyToMaterial)
G4VDNAModel Constructeur of the G4VDNAModel class.
Definition: G4VDNAModel.cc:35
const G4String fStringOfMaterials
fStringOfMaterials The user can decide to specify by hand which are the materials the be activated am...
Definition: G4VDNAModel.hh:286
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
TableMapData fTableData
fTableData It contains the cross section data and can be used like: dataTable=fTableData[material][pa...
Definition: G4VDNAModel.hh:292
std::vector< G4double > fModelScaleFactors
List the model scale factors (they could change with material)
Definition: G4VDNAModel.hh:298
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, std::map< G4String, G4double > > fHighEnergyLimits
List the high energy limits.
Definition: G4VDNAModel.hh:301
std::vector< G4String > fModelCSFiles
List the cross section data files.
Definition: G4VDNAModel.hh:296
std::vector< G4String > fModelMaterials
List the materials that can be activated (and will be by default) within the model.
Definition: G4VDNAModel.hh:294
std::vector< G4String > fModelParticles
List the particles that can be activated within the model.
Definition: G4VDNAModel.hh:295
std::map< G4String, G4double >::const_iterator ItCompoMapData
Definition: G4VDNAModel.hh:185
std::map< G4String, std::map< G4String, G4double > > fLowEnergyLimits
List the low energy limits.
Definition: G4VDNAModel.hh:300
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
string material
Definition: eplot.py:19
#define DBL_MAX
Definition: templates.hh:62