Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VDNAModel Class Referenceabstract

The G4VDNAModel class. More...

#include <G4VDNAModel.hh>

+ Inheritance diagram for G4VDNAModel:

Public Member Functions

 G4VDNAModel (const G4String &nam, const G4String &applyToMaterial)
 G4VDNAModel Constructeur of the G4VDNAModel class.
 
virtual ~G4VDNAModel ()
 ~G4VDNAModel
 
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. It is used by the process to determine the step path and must return a cross section times a number of molecules per volume unit.
 
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 after the ModelInterface or if any charateristic of the incident particle will change.
 
G4bool IsMaterialDefine (const G4String &materialName)
 IsMaterialDefine Check if the given material is defined in the simulation.
 
G4bool IsMaterialExistingInModel (const G4String &materialName)
 IsMaterialExistingInModel Check if the given material is defined in the current model class.
 
G4bool IsParticleExistingInModelForMaterial (const G4String &particleName, const G4String &materialName)
 IsParticleExistingInModelForMaterial To check two things: 1- is the material existing in model ? 2- if yes, is the particle defined for that material ?
 
G4String GetName ()
 GetName.
 
G4double GetHighELimit (const G4String &material, const G4String &particle)
 GetHighEnergyLimit.
 
G4double GetLowELimit (const G4String &material, const G4String &particle)
 GetLowEnergyLimit.
 
void SetHighELimit (const G4String &material, const G4String &particle, G4double lim)
 SetHighEnergyLimit.
 
void SetLowELimit (const G4String &material, const G4String &particle, G4double lim)
 SetLowEnergyLimit.
 

Protected Types

typedef std::map< G4String, std::map< G4String, G4DNACrossSectionDataSet *, std::less< G4String > > > TableMapData
 
typedef std::map< G4String, std::map< G4String, G4double > > RatioMapData
 
typedef std::map< G4String, G4double >::const_iterator ItCompoMapData
 

Protected Member Functions

TableMapDataGetTableData ()
 GetTableData.
 
std::vector< G4StringBuildApplyToMatVect (const G4String &materials)
 BuildApplyToMatVect Build the material name vector which is used to know the materials the user want to include in the model.
 
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->loadData()
 
G4int RandomSelectShell (G4double k, const G4String &particle, const G4String &materialName)
 RandomSelectShell Method to randomely select a shell from the data table uploaded. The size of the table (number of columns) is used to determine the total number of possible shells.
 
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. It adds a material to the model and fills vectors with informations.
 
void AddCrossSectionData (G4String materialName, G4String particleName, G4String fileCS, G4double scaleFactor)
 AddCrossSectionData 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. Not every model needs differential cross sections.
 
void LoadCrossSectionData (const G4String &particleName)
 LoadCrossSectionData Method to loop on all the registered materials in the model and load the corresponding data.
 
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 sections. The read method for that kind of information is not standardized yet.
 
void EnableForMaterialAndParticle (const G4String &materialName, const G4String &particleName)
 EnableMaterialAndParticle.
 

Detailed Description

The G4VDNAModel class.

All the models using the DNA material management should inherit from that class. The goal is to allow the use of the material management system with little code interferences within the model classes.

Definition at line 49 of file G4VDNAModel.hh.

Member Typedef Documentation

◆ ItCompoMapData

typedef std::map<G4String,G4double>::const_iterator G4VDNAModel::ItCompoMapData
protected

Definition at line 185 of file G4VDNAModel.hh.

◆ RatioMapData

typedef std::map<G4String,std::map<G4String, G4double> > G4VDNAModel::RatioMapData
protected

Definition at line 184 of file G4VDNAModel.hh.

◆ TableMapData

typedef std::map<G4String, std::map<G4String,G4DNACrossSectionDataSet*,std::less<G4String> > > G4VDNAModel::TableMapData
protected

Definition at line 183 of file G4VDNAModel.hh.

Constructor & Destructor Documentation

◆ G4VDNAModel()

G4VDNAModel::G4VDNAModel ( const G4String nam,
const G4String applyToMaterial 
)

G4VDNAModel Constructeur of the G4VDNAModel class.

Parameters
nam
applyToMaterial

Definition at line 35 of file G4VDNAModel.cc.

36 : fStringOfMaterials(applyToMaterial), fName(nam)
37{
38
39}

◆ ~G4VDNAModel()

G4VDNAModel::~G4VDNAModel ( )
virtual

~G4VDNAModel

Definition at line 41 of file G4VDNAModel.cc.

42{
43 // Clean fTableData
44 std::map<G4String, std::map<G4String,G4DNACrossSectionDataSet*,std::less<G4String> > >::iterator posOuter;
45 std::map<G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator posInner;
46 // iterate on each material
47 for (posOuter = fTableData.begin(); posOuter != fTableData.end(); ++posOuter)
48 {
49 // iterate on each particle
50 for(posInner = posOuter->second.begin(); posInner != posOuter->second.end(); ++posInner)
51 {
52 G4DNACrossSectionDataSet* table = posInner->second;
53 if(table != 0) delete table;
54 }
55 }
56}

Member Function Documentation

◆ AddCrossSectionData() [1/2]

void G4VDNAModel::AddCrossSectionData ( G4String  materialName,
G4String  particleName,
G4String  fileCS,
G4double  scaleFactor 
)
protected

AddCrossSectionData 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. Not every model needs differential cross sections.

Parameters
materialName
particleName
fileCS
scaleFactor

Definition at line 67 of file G4VDNAModel.cc.

68{
69 fModelMaterials.push_back(materialName);
70 fModelParticles.push_back(particleName);
71 fModelCSFiles.push_back(fileCS);
72 fModelScaleFactors.push_back(scaleFactor);
73}

◆ AddCrossSectionData() [2/2]

void G4VDNAModel::AddCrossSectionData ( G4String  materialName,
G4String  particleName,
G4String  fileCS,
G4String  fileDiffCS,
G4double  scaleFactor 
)
protected

AddCrossSectionData 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.

Parameters
materialName
particleName
fileCS
fileDiffCS
scaleFactor

Definition at line 58 of file G4VDNAModel.cc.

59{
60 fModelMaterials.push_back(materialName);
61 fModelParticles.push_back(particleName);
62 fModelCSFiles.push_back(fileCS);
63 fModelDiffCSFiles.push_back(fileDiffCS);
64 fModelScaleFactors.push_back(scaleFactor);
65}

Referenced by G4DNAPTBElasticModel::Initialise(), G4DNAPTBExcitationModel::Initialise(), and G4DNAPTBIonisationModel::Initialise().

◆ BuildApplyToMatVect()

std::vector< G4String > G4VDNAModel::BuildApplyToMatVect ( const G4String materials)
protected

BuildApplyToMatVect Build the material name vector which is used to know the materials the user want to include in the model.

Parameters
materials
Returns
a vector with all the material names

Definition at line 139 of file G4VDNAModel.cc.

140{
141 // output material vector
142 std::vector<G4String> materialVect;
143
144 // if we don't find any "/" then it means we only have one "material" (could be the "all" option)
145 if(materials.find("/")==std::string::npos)
146 {
147 // we add the material to the output vector
148 materialVect.push_back(materials);
149 }
150 // if we have several materials listed in the string then we must retrieve them
151 else
152 {
153 G4String materialsNonIdentified = materials;
154
155 while(materialsNonIdentified.find_first_of("/") != std::string::npos)
156 {
157 // we select the first material and stop at the "/" caracter
158 G4String mat = materialsNonIdentified.substr(0, materialsNonIdentified.find_first_of("/"));
159 materialVect.push_back(mat);
160
161 // we remove the previous material from the materialsNonIdentified string
162 materialsNonIdentified = materialsNonIdentified.substr(materialsNonIdentified.find_first_of("/")+1,
163 materialsNonIdentified.size()-materialsNonIdentified.find_first_of("/"));
164 }
165
166 // we don't find "/" anymore, it means we only have one material string left
167 // we get it
168 materialVect.push_back(materialsNonIdentified);
169 }
170
171 return materialVect;
172}

Referenced by LoadCrossSectionData().

◆ CrossSectionPerVolume()

virtual G4double G4VDNAModel::CrossSectionPerVolume ( const G4Material material,
const G4String materialName,
const G4ParticleDefinition p,
G4double  ekin,
G4double  emin,
G4double  emax 
)
pure virtual

CrossSectionPerVolume Every model must implement its own CrossSectionPerVolume method. It is used by the process to determine the step path and must return a cross section times a number of molecules per volume unit.

Parameters
material
materialName
p
ekin
emin
emax
Returns
crossSection*numberOfMoleculesPerVolumeUnit

Implemented in G4DNADummyModel, G4DNAPTBElasticModel, G4DNAPTBExcitationModel, G4DNAPTBIonisationModel, and G4DNAVacuumModel.

Referenced by G4DNAModelInterface::CrossSectionPerVolume().

◆ EnableForMaterialAndParticle()

void G4VDNAModel::EnableForMaterialAndParticle ( const G4String materialName,
const G4String particleName 
)
protected

EnableMaterialAndParticle.

Parameters
materialName
particleNameMeant 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 are registered here. The data should obviously be registered somewhere in the child class. This method is here to allow an easy use of the no-ModelInterface dna models within the ModelInterface system.

Definition at line 134 of file G4VDNAModel.cc.

135{
136 fTableData[materialName][particleName] = 0;
137}

Referenced by G4DNAVacuumModel::Initialise(), and G4DNADummyModel::Initialise().

◆ GetHighELimit()

G4double G4VDNAModel::GetHighELimit ( const G4String material,
const G4String particle 
)
inline

GetHighEnergyLimit.

Parameters
material
particle
Returns
fHighEnergyLimits[material][particle]

Definition at line 153 of file G4VDNAModel.hh.

153{return fHighEnergyLimits[material][particle];}

Referenced by G4DNAPTBElasticModel::CrossSectionPerVolume(), G4DNAPTBExcitationModel::CrossSectionPerVolume(), G4DNAPTBIonisationModel::CrossSectionPerVolume(), G4DNAPTBElasticModel::SampleSecondaries(), G4DNAPTBExcitationModel::SampleSecondaries(), and G4DNAPTBIonisationModel::SampleSecondaries().

◆ GetLowELimit()

G4double G4VDNAModel::GetLowELimit ( const G4String material,
const G4String particle 
)
inline

GetLowEnergyLimit.

Parameters
material
particle
Returns
fLowEnergyLimits[material][particle]

Definition at line 161 of file G4VDNAModel.hh.

161{return fLowEnergyLimits[material][particle];}

Referenced by G4DNAPTBElasticModel::CrossSectionPerVolume(), G4DNAPTBExcitationModel::CrossSectionPerVolume(), G4DNAPTBIonisationModel::CrossSectionPerVolume(), G4DNAPTBElasticModel::SampleSecondaries(), G4DNAPTBExcitationModel::SampleSecondaries(), and G4DNAPTBIonisationModel::SampleSecondaries().

◆ GetName()

G4String G4VDNAModel::GetName ( )
inline

GetName.

Returns
the name of the model

Definition at line 145 of file G4VDNAModel.hh.

145{return fName;}

Referenced by IsMaterialDefine().

◆ GetTableData()

TableMapData * G4VDNAModel::GetTableData ( )
inlineprotected

GetTableData.

Returns
a pointer to a map with the following structure: [materialName][particleName]=G4DNACrossSectionDataSet*

Definition at line 193 of file G4VDNAModel.hh.

193{return &fTableData;}

Referenced by G4DNAPTBElasticModel::CrossSectionPerVolume(), G4DNAPTBExcitationModel::CrossSectionPerVolume(), G4DNAPTBIonisationModel::CrossSectionPerVolume(), and RandomSelectShell().

◆ Initialise()

virtual void G4VDNAModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector cuts,
G4ParticleChangeForGamma fpChangeForGamme = nullptr 
)
pure virtual

Initialise Each model must implement an Initialize method.

Parameters
particle
cuts

Implemented in G4DNAVacuumModel, G4DNAPTBElasticModel, G4DNADummyModel, G4DNAPTBExcitationModel, and G4DNAPTBIonisationModel.

◆ IsMaterialDefine()

G4bool G4VDNAModel::IsMaterialDefine ( const G4String materialName)

IsMaterialDefine Check if the given material is defined in the simulation.

Parameters
materialName
Returns
true if the material is defined in the simulation

Definition at line 237 of file G4VDNAModel.cc.

238{
239 // Check if the given material is defined in the simulation
240
241 G4bool exist (false);
242
243 double matTableSize = G4Material::GetMaterialTable()->size();
244
245 for(int i=0;i<matTableSize;i++)
246 {
247 if(materialName == G4Material::GetMaterialTable()->at(i)->GetName())
248 {
249 exist = true;
250 return exist;
251 }
252 }
253
254 return exist;
255}
bool G4bool
Definition: G4Types.hh:86
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:677
G4String GetName()
GetName.
Definition: G4VDNAModel.hh:145

◆ IsMaterialExistingInModel()

G4bool G4VDNAModel::IsMaterialExistingInModel ( const G4String materialName)

IsMaterialExistingInModel Check if the given material is defined in the current model class.

Parameters
materialName
Returns
true if the material is defined in the model

Definition at line 257 of file G4VDNAModel.cc.

258{
259 // Check if the given material is defined in the current model class
260
261 if (fTableData.find(materialName) == fTableData.end())
262 {
263 return false;
264 }
265 else
266 {
267 return true;
268 }
269}

Referenced by IsParticleExistingInModelForMaterial().

◆ IsParticleExistingInModelForMaterial()

G4bool G4VDNAModel::IsParticleExistingInModelForMaterial ( const G4String particleName,
const G4String materialName 
)

IsParticleExistingInModelForMaterial To check two things: 1- is the material existing in model ? 2- if yes, is the particle defined for that material ?

Parameters
particleName
materialName
Returns
true if the particle/material couple is defined in the model

Definition at line 271 of file G4VDNAModel.cc.

272{
273 // To check two things:
274 // 1- is the material existing in model ?
275 // 2- if yes, is the particle defined for that material ?
276
277 if(IsMaterialExistingInModel(materialName))
278 {
279 if (fTableData[materialName].find(particleName) == fTableData[materialName].end())
280 {
281 return false;
282 }
283 else return true;
284 }
285 else return false;
286}
G4bool IsMaterialExistingInModel(const G4String &materialName)
IsMaterialExistingInModel Check if the given material is defined in the current model class.
Definition: G4VDNAModel.cc:257

◆ LoadCrossSectionData()

void G4VDNAModel::LoadCrossSectionData ( const G4String particleName)
protected

LoadCrossSectionData Method to loop on all the registered materials in the model and load the corresponding data.

Definition at line 75 of file G4VDNAModel.cc.

76{
77 G4String fileElectron, fileDiffElectron;
78 G4String materialName, modelParticleName;
79 G4double scaleFactor;
80
81 // construct applyToMatVect with materials specified by the user
82 std::vector<G4String> applyToMatVect = BuildApplyToMatVect(fStringOfMaterials);
83
84 // iterate on each material contained into the fStringOfMaterials variable (through applyToMatVect)
85 for(unsigned int i=0;i<applyToMatVect.size();++i)
86 {
87 // We have selected a material coming from applyToMatVect
88 // We try to find if this material correspond to a model registered material
89 // If it is, then isMatFound becomes true
90 G4bool isMatFound = false;
91
92 // We iterate on each model registered materials to load the CS data
93 // We have to do a for loop because of the "all" option
94 // applyToMatVect[i] == "all" implies applyToMatVect.size()=1 and we want to iterate on all registered materials
95 for(unsigned int j=0;j<fModelMaterials.size();++j)
96 {
97 if(applyToMatVect[i] == fModelMaterials[j] || applyToMatVect[i] == "all")
98 {
99 isMatFound = true;
100 materialName = fModelMaterials[j];
101 modelParticleName = fModelParticles[j];
102 fileElectron = fModelCSFiles[j];
103 if(!fModelDiffCSFiles.empty()) fileDiffElectron = fModelDiffCSFiles[j];
104 scaleFactor = fModelScaleFactors[j];
105
106 ReadAndSaveCSFile(materialName, modelParticleName, fileElectron, scaleFactor);
107
108 if(!fModelDiffCSFiles.empty()) ReadDiffCSFile(materialName, modelParticleName, fileDiffElectron, scaleFactor);
109
110 }
111 }
112
113 // check if we found a correspondance, if not: fatal error
114 if(!isMatFound)
115 {
116 std::ostringstream oss;
117 oss << applyToMatVect[i] << " material was not found. It means the material specified in the UserPhysicsList is not a model material for ";
118 oss << particleName;
119 G4Exception("G4VDNAModel::LoadCrossSectionData","em0003",
120 FatalException, oss.str().c_str());
121 return;
122 }
123 }
124}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
double G4double
Definition: G4Types.hh:83
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
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
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

Referenced by G4DNAPTBElasticModel::Initialise(), G4DNAPTBExcitationModel::Initialise(), and G4DNAPTBIonisationModel::Initialise().

◆ RandomSelectShell()

G4int G4VDNAModel::RandomSelectShell ( G4double  k,
const G4String particle,
const G4String materialName 
)
protected

RandomSelectShell Method to randomely select a shell from the data table uploaded. The size of the table (number of columns) is used to determine the total number of possible shells.

Parameters
k
particle
materialName
Returns
the selected shell

Definition at line 182 of file G4VDNAModel.cc.

183{
184 G4int level = 0;
185
186 TableMapData* tableData = GetTableData();
187
188 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
189 pos = (*tableData)[materialName].find(particle);
190
191 if (pos != (*tableData)[materialName].end())
192 {
193 G4DNACrossSectionDataSet* table = pos->second;
194
195 if (table != 0)
196 {
197 G4double* valuesBuffer = new G4double[table->NumberOfComponents()];
198 const G4int n = (G4int)table->NumberOfComponents();
199 G4int i(n);
200 G4double value = 0.;
201
202 while (i>0)
203 {
204 --i;
205 valuesBuffer[i] = table->GetComponent(i)->FindValue(k);
206 value += valuesBuffer[i];
207 }
208
209 value *= G4UniformRand();
210
211 i = n;
212
213 while (i > 0)
214 {
215 --i;
216
217 if (valuesBuffer[i] > value)
218 {
219 delete[] valuesBuffer;
220 return i;
221 }
222 value -= valuesBuffer[i];
223 }
224
225 if (valuesBuffer) delete[] valuesBuffer;
226
227 }
228 }
229 else
230 {
231 G4Exception("G4VDNAModel::RandomSelectShell","em0002",
232 FatalException,"Model not applicable to particle type.");
233 }
234 return level;
235}
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
virtual size_t NumberOfComponents(void) const
virtual const G4VEMDataSet * GetComponent(G4int componentId) const
TableMapData * GetTableData()
GetTableData.
Definition: G4VDNAModel.hh:193
std::map< G4String, std::map< G4String, G4DNACrossSectionDataSet *, std::less< G4String > > > TableMapData
Definition: G4VDNAModel.hh:183
virtual G4double FindValue(G4double x, G4int componentId=0) const =0

Referenced by G4DNAPTBExcitationModel::SampleSecondaries(), and G4DNAPTBIonisationModel::SampleSecondaries().

◆ ReadAndSaveCSFile()

void G4VDNAModel::ReadAndSaveCSFile ( const G4String materialName,
const G4String particleName,
const G4String file,
G4double  scaleFactor 
)
protected

ReadAndSaveCSFile Read and save a "simple" cross section file : use of G4DNACrossSectionDataSet->loadData()

Parameters
materialName
particleName
file
scaleFactor

Definition at line 174 of file G4VDNAModel.cc.

177{
178 fTableData[materialName][particleName] = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV, scaleFactor);
179 fTableData[materialName][particleName]->LoadData(file);
180}

Referenced by LoadCrossSectionData().

◆ ReadDiffCSFile()

void G4VDNAModel::ReadDiffCSFile ( const G4String materialName,
const G4String particleName,
const G4String path,
const G4double  scaleFactor 
)
protectedvirtual

ReadDiffCSFile Virtual method that need to be implemented if one wish to use the differential cross sections. The read method for that kind of information is not standardized yet.

Parameters
materialName
particleName
path
scaleFactor

Definition at line 126 of file G4VDNAModel.cc.

127{
128 G4String text("ReadDiffCSFile must be implemented in the model class using a differential cross section data file");
129
130 G4Exception("G4VDNAModel::ReadDiffCSFile","em0003",
131 FatalException, text);
132}

Referenced by LoadCrossSectionData().

◆ SampleSecondaries()

virtual void G4VDNAModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple ,
const G4String materialName,
const G4DynamicParticle ,
G4ParticleChangeForGamma particleChangeForGamma,
G4double  tmin = 0,
G4double  tmax = DBL_MAX 
)
pure virtual

SampleSecondaries Each model must implement SampleSecondaries to decide if a particle will be created after the ModelInterface or if any charateristic of the incident particle will change.

Parameters
materialName
particleChangeForGamma
tmin
tmax

Implemented in G4DNADummyModel, G4DNAPTBElasticModel, G4DNAPTBExcitationModel, G4DNAPTBIonisationModel, and G4DNAVacuumModel.

Referenced by G4DNAModelInterface::SampleSecondaries().

◆ SetHighELimit()

void G4VDNAModel::SetHighELimit ( const G4String material,
const G4String particle,
G4double  lim 
)
inline

SetHighEnergyLimit.

Parameters
material
particle
lim

Definition at line 169 of file G4VDNAModel.hh.

169{fHighEnergyLimits[material][particle]=lim;}

Referenced by G4DNAPTBElasticModel::Initialise(), G4DNADummyModel::Initialise(), G4DNAPTBExcitationModel::Initialise(), and G4DNAPTBIonisationModel::Initialise().

◆ SetLowELimit()

void G4VDNAModel::SetLowELimit ( const G4String material,
const G4String particle,
G4double  lim 
)
inline

SetLowEnergyLimit.

Parameters
material
particle
lim

Definition at line 177 of file G4VDNAModel.hh.

177{fLowEnergyLimits[material][particle]=lim;}

Referenced by G4DNAPTBElasticModel::Initialise(), G4DNADummyModel::Initialise(), G4DNAPTBExcitationModel::Initialise(), and G4DNAPTBIonisationModel::Initialise().


The documentation for this class was generated from the following files: