Geant4-11
Public Member Functions | Protected Types | Protected Member Functions | Private Types | Private Member Functions | Private Attributes
G4DNAPTBExcitationModel Class Reference

The G4DNAPTBExcitationModel class This class implements the PTB excitation model. More...

#include <G4DNAPTBExcitationModel.hh>

Inheritance diagram for G4DNAPTBExcitationModel:
G4VDNAModel

Public Member Functions

virtual G4double CrossSectionPerVolume (const G4Material *material, const G4String &materialName, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
 CrossSectionPerVolume Retrieve the cross section corresponding to the current material, particle and energy. More...
 
 G4DNAPTBExcitationModel (const G4String &applyToMaterial="all", const G4ParticleDefinition *p=0, const G4String &nam="DNAPTBExcitationModel")
 G4DNAPTBExcitationModel Constructor. More...
 
G4double GetHighELimit (const G4String &material, const G4String &particle)
 GetHighEnergyLimit. More...
 
G4double GetLowELimit (const G4String &material, const G4String &particle)
 GetLowEnergyLimit. More...
 
G4String GetName ()
 GetName. More...
 
virtual void Initialise (const G4ParticleDefinition *particle, const G4DataVector &= *(new G4DataVector()), G4ParticleChangeForGamma *fpChangeForGamme=nullptr)
 Initialise Set the materials for which the model can be used and defined the energy limits. More...
 
G4bool IsMaterialDefine (const G4String &materialName)
 IsMaterialDefine Check if the given material is defined in the simulation. More...
 
G4bool IsMaterialExistingInModel (const G4String &materialName)
 IsMaterialExistingInModel Check if the given material is defined in the current model class. More...
 
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 ? More...
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4String &materialName, const G4DynamicParticle *, G4ParticleChangeForGamma *particleChangeForGamma, G4double tmin, G4double tmax)
 SampleSecondaries If the model is selected for the ModelInterface then the SampleSecondaries method will be called. The method sets the incident particle characteristics after the ModelInterface. More...
 
void SetHighELimit (const G4String &material, const G4String &particle, G4double lim)
 SetHighEnergyLimit. More...
 
void SetLowELimit (const G4String &material, const G4String &particle, G4double lim)
 SetLowEnergyLimit. More...
 
virtual ~G4DNAPTBExcitationModel ()
 ~G4DNAPTBExcitationModel Destructor More...
 

Protected Types

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

Protected Member Functions

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. More...
 
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. More...
 
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. More...
 
void EnableForMaterialAndParticle (const G4String &materialName, const G4String &particleName)
 EnableMaterialAndParticle. More...
 
TableMapDataGetTableData ()
 GetTableData. More...
 
void LoadCrossSectionData (const G4String &particleName)
 LoadCrossSectionData Method to loop on all the registered materials in the model and load the corresponding data. More...
 
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. More...
 
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() More...
 
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. More...
 

Private Types

typedef std::map< G4String, G4double, std::less< G4String > > MapMeanEnergy
 

Private Member Functions

 G4DNAPTBExcitationModel (const G4DNAPTBExcitationModel &)
 
G4DNAPTBExcitationModeloperator= (const G4DNAPTBExcitationModel &right)
 

Private Attributes

std::map< G4String, std::map< G4String, G4double > > fHighEnergyLimits
 List the high energy limits. More...
 
std::map< G4String, std::map< G4String, G4double > > fLowEnergyLimits
 List the low energy limits. More...
 
std::vector< G4StringfModelCSFiles
 List the cross section data files. More...
 
std::vector< G4StringfModelDiffCSFiles
 List the differential corss section data files. More...
 
std::vector< G4StringfModelMaterials
 List the materials that can be activated (and will be by default) within the model. More...
 
std::vector< G4StringfModelParticles
 List the particles that can be activated within the model. More...
 
std::vector< G4doublefModelScaleFactors
 List the model scale factors (they could change with material) More...
 
G4String fName
 model name More...
 
const G4String fStringOfMaterials
 fStringOfMaterials The user can decide to specify by hand which are the materials the be activated among those implemented in the model. If the user does then only the specified materials contained in this string variable will be activated. The string is like: mat1/mat2/mat3/mat4 More...
 
TableMapData fTableData
 fTableData It contains the cross section data and can be used like: dataTable=fTableData[material][particle] More...
 
MapMeanEnergy tableMeanEnergyPTB
 map: [materialName]=energyValue More...
 
G4int verboseLevel
 verbose level More...
 
G4DNAWaterExcitationStructure waterStructure
 

Detailed Description

The G4DNAPTBExcitationModel class This class implements the PTB excitation model.

Definition at line 50 of file G4DNAPTBExcitationModel.hh.

Member Typedef Documentation

◆ ItCompoMapData

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

Definition at line 185 of file G4VDNAModel.hh.

◆ MapMeanEnergy

typedef std::map<G4String,G4double,std::less<G4String> > G4DNAPTBExcitationModel::MapMeanEnergy
private

Definition at line 120 of file G4DNAPTBExcitationModel.hh.

◆ RatioMapData

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

Definition at line 184 of file G4VDNAModel.hh.

◆ TableMapData

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

Definition at line 183 of file G4VDNAModel.hh.

Constructor & Destructor Documentation

◆ G4DNAPTBExcitationModel() [1/2]

G4DNAPTBExcitationModel::G4DNAPTBExcitationModel ( const G4String applyToMaterial = "all",
const G4ParticleDefinition p = 0,
const G4String nam = "DNAPTBExcitationModel" 
)

G4DNAPTBExcitationModel Constructor.

Parameters
applyToMaterial
p
nam

Definition at line 36 of file G4DNAPTBExcitationModel.cc.

38 : G4VDNAModel(nam, applyToMaterial)
39{
40 verboseLevel= 0;
41 // Verbosity scale:
42 // 0 = nothing
43 // 1 = warning for energy non-conservation
44 // 2 = details of energy budget
45 // 3 = calculation of cross sections, file openings, sampling of atoms
46 // 4 = entering in methods
47
48 // initialisation of mean energy loss for each material
49 tableMeanEnergyPTB["THF"] = 8.01*eV;
50 tableMeanEnergyPTB["PY"] = 7.61*eV;
51 tableMeanEnergyPTB["PU"] = 7.61*eV;
52 tableMeanEnergyPTB["TMP"] = 8.01*eV;
53
54 if( verboseLevel>0 )
55 {
56 G4cout << "PTB excitation model is constructed " << G4endl;
57 }
58}
static constexpr double eV
Definition: G4SIunits.hh:201
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
MapMeanEnergy tableMeanEnergyPTB
map: [materialName]=energyValue
G4VDNAModel(const G4String &nam, const G4String &applyToMaterial)
G4VDNAModel Constructeur of the G4VDNAModel class.
Definition: G4VDNAModel.cc:35

References eV, G4cout, G4endl, tableMeanEnergyPTB, and verboseLevel.

◆ ~G4DNAPTBExcitationModel()

G4DNAPTBExcitationModel::~G4DNAPTBExcitationModel ( )
virtual

~G4DNAPTBExcitationModel Destructor

Definition at line 62 of file G4DNAPTBExcitationModel.cc.

63{
64
65}

◆ G4DNAPTBExcitationModel() [2/2]

G4DNAPTBExcitationModel::G4DNAPTBExcitationModel ( const G4DNAPTBExcitationModel )
private

Member Function Documentation

◆ AddCrossSectionData() [1/2]

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

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}
std::vector< G4double > fModelScaleFactors
List the model scale factors (they could change with material)
Definition: G4VDNAModel.hh:298
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

References G4VDNAModel::fModelCSFiles, G4VDNAModel::fModelMaterials, G4VDNAModel::fModelParticles, and G4VDNAModel::fModelScaleFactors.

◆ AddCrossSectionData() [2/2]

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

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}
std::vector< G4String > fModelDiffCSFiles
List the differential corss section data files.
Definition: G4VDNAModel.hh:297

References G4VDNAModel::fModelCSFiles, G4VDNAModel::fModelDiffCSFiles, G4VDNAModel::fModelMaterials, G4VDNAModel::fModelParticles, and G4VDNAModel::fModelScaleFactors.

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

◆ BuildApplyToMatVect()

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

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 G4VDNAModel::LoadCrossSectionData().

◆ CrossSectionPerVolume()

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

CrossSectionPerVolume Retrieve the cross section corresponding to the current material, particle and energy.

Parameters
material
materialName
p
ekin
emin
emax
Returns
the cross section value

Implements G4VDNAModel.

Definition at line 186 of file G4DNAPTBExcitationModel.cc.

192{
193 if (verboseLevel > 3)
194 G4cout << "Calling CrossSectionPerVolume() of G4DNAPTBExcitationModel" << G4endl;
195
196 // Get the name of the current particle
197 G4String particleName = particleDefinition->GetParticleName();
198
199 // initialise variables
200 G4double lowLim = 0;
201 G4double highLim = 0;
202 G4double sigma=0;
203
204 // Get the low energy limit for the current particle
205 lowLim = GetLowELimit(materialName, particleName);
206
207 // Get the high energy limit for the current particle
208 highLim = GetHighELimit(materialName, particleName);
209
210 // Check that we are in the correct energy range
211 if (ekin >= lowLim && ekin < highLim)
212 {
213 // Get the map with all the data tables
214 TableMapData* tableData = GetTableData();
215
216 // Retrieve the cross section value
217 sigma = (*tableData)[materialName][particleName]->FindValue(ekin);
218
219 if (verboseLevel > 2)
220 {
221 G4cout << "__________________________________" << G4endl;
222 G4cout << "°°° G4DNAPTBExcitationModel - XS INFO START" << G4endl;
223 G4cout << "°°° Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
224 G4cout << "°°° Cross section per "<< materialName <<" molecule (cm^2)=" << sigma/cm/cm << G4endl;
225 G4cout << "°°° G4DNAPTBExcitationModel - XS INFO END" << G4endl;
226 }
227
228 }
229
230 // Return the cross section value
231 return sigma;
232}
static constexpr double cm
Definition: G4SIunits.hh:99
double G4double
Definition: G4Types.hh:83
TableMapData * GetTableData()
GetTableData.
Definition: G4VDNAModel.hh:193
std::map< G4String, std::map< G4String, G4DNACrossSectionDataSet *, std::less< G4String > > > TableMapData
Definition: G4VDNAModel.hh:183
G4double GetHighELimit(const G4String &material, const G4String &particle)
GetHighEnergyLimit.
Definition: G4VDNAModel.hh:153
G4double GetLowELimit(const G4String &material, const G4String &particle)
GetLowEnergyLimit.
Definition: G4VDNAModel.hh:161

References cm, eV, G4cout, G4endl, G4VDNAModel::GetHighELimit(), G4VDNAModel::GetLowELimit(), G4ParticleDefinition::GetParticleName(), G4VDNAModel::GetTableData(), and verboseLevel.

◆ EnableForMaterialAndParticle()

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

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}
TableMapData fTableData
fTableData It contains the cross section data and can be used like: dataTable=fTableData[material][pa...
Definition: G4VDNAModel.hh:292

References G4VDNAModel::fTableData.

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

◆ GetHighELimit()

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

GetHighEnergyLimit.

Parameters
material
particle
Returns
fHighEnergyLimits[material][particle]

Definition at line 153 of file G4VDNAModel.hh.

153{return fHighEnergyLimits[material][particle];}
std::map< G4String, std::map< G4String, G4double > > fHighEnergyLimits
List the high energy limits.
Definition: G4VDNAModel.hh:301
string material
Definition: eplot.py:19

References G4VDNAModel::fHighEnergyLimits, and eplot::material.

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

◆ GetLowELimit()

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

GetLowEnergyLimit.

Parameters
material
particle
Returns
fLowEnergyLimits[material][particle]

Definition at line 161 of file G4VDNAModel.hh.

161{return fLowEnergyLimits[material][particle];}
std::map< G4String, std::map< G4String, G4double > > fLowEnergyLimits
List the low energy limits.
Definition: G4VDNAModel.hh:300

References G4VDNAModel::fLowEnergyLimits, and eplot::material.

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

◆ GetName()

G4String G4VDNAModel::GetName ( )
inlineinherited

GetName.

Returns
the name of the model

Definition at line 145 of file G4VDNAModel.hh.

145{return fName;}
G4String fName
model name
Definition: G4VDNAModel.hh:303

References G4VDNAModel::fName.

Referenced by G4VDNAModel::IsMaterialDefine().

◆ GetTableData()

TableMapData * G4VDNAModel::GetTableData ( )
inlineprotectedinherited

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;}

References G4VDNAModel::fTableData.

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

◆ Initialise()

void G4DNAPTBExcitationModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector = *(new G4DataVector()),
G4ParticleChangeForGamma fpChangeForGamme = nullptr 
)
virtual

Initialise Set the materials for which the model can be used and defined the energy limits.

Implements G4VDNAModel.

Definition at line 69 of file G4DNAPTBExcitationModel.cc.

71{
72 if (verboseLevel > 3)
73 G4cout << "Calling G4DNAPTBExcitationModel::Initialise()" << G4endl;
74
75 G4double scaleFactor = 1e-16*cm*cm;
76 G4double scaleFactorBorn = (1.e-22 / 3.343) * m*m;
77
79
80 //*******************************************************
81 // Cross section data
82 //*******************************************************
83
84 if(particle == electronDef)
85 {
86 G4String particleName = particle->GetParticleName();
87
89 particleName,
90 "dna/sigma_excitation_e-_PTB_THF",
91 scaleFactor);
92 SetLowELimit("THF", particleName, 9.*eV);
93 SetHighELimit("THF", particleName, 1.*keV);
94
96 particleName,
97 "dna/sigma_excitation_e-_PTB_PY",
98 scaleFactor);
99 SetLowELimit("PY", particleName, 9.*eV);
100 SetHighELimit("PY", particleName, 1.*keV);
101
103 particleName,
104 "dna/sigma_excitation_e-_PTB_PU",
105 scaleFactor);
106 SetLowELimit("PU", particleName, 9.*eV);
107 SetHighELimit("PU", particleName, 1.*keV);
108
110 particleName,
111 "dna/sigma_excitation_e-_PTB_TMP",
112 scaleFactor);
113 SetLowELimit("TMP", particleName, 9.*eV);
114 SetHighELimit("TMP", particleName, 1.*keV);
115
116 AddCrossSectionData("G4_WATER",
117 particleName,
118 "dna/sigma_excitation_e_born",
119 scaleFactorBorn);
120 SetLowELimit("G4_WATER", particleName, 9.*eV);
121 SetHighELimit("G4_WATER", particleName, 1.*keV);
122
123 // DNA materials
124 //
125 AddCrossSectionData("backbone_THF",
126 particleName,
127 "dna/sigma_excitation_e-_PTB_THF",
128 scaleFactor*33./30);
129 SetLowELimit("backbone_THF", particleName, 9.*eV);
130 SetHighELimit("backbone_THF", particleName, 1.*keV);
131
132 AddCrossSectionData("cytosine_PY",
133 particleName,
134 "dna/sigma_excitation_e-_PTB_PY",
135 scaleFactor*42./30);
136 SetLowELimit("cytosine_PY", particleName, 9.*eV);
137 SetHighELimit("cytosine_PY", particleName, 1.*keV);
138
139 AddCrossSectionData("thymine_PY",
140 particleName,
141 "dna/sigma_excitation_e-_PTB_PY",
142 scaleFactor*48./30);
143 SetLowELimit("thymine_PY", particleName, 9.*eV);
144 SetHighELimit("thymine_PY", particleName, 1.*keV);
145
146 AddCrossSectionData("adenine_PU",
147 particleName,
148 "dna/sigma_excitation_e-_PTB_PU",
149 scaleFactor*50./44);
150 SetLowELimit("adenine_PU", particleName, 9.*eV);
151 SetHighELimit("adenine_PU", particleName, 1.*keV);
152
153 AddCrossSectionData("guanine_PU",
154 particleName,
155 "dna/sigma_excitation_e-_PTB_PU",
156 scaleFactor*56./44);
157 SetLowELimit("guanine_PU", particleName, 9.*eV);
158 SetHighELimit("guanine_PU", particleName, 1.*keV);
159
160 AddCrossSectionData("backbone_TMP",
161 particleName,
162 "dna/sigma_excitation_e-_PTB_TMP",
163 scaleFactor*33./50);
164 SetLowELimit("backbone_TMP", particleName, 9.*eV);
165 SetHighELimit("backbone_TMP", particleName, 1.*keV);
166 }
167
168 //*******************************************************
169 // Load data
170 //*******************************************************
171
173
174 //*******************************************************
175 // Verbose
176 //*******************************************************
177
178 if( verboseLevel>0 )
179 {
180 G4cout << "PTB excitation model is initialized " << G4endl;
181 }
182}
static constexpr double m
Definition: G4SIunits.hh:109
static constexpr double keV
Definition: G4SIunits.hh:202
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:88
const G4String & GetParticleName() const
void SetHighELimit(const G4String &material, const G4String &particle, G4double lim)
SetHighEnergyLimit.
Definition: G4VDNAModel.hh:169
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
void SetLowELimit(const G4String &material, const G4String &particle, G4double lim)
SetLowEnergyLimit.
Definition: G4VDNAModel.hh:177
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

References G4VDNAModel::AddCrossSectionData(), cm, G4Electron::ElectronDefinition(), eV, G4cout, G4endl, G4ParticleDefinition::GetParticleName(), keV, G4VDNAModel::LoadCrossSectionData(), m, G4VDNAModel::SetHighELimit(), G4VDNAModel::SetLowELimit(), and verboseLevel.

◆ IsMaterialDefine()

G4bool G4VDNAModel::IsMaterialDefine ( const G4String materialName)
inherited

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:672
G4String GetName()
GetName.
Definition: G4VDNAModel.hh:145

References G4Material::GetMaterialTable(), and G4VDNAModel::GetName().

◆ IsMaterialExistingInModel()

G4bool G4VDNAModel::IsMaterialExistingInModel ( const G4String materialName)
inherited

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}

References G4VDNAModel::fTableData.

Referenced by G4VDNAModel::IsParticleExistingInModelForMaterial().

◆ IsParticleExistingInModelForMaterial()

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

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

References G4VDNAModel::fTableData, and G4VDNAModel::IsMaterialExistingInModel().

◆ LoadCrossSectionData()

void G4VDNAModel::LoadCrossSectionData ( const G4String particleName)
protectedinherited

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: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 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

References G4VDNAModel::BuildApplyToMatVect(), FatalException, G4VDNAModel::fModelCSFiles, G4VDNAModel::fModelDiffCSFiles, G4VDNAModel::fModelMaterials, G4VDNAModel::fModelParticles, G4VDNAModel::fModelScaleFactors, G4VDNAModel::fStringOfMaterials, G4Exception(), G4VDNAModel::ReadAndSaveCSFile(), and G4VDNAModel::ReadDiffCSFile().

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

◆ operator=()

G4DNAPTBExcitationModel & G4DNAPTBExcitationModel::operator= ( const G4DNAPTBExcitationModel right)
private

◆ RandomSelectShell()

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

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 size_t n(table->NumberOfComponents());
199 size_t 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}
static const G4double pos
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
virtual G4double FindValue(G4double x, G4int componentId=0) const =0

References FatalException, G4VEMDataSet::FindValue(), G4Exception(), G4UniformRand, G4DNACrossSectionDataSet::GetComponent(), G4VDNAModel::GetTableData(), CLHEP::detail::n, G4DNACrossSectionDataSet::NumberOfComponents(), and pos.

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

◆ ReadAndSaveCSFile()

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

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}

References eV, geant4_check_module_cycles::file, and G4VDNAModel::fTableData.

Referenced by G4VDNAModel::LoadCrossSectionData().

◆ ReadDiffCSFile()

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

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

Reimplemented in G4DNAPTBIonisationModel, and G4DNAPTBElasticModel.

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}

References FatalException, and G4Exception().

Referenced by G4VDNAModel::LoadCrossSectionData().

◆ SampleSecondaries()

void G4DNAPTBExcitationModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple ,
const G4String materialName,
const G4DynamicParticle aDynamicParticle,
G4ParticleChangeForGamma particleChangeForGamma,
G4double  tmin,
G4double  tmax 
)
virtual

SampleSecondaries If the model is selected for the ModelInterface then the SampleSecondaries method will be called. The method sets the incident particle characteristics after the ModelInterface.

Parameters
materialName
particleChangeForGamma
tmin
tmax

Implements G4VDNAModel.

Definition at line 236 of file G4DNAPTBExcitationModel.cc.

243{
244 if (verboseLevel > 3)
245 G4cout << "Calling SampleSecondaries() of G4DNAPTBExcitationModel" << G4endl;
246
247 // Get the incident particle kinetic energy
248 G4double k = aDynamicParticle->GetKineticEnergy();
249
250 if(materialName!="G4_WATER")
251 {
252 // Retrieve the excitation energy for the current material
253 G4double excitationEnergy = tableMeanEnergyPTB[materialName];
254
255 // Calculate the new energy of the particle
256 G4double newEnergy = k - excitationEnergy;
257
258 // Check that the new energy is above zero before applying it the particle.
259 // Otherwise, do nothing.
260 if (newEnergy > 0)
261 {
262 particleChangeForGamma->ProposeMomentumDirection(aDynamicParticle->GetMomentumDirection());
263 particleChangeForGamma->SetProposedKineticEnergy(newEnergy);
264 particleChangeForGamma->ProposeLocalEnergyDeposit(excitationEnergy);
265 }
266 }
267 else
268 {
269 const G4String& particleName = aDynamicParticle->GetDefinition()->GetParticleName();
270
271 G4int level = RandomSelectShell(k,particleName, materialName);
272 G4double excitationEnergy = waterStructure.ExcitationEnergy(level);
273 G4double newEnergy = k - excitationEnergy;
274
275 if (newEnergy > 0)
276 {
277 particleChangeForGamma->ProposeMomentumDirection(aDynamicParticle->GetMomentumDirection());
278 particleChangeForGamma->SetProposedKineticEnergy(newEnergy);
279 particleChangeForGamma->ProposeLocalEnergyDeposit(excitationEnergy);
280 }
281
282 const G4Track * theIncomingTrack = particleChangeForGamma->GetCurrentTrack();
284 level,
285 theIncomingTrack);
286 }
287}
@ eExcitedMolecule
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
G4DNAWaterExcitationStructure waterStructure
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
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
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

References G4DNAChemistryManager::CreateWaterMolecule(), eExcitedMolecule, G4DNAWaterExcitationStructure::ExcitationEnergy(), G4cout, G4endl, G4ParticleChangeForGamma::GetCurrentTrack(), G4DynamicParticle::GetDefinition(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), G4ParticleDefinition::GetParticleName(), G4DNAChemistryManager::Instance(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChangeForGamma::ProposeMomentumDirection(), G4VDNAModel::RandomSelectShell(), G4ParticleChangeForGamma::SetProposedKineticEnergy(), tableMeanEnergyPTB, verboseLevel, and waterStructure.

◆ SetHighELimit()

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

SetHighEnergyLimit.

Parameters
material
particle
lim

Definition at line 169 of file G4VDNAModel.hh.

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

References G4VDNAModel::fHighEnergyLimits, and eplot::material.

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

◆ SetLowELimit()

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

SetLowEnergyLimit.

Parameters
material
particle
lim

Definition at line 177 of file G4VDNAModel.hh.

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

References G4VDNAModel::fLowEnergyLimits, and eplot::material.

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

Field Documentation

◆ fHighEnergyLimits

std::map<G4String, std::map<G4String, G4double> > G4VDNAModel::fHighEnergyLimits
privateinherited

List the high energy limits.

Definition at line 301 of file G4VDNAModel.hh.

Referenced by G4VDNAModel::GetHighELimit(), and G4VDNAModel::SetHighELimit().

◆ fLowEnergyLimits

std::map<G4String, std::map<G4String, G4double> > G4VDNAModel::fLowEnergyLimits
privateinherited

List the low energy limits.

Definition at line 300 of file G4VDNAModel.hh.

Referenced by G4VDNAModel::GetLowELimit(), and G4VDNAModel::SetLowELimit().

◆ fModelCSFiles

std::vector<G4String> G4VDNAModel::fModelCSFiles
privateinherited

List the cross section data files.

Definition at line 296 of file G4VDNAModel.hh.

Referenced by G4VDNAModel::AddCrossSectionData(), and G4VDNAModel::LoadCrossSectionData().

◆ fModelDiffCSFiles

std::vector<G4String> G4VDNAModel::fModelDiffCSFiles
privateinherited

List the differential corss section data files.

Definition at line 297 of file G4VDNAModel.hh.

Referenced by G4VDNAModel::AddCrossSectionData(), and G4VDNAModel::LoadCrossSectionData().

◆ fModelMaterials

std::vector<G4String> G4VDNAModel::fModelMaterials
privateinherited

List the materials that can be activated (and will be by default) within the model.

Definition at line 294 of file G4VDNAModel.hh.

Referenced by G4VDNAModel::AddCrossSectionData(), and G4VDNAModel::LoadCrossSectionData().

◆ fModelParticles

std::vector<G4String> G4VDNAModel::fModelParticles
privateinherited

List the particles that can be activated within the model.

Definition at line 295 of file G4VDNAModel.hh.

Referenced by G4VDNAModel::AddCrossSectionData(), and G4VDNAModel::LoadCrossSectionData().

◆ fModelScaleFactors

std::vector<G4double> G4VDNAModel::fModelScaleFactors
privateinherited

List the model scale factors (they could change with material)

Definition at line 298 of file G4VDNAModel.hh.

Referenced by G4VDNAModel::AddCrossSectionData(), and G4VDNAModel::LoadCrossSectionData().

◆ fName

G4String G4VDNAModel::fName
privateinherited

model name

Definition at line 303 of file G4VDNAModel.hh.

Referenced by G4VDNAModel::GetName().

◆ fStringOfMaterials

const G4String G4VDNAModel::fStringOfMaterials
privateinherited

fStringOfMaterials The user can decide to specify by hand which are the materials the be activated among those implemented in the model. If the user does then only the specified materials contained in this string variable will be activated. The string is like: mat1/mat2/mat3/mat4

Definition at line 286 of file G4VDNAModel.hh.

Referenced by G4VDNAModel::LoadCrossSectionData().

◆ fTableData

TableMapData G4VDNAModel::fTableData
privateinherited

fTableData It contains the cross section data and can be used like: dataTable=fTableData[material][particle]

Definition at line 292 of file G4VDNAModel.hh.

Referenced by G4VDNAModel::EnableForMaterialAndParticle(), G4VDNAModel::GetTableData(), G4VDNAModel::IsMaterialExistingInModel(), G4VDNAModel::IsParticleExistingInModelForMaterial(), G4VDNAModel::ReadAndSaveCSFile(), and G4VDNAModel::~G4VDNAModel().

◆ tableMeanEnergyPTB

MapMeanEnergy G4DNAPTBExcitationModel::tableMeanEnergyPTB
private

map: [materialName]=energyValue

Definition at line 121 of file G4DNAPTBExcitationModel.hh.

Referenced by G4DNAPTBExcitationModel(), and SampleSecondaries().

◆ verboseLevel

G4int G4DNAPTBExcitationModel::verboseLevel
private

◆ waterStructure

G4DNAWaterExcitationStructure G4DNAPTBExcitationModel::waterStructure
private

Definition at line 118 of file G4DNAPTBExcitationModel.hh.

Referenced by SampleSecondaries().


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