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

The G4DNAPTBElasticModel class This class implements the elastic model for the DNA materials and precursors. More...

#include <G4DNAPTBElasticModel.hh>

Inheritance diagram for G4DNAPTBElasticModel:
G4VDNAModel

Public Member Functions

virtual G4double CrossSectionPerVolume (const G4Material *material, const G4String &materialName, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
 CrossSectionPerVolume This method is mandatory for any model class. It finds and return the cross section value for the current material, particle and energy values. The number of molecule per volume is not used here but in the G4DNAModelInterface class. More...
 
 G4DNAPTBElasticModel (const G4String &applyToMaterial="all", const G4ParticleDefinition *p=0, const G4String &nam="DNAPTBElasticModel")
 G4DNAPTBElasticModel 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 &, G4ParticleChangeForGamma *fpChangeForGamme=nullptr)
 Initialise Mandatory method for every model class. The material/particle for which the model can be used have to be added here through the AddCrossSectionData method. Then the LoadCrossSectionData method must be called to trigger the load process. Scale factors to be applied to the cross section can be defined here. 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 Method called after CrossSectionPerVolume if the process is the one which is selected (according to the sampling on the calculated path length). Here, the characteristics of the incident and created (if any) particle(s) are set (energy, momentum ...). 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 ~G4DNAPTBElasticModel ()
 ~G4DNAPTBElasticModel 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...
 

Private Types

typedef std::map< G4String, std::map< G4String, std::map< double, std::map< double, double > > > > TriDimensionMap
 
typedef std::map< G4String, std::map< G4String, std::map< double, std::vector< double > > > > VecMap
 

Private Member Functions

 G4DNAPTBElasticModel (G4DNAPTBElasticModel &)
 
G4double LinLinInterpolate (G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
 LinLinInterpolate. More...
 
G4double LinLogInterpolate (G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
 LinLogInterpolate. More...
 
G4double LogLogInterpolate (G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
 LogLogInterpolate. More...
 
G4DNAPTBElasticModeloperator= (const G4DNAPTBElasticModel &right)
 
G4double QuadInterpolator (G4double e11, G4double e12, G4double e21, G4double e22, G4double x11, G4double x12, G4double x21, G4double x22, G4double t1, G4double t2, G4double t, G4double e)
 QuadInterpolator. More...
 
G4double RandomizeCosTheta (G4double k, const G4String &materialName)
 RandomizeCosTheta. More...
 
void ReadDiffCSFile (const G4String &materialName, const G4String &particleName, const G4String &file, const G4double)
 ReadDiffCSFile Method to read the differential cross section files. This method is not standard yet so every model must implement its own. More...
 
G4double Theta (G4ParticleDefinition *fParticleDefinition, G4double k, G4double integrDiff, const G4String &materialName)
 Theta To return an angular theta value from the differential file. This method uses interpolations to calculate the theta value. More...
 

Private Attributes

TriDimensionMap diffCrossSectionData
 A map: [materialName][particleName]=DiffCrossSectionTable. More...
 
VecMap eValuesVect
 
std::map< G4String, std::map< G4String, G4double > > fHighEnergyLimits
 List the high energy limits. More...
 
G4double fKillBelowEnergy
 energy kill limit 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...
 
std::map< G4String, double > killBelowEnergyTable
 map to save the different energy kill limits for the materials More...
 
std::map< G4String, std::map< G4String, std::vector< double > > > tValuesVec
 map with vectors containing all the incident (T) energy of the differential file More...
 
G4int verboseLevel
 verbose level More...
 

Detailed Description

The G4DNAPTBElasticModel class This class implements the elastic model for the DNA materials and precursors.

Definition at line 48 of file G4DNAPTBElasticModel.hh.

Member Typedef Documentation

◆ ItCompoMapData

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

Definition at line 185 of file G4VDNAModel.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.

◆ TriDimensionMap

typedef std::map<G4String, std::map<G4String, std::map<double, std::map<double, double> > > > G4DNAPTBElasticModel::TriDimensionMap
private

Definition at line 125 of file G4DNAPTBElasticModel.hh.

◆ VecMap

typedef std::map<G4String, std::map<G4String, std::map<double, std::vector<double> > > > G4DNAPTBElasticModel::VecMap
private

Definition at line 128 of file G4DNAPTBElasticModel.hh.

Constructor & Destructor Documentation

◆ G4DNAPTBElasticModel() [1/2]

G4DNAPTBElasticModel::G4DNAPTBElasticModel ( const G4String applyToMaterial = "all",
const G4ParticleDefinition p = 0,
const G4String nam = "DNAPTBElasticModel" 
)

G4DNAPTBElasticModel Constructor.

Parameters
applyToMaterial
p
nam

Definition at line 38 of file G4DNAPTBElasticModel.cc.

40 : G4VDNAModel(nam, applyToMaterial)
41{
42 fKillBelowEnergy = 10*eV; // will be override by the limits defined for each material
43
44 verboseLevel= 0;
45 // Verbosity scale:
46 // 0 = nothing
47 // 1 = warning for energy non-conservation
48 // 2 = details of energy budget
49 // 3 = calculation of cross sections, file openings, sampling of atoms
50 // 4 = entering in methods
51
52 if( verboseLevel>0 )
53 {
54 G4cout << "PTB Elastic model is constructed " << G4endl;
55 }
56}
static constexpr double eV
Definition: G4SIunits.hh:201
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double fKillBelowEnergy
energy kill limit
G4int verboseLevel
verbose level
G4VDNAModel(const G4String &nam, const G4String &applyToMaterial)
G4VDNAModel Constructeur of the G4VDNAModel class.
Definition: G4VDNAModel.cc:35

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

◆ ~G4DNAPTBElasticModel()

G4DNAPTBElasticModel::~G4DNAPTBElasticModel ( )
virtual

~G4DNAPTBElasticModel Destructor

Definition at line 60 of file G4DNAPTBElasticModel.cc.

61{
62
63}

◆ G4DNAPTBElasticModel() [2/2]

G4DNAPTBElasticModel::G4DNAPTBElasticModel ( G4DNAPTBElasticModel )
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 Initialise(), G4DNAPTBExcitationModel::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 G4DNAPTBElasticModel::CrossSectionPerVolume ( const G4Material material,
const G4String materialName,
const G4ParticleDefinition p,
G4double  ekin,
G4double  emin,
G4double  emax 
)
virtual

CrossSectionPerVolume This method is mandatory for any model class. It finds and return the cross section value for the current material, particle and energy values. The number of molecule per volume is not used here but in the G4DNAModelInterface class.

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

Implements G4VDNAModel.

Definition at line 295 of file G4DNAPTBElasticModel.cc.

301{
302 if (verboseLevel > 3)
303 G4cout << "Calling CrossSectionPerVolume() of G4DNAPTBElasticModel" << G4endl;
304
305 // Get the name of the current particle
306 const G4String& particleName = p->GetParticleName();
307
308 // set killBelowEnergy value for current material
309 fKillBelowEnergy = GetLowELimit(materialName, particleName);
310
311 // initialise the return value (cross section) to zero
312 G4double sigma(0);
313
314 // check if we are below the high energy limit
315 if (ekin < GetHighELimit(materialName, particleName) )
316 {
317 // This is used to kill the particle if its kinetic energy is below fKillBelowEnergy.
318 // If the energy is lower then we return a maximum cross section and thus the SampleSecondaries method will be called for sure.
319 // SampleSecondaries will remove the particle from the simulation.
320 //
321 //SI : XS must not be zero otherwise sampling of secondaries method ignored
322 if (ekin < fKillBelowEnergy) return DBL_MAX;
323
324 // Get the tables with the cross section data
325 TableMapData* tableData = GetTableData();
326
327 // Retrieve the cross section value
328 sigma = (*tableData)[materialName][particleName]->FindValue(ekin);
329 }
330
331 if (verboseLevel > 2)
332 {
333 G4cout << "__________________________________" << G4endl;
334 G4cout << "°°° G4DNAPTBElasticModel - XS INFO START" << G4endl;
335 G4cout << "°°° Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
336 G4cout << "°°° Cross section per molecule (cm^2)=" << sigma/cm/cm << G4endl;
337 G4cout << "°°° G4DNAPTBElasticModel - XS INFO END" << G4endl;
338 }
339
340 // Return the cross section
341 return sigma;
342}
static constexpr double cm
Definition: G4SIunits.hh:99
double G4double
Definition: G4Types.hh:83
const G4String & GetParticleName() const
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
#define DBL_MAX
Definition: templates.hh:62

References cm, DBL_MAX, eV, fKillBelowEnergy, 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 CrossSectionPerVolume(), G4DNAPTBExcitationModel::CrossSectionPerVolume(), G4DNAPTBIonisationModel::CrossSectionPerVolume(), 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 CrossSectionPerVolume(), G4DNAPTBExcitationModel::CrossSectionPerVolume(), G4DNAPTBIonisationModel::CrossSectionPerVolume(), 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 CrossSectionPerVolume(), G4DNAPTBExcitationModel::CrossSectionPerVolume(), G4DNAPTBIonisationModel::CrossSectionPerVolume(), and G4VDNAModel::RandomSelectShell().

◆ Initialise()

void G4DNAPTBElasticModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector ,
G4ParticleChangeForGamma fpChangeForGamme = nullptr 
)
virtual

Initialise Mandatory method for every model class. The material/particle for which the model can be used have to be added here through the AddCrossSectionData method. Then the LoadCrossSectionData method must be called to trigger the load process. Scale factors to be applied to the cross section can be defined here.

Implements G4VDNAModel.

Definition at line 67 of file G4DNAPTBElasticModel.cc.

69{
70 if (verboseLevel > 3)
71 G4cout << "Calling G4DNAPTBElasticModel::Initialise()" << G4endl;
72
73 G4double scaleFactor = 1e-16*cm*cm;
74
76
77 //*******************************************************
78 // Cross section data
79 //*******************************************************
80
81 if(particle == electronDef)
82 {
83 G4String particleName = particle->GetParticleName();
84
86 particleName,
87 "dna/sigma_elastic_e-_PTB_THF",
88 "dna/sigmadiff_cumulated_elastic_e-_PTB_THF",
89 scaleFactor);
90 SetLowELimit("THF", particleName, 10*eV);
91 SetHighELimit("THF", particleName, 1*keV);
92
94 particleName,
95 "dna/sigma_elastic_e-_PTB_PY",
96 "dna/sigmadiff_cumulated_elastic_e-_PTB_PY",
97 scaleFactor);
98 SetLowELimit("PY", particleName, 10*eV);
99 SetHighELimit("PY", particleName, 1*keV);
100
102 particleName,
103 "dna/sigma_elastic_e-_PTB_PU",
104 "dna/sigmadiff_cumulated_elastic_e-_PTB_PU",
105 scaleFactor);
106 SetLowELimit("PU", particleName, 10*eV);
107 SetHighELimit("PU", particleName, 1*keV);
108
110 particleName,
111 "dna/sigma_elastic_e-_PTB_TMP",
112 "dna/sigmadiff_cumulated_elastic_e-_PTB_TMP",
113 scaleFactor);
114 SetLowELimit("TMP", particleName, 10*eV);
115 SetHighELimit("TMP", particleName, 1*keV);
116
117 AddCrossSectionData("G4_WATER",
118 particleName,
119 "dna/sigma_elastic_e_champion",
120 "dna/sigmadiff_cumulated_elastic_e_champion",
121 scaleFactor);
122 SetLowELimit("G4_WATER", particleName, 10*eV);
123 SetHighELimit("G4_WATER", particleName, 1*keV);
124
125 // DNA materials
126 //
127 AddCrossSectionData("backbone_THF",
128 particleName,
129 "dna/sigma_elastic_e-_PTB_THF",
130 "dna/sigmadiff_cumulated_elastic_e-_PTB_THF",
131 scaleFactor*33./30);
132 SetLowELimit("backbone_THF", particleName, 10*eV);
133 SetHighELimit("backbone_THF", particleName, 1*keV);
134
135 AddCrossSectionData("cytosine_PY",
136 particleName,
137 "dna/sigma_elastic_e-_PTB_PY",
138 "dna/sigmadiff_cumulated_elastic_e-_PTB_PY",
139 scaleFactor*42./30);
140 SetLowELimit("cytosine_PY", particleName, 10*eV);
141 SetHighELimit("cytosine_PY", particleName, 1*keV);
142
143 AddCrossSectionData("thymine_PY",
144 particleName,
145 "dna/sigma_elastic_e-_PTB_PY",
146 "dna/sigmadiff_cumulated_elastic_e-_PTB_PY",
147 scaleFactor*48./30);
148 SetLowELimit("thymine_PY", particleName, 10*eV);
149 SetHighELimit("thymine_PY", particleName, 1*keV);
150
151 AddCrossSectionData("adenine_PU",
152 particleName,
153 "dna/sigma_elastic_e-_PTB_PU",
154 "dna/sigmadiff_cumulated_elastic_e-_PTB_PU",
155 scaleFactor*50./44);
156 SetLowELimit("adenine_PU", particleName, 10*eV);
157 SetHighELimit("adenine_PU", particleName, 1*keV);
158
159 AddCrossSectionData("guanine_PU",
160 particleName,
161 "dna/sigma_elastic_e-_PTB_PU",
162 "dna/sigmadiff_cumulated_elastic_e-_PTB_PU",
163 scaleFactor*56./44);
164 SetLowELimit("guanine_PU", particleName, 10*eV);
165 SetHighELimit("guanine_PU", particleName, 1*keV);
166
167 AddCrossSectionData("backbone_TMP",
168 particleName,
169 "dna/sigma_elastic_e-_PTB_TMP",
170 "dna/sigmadiff_cumulated_elastic_e-_PTB_TMP",
171 scaleFactor*33./50);
172 SetLowELimit("backbone_TMP", particleName, 10*eV);
173 SetHighELimit("backbone_TMP", particleName, 1*keV);
174 }
175
176 //*******************************************************
177 // Load the data
178 //*******************************************************
179
181
182 //*******************************************************
183 // Verbose output
184 //*******************************************************
185
186 if (verboseLevel > 2)
187 G4cout << "Loaded cross section files for PTB Elastic model" << G4endl;
188
189 if( verboseLevel>0 )
190 {
191 G4cout << "PTB Elastic model is initialized " << G4endl;
192 }
193}
static constexpr double keV
Definition: G4SIunits.hh:202
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:88
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(), 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().

◆ LinLinInterpolate()

G4double G4DNAPTBElasticModel::LinLinInterpolate ( G4double  e1,
G4double  e2,
G4double  e,
G4double  xs1,
G4double  xs2 
)
private

LinLinInterpolate.

Parameters
e1
e2
e
xs1
xs2
Returns

Definition at line 470 of file G4DNAPTBElasticModel.cc.

475{
476 G4double d1 = xs1;
477 G4double d2 = xs2;
478 G4double value = (d1 + (d2 - d1)*(e - e1)/ (e2 - e1));
479 return value;
480}
static const G4double e1[44]
static const G4double e2[44]
static const G4double d1
static const G4double d2

References d1, d2, e1, and e2.

Referenced by QuadInterpolator().

◆ LinLogInterpolate()

G4double G4DNAPTBElasticModel::LinLogInterpolate ( G4double  e1,
G4double  e2,
G4double  e,
G4double  xs1,
G4double  xs2 
)
private

LinLogInterpolate.

Parameters
e1
e2
e
xs1
xs2
Returns

Definition at line 456 of file G4DNAPTBElasticModel.cc.

461{
462 G4double d1 = std::log(xs1);
463 G4double d2 = std::log(xs2);
464 G4double value = std::exp(d1 + (d2 - d1)*(e - e1)/ (e2 - e1));
465 return value;
466}

References d1, d2, e1, and e2.

◆ 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 Initialise(), G4DNAPTBExcitationModel::Initialise(), and G4DNAPTBIonisationModel::Initialise().

◆ LogLogInterpolate()

G4double G4DNAPTBElasticModel::LogLogInterpolate ( G4double  e1,
G4double  e2,
G4double  e,
G4double  xs1,
G4double  xs2 
)
private

LogLogInterpolate.

Parameters
e1
e2
e
xs1
xs2
Returns

Definition at line 484 of file G4DNAPTBElasticModel.cc.

489{
490 G4double a = (std::log10(xs2)-std::log10(xs1)) / (std::log10(e2)-std::log10(e1));
491 G4double b = std::log10(xs2) - a*std::log10(e2);
492 G4double sigma = a*std::log10(e) + b;
493 G4double value = (std::pow(10.,sigma));
494 return value;
495}

References e1, and e2.

◆ operator=()

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

◆ QuadInterpolator()

G4double G4DNAPTBElasticModel::QuadInterpolator ( G4double  e11,
G4double  e12,
G4double  e21,
G4double  e22,
G4double  x11,
G4double  x12,
G4double  x21,
G4double  x22,
G4double  t1,
G4double  t2,
G4double  t,
G4double  e 
)
private

QuadInterpolator.

Parameters
e11
e12
e21
e22
x11
x12
x21
x22
t1
t2
t
e
Returns

Definition at line 499 of file G4DNAPTBElasticModel.cc.

505{
506 // Log-Log
507 /*
508 G4double interpolatedvalue1 = LogLogInterpolate(e11, e12, e, xs11, xs12);
509 G4double interpolatedvalue2 = LogLogInterpolate(e21, e22, e, xs21, xs22);
510 G4double value = LogLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
511
512
513 // Lin-Log
514 G4double interpolatedvalue1 = LinLogInterpolate(e11, e12, e, xs11, xs12);
515 G4double interpolatedvalue2 = LinLogInterpolate(e21, e22, e, xs21, xs22);
516 G4double value = LinLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
517*/
518
519 // Lin-Lin
520 G4double interpolatedvalue1 = LinLinInterpolate(e11, e12, e, xs11, xs12);
521 G4double interpolatedvalue2 = LinLinInterpolate(e21, e22, e, xs21, xs22);
522 G4double value = LinLinInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
523
524 return value;
525}
G4double LinLinInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
LinLinInterpolate.

References LinLinInterpolate().

Referenced by Theta().

◆ RandomizeCosTheta()

G4double G4DNAPTBElasticModel::RandomizeCosTheta ( G4double  k,
const G4String materialName 
)
private

RandomizeCosTheta.

Parameters
k
materialName
Returns

Definition at line 529 of file G4DNAPTBElasticModel.cc.

530{
531 G4double integrdiff=0;
533 integrdiff = uniformRand;
534
535 G4double theta=0.;
536 G4double cosTheta=0.;
537 theta = Theta(G4Electron::ElectronDefinition(),k/eV,integrdiff, materialName);
538
539 cosTheta= std::cos(theta*pi/180);
540
541 return cosTheta;
542}
static constexpr double pi
Definition: G4SIunits.hh:55
#define G4UniformRand()
Definition: Randomize.hh:52
G4double Theta(G4ParticleDefinition *fParticleDefinition, G4double k, G4double integrDiff, const G4String &materialName)
Theta To return an angular theta value from the differential file. This method uses interpolations to...

References G4Electron::ElectronDefinition(), eV, G4UniformRand, pi, Theta(), and anonymous_namespace{G4HadPhaseSpaceGenbod.cc}::uniformRand().

Referenced by SampleSecondaries().

◆ 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
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 G4DNAPTBExcitationModel::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 G4DNAPTBElasticModel::ReadDiffCSFile ( const G4String materialName,
const G4String particleName,
const G4String file,
const  G4double 
)
privatevirtual

ReadDiffCSFile Method to read the differential cross section files. This method is not standard yet so every model must implement its own.

Parameters
materialName
particleName
file

Reimplemented from G4VDNAModel.

Definition at line 197 of file G4DNAPTBElasticModel.cc.

201{
202 // Method to read and save the information contained within the differential cross section files.
203 // This method is not yet standard.
204
205 // get the path of the G4LEDATA data folder
206 char *path = std::getenv("G4LEDATA");
207 // if it is not found then quit and print error message
208 if(!path)
209 {
210 G4Exception("G4DNAPTBElasticModel::ReadAllDiffCSFiles","em0006",
211 FatalException,"G4LEDATA environment variable not set.");
212 return;
213 }
214
215 // build the fullFileName path of the data file
216 std::ostringstream fullFileName;
217 fullFileName << path <<"/"<< file<<".dat";
218
219 // open the data file
220 std::ifstream diffCrossSection (fullFileName.str().c_str());
221 // error if file is not there
222 std::stringstream endPath;
223 if (!diffCrossSection)
224 {
225 endPath << "Missing data file: "<<file;
226 G4Exception("G4DNAPTBElasticModel::Initialise","em0003",
227 FatalException, endPath.str().c_str());
228 }
229
230 tValuesVec[materialName][particleName].push_back(0.);
231
232 G4String line;
233
234 // read the file line by line until we reach the end of file point
235 while(std::getline(diffCrossSection, line))
236 {
237 // check if the line is comment or empty
238 //
239 std::istringstream testIss(line);
241 testIss >> test;
242 // check first caracter to determine if following information is data or comments
243 if(test=="#")
244 {
245 // skip the line by beginning a new while loop.
246 continue;
247 }
248 // check if line is empty
249 else if(line.empty())
250 {
251 // skip the line by beginning a new while loop.
252 continue;
253 }
254 //
255 // end of the check
256
257 // transform the line into a iss
258 std::istringstream iss(line);
259
260 // Variables to be filled by the input file
261 double tDummy;
262 double eDummy;
263
264 // fill the variables with the content of the line
265 iss>>tDummy>>eDummy;
266
267 // SI : mandatory Vecm initialization
268
269 // Fill two vectors contained in maps of types:
270 // [materialName][particleName]=vector
271 // [materialName][particleName][T]=vector
272 // to list all the incident energies (tValues) and all the output energies (eValues) within the file
273 //
274 // Check if we already have the current T value in the vector.
275 // If not then add it
276 if (tDummy != tValuesVec[materialName][particleName].back())
277 {
278 // Add the current T value
279 tValuesVec[materialName][particleName].push_back(tDummy);
280
281 // Make it correspond to a default zero E value
282 eValuesVect[materialName][particleName][tDummy].push_back(0.);
283 }
284
285 // Put the differential cross section value of the input file within the diffCrossSectionData map
286 iss>>diffCrossSectionData[materialName][particleName][tDummy][eDummy];
287
288 // If the current E value (eDummy) is different from the one already registered in the eVector then add it to the vector
289 if (eDummy != eValuesVect[materialName][particleName][tDummy].back()) eValuesVect[materialName][particleName][tDummy].push_back(eDummy);
290 }
291}
TriDimensionMap diffCrossSectionData
A map: [materialName][particleName]=DiffCrossSectionTable.
std::map< G4String, std::map< G4String, std::vector< double > > > tValuesVec
map with vectors containing all the incident (T) energy of the differential file
def test()
Definition: mcscore.py:117
Definition: test.py:1

References diffCrossSectionData, eValuesVect, FatalException, geant4_check_module_cycles::file, G4Exception(), mcscore::test(), and tValuesVec.

◆ SampleSecondaries()

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

SampleSecondaries Method called after CrossSectionPerVolume if the process is the one which is selected (according to the sampling on the calculated path length). Here, the characteristics of the incident and created (if any) particle(s) are set (energy, momentum ...).

Parameters
materialName
particleChangeForGamma
tmin
tmax

Implements G4VDNAModel.

Definition at line 346 of file G4DNAPTBElasticModel.cc.

353{
354 if (verboseLevel > 3)
355 G4cout << "Calling SampleSecondaries() of G4DNAPTBElasticModel" << G4endl;
356
357 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
358
359 const G4String& particleName = aDynamicElectron->GetParticleDefinition()->GetParticleName();
360
361 // set killBelowEnergy value for material
362 fKillBelowEnergy = GetLowELimit(materialName, particleName);
363
364 // If the particle (electron here) energy is below the kill limit then we remove it from the simulation
365 if (electronEnergy0 < fKillBelowEnergy)
366 {
367 particleChangeForGamma->SetProposedKineticEnergy(0.);
368 particleChangeForGamma->ProposeTrackStatus(fStopAndKill);
369 particleChangeForGamma->ProposeLocalEnergyDeposit(electronEnergy0);
370 }
371 // If we are above the kill limite and below the high limit then we proceed
372 else if (electronEnergy0>= fKillBelowEnergy && electronEnergy0 < GetHighELimit(materialName, particleName) )
373 {
374 // Random sampling of the cosTheta
375 G4double cosTheta = RandomizeCosTheta(electronEnergy0, materialName);
376
377 // Random sampling of phi
378 G4double phi = 2. * pi * G4UniformRand();
379
380 G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
381 G4ThreeVector xVers = zVers.orthogonal();
382 G4ThreeVector yVers = zVers.cross(xVers);
383
384 G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
385 G4double yDir = xDir;
386 xDir *= std::cos(phi);
387 yDir *= std::sin(phi);
388
389 // Particle direction after ModelInterface
390 G4ThreeVector zPrikeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
391
392 // Give the new direction
393 particleChangeForGamma->ProposeMomentumDirection(zPrikeVers.unit()) ;
394
395 // Update the energy which does not change here
396 particleChangeForGamma->SetProposedKineticEnergy(electronEnergy0);
397 }
398}
@ fStopAndKill
Hep3Vector orthogonal() const
Hep3Vector cross(const Hep3Vector &) const
G4double RandomizeCosTheta(G4double k, const G4String &materialName)
RandomizeCosTheta.
const G4ThreeVector & GetMomentumDirection() const
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

References CLHEP::Hep3Vector::cross(), fKillBelowEnergy, fStopAndKill, G4cout, G4endl, G4UniformRand, G4VDNAModel::GetHighELimit(), G4DynamicParticle::GetKineticEnergy(), G4VDNAModel::GetLowELimit(), G4DynamicParticle::GetMomentumDirection(), G4DynamicParticle::GetParticleDefinition(), G4ParticleDefinition::GetParticleName(), CLHEP::Hep3Vector::orthogonal(), pi, G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChangeForGamma::ProposeMomentumDirection(), G4VParticleChange::ProposeTrackStatus(), RandomizeCosTheta(), G4ParticleChangeForGamma::SetProposedKineticEnergy(), CLHEP::Hep3Vector::unit(), and verboseLevel.

◆ 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 Initialise(), G4DNADummyModel::Initialise(), G4DNAPTBExcitationModel::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 Initialise(), G4DNADummyModel::Initialise(), G4DNAPTBExcitationModel::Initialise(), and G4DNAPTBIonisationModel::Initialise().

◆ Theta()

G4double G4DNAPTBElasticModel::Theta ( G4ParticleDefinition fParticleDefinition,
G4double  k,
G4double  integrDiff,
const G4String materialName 
)
private

Theta To return an angular theta value from the differential file. This method uses interpolations to calculate the theta value.

Parameters
fParticleDefinition
k
integrDiff
materialName
Returns
a theta value

Definition at line 402 of file G4DNAPTBElasticModel.cc.

404{
405 G4double theta = 0.;
406 G4double valueT1 = 0;
407 G4double valueT2 = 0;
408 G4double valueE21 = 0;
409 G4double valueE22 = 0;
410 G4double valueE12 = 0;
411 G4double valueE11 = 0;
412 G4double xs11 = 0;
413 G4double xs12 = 0;
414 G4double xs21 = 0;
415 G4double xs22 = 0;
416 G4String particleName = particleDefinition->GetParticleName();
417
418 if (particleDefinition == G4Electron::ElectronDefinition())
419 {
420 std::vector<double>::iterator t2 = std::upper_bound(tValuesVec[materialName][particleName].begin(),tValuesVec[materialName][particleName].end(), k);
421 std::vector<double>::iterator t1 = t2-1;
422
423 std::vector<double>::iterator e12 = std::upper_bound(eValuesVect[materialName][particleName][(*t1)].begin(),eValuesVect[materialName][particleName][(*t1)].end(), integrDiff);
424 std::vector<double>::iterator e11 = e12-1;
425
426 std::vector<double>::iterator e22 = std::upper_bound(eValuesVect[materialName][particleName][(*t2)].begin(),eValuesVect[materialName][particleName][(*t2)].end(), integrDiff);
427 std::vector<double>::iterator e21 = e22-1;
428
429 valueT1 =*t1;
430 valueT2 =*t2;
431 valueE21 =*e21;
432 valueE22 =*e22;
433 valueE12 =*e12;
434 valueE11 =*e11;
435
436 xs11 = diffCrossSectionData[materialName][particleName][valueT1][valueE11];
437 xs12 = diffCrossSectionData[materialName][particleName][valueT1][valueE12];
438 xs21 = diffCrossSectionData[materialName][particleName][valueT2][valueE21];
439 xs22 = diffCrossSectionData[materialName][particleName][valueT2][valueE22];
440 }
441
442 if (xs11==0 && xs12==0 && xs21==0 && xs22==0) return (0.);
443
444 theta = QuadInterpolator ( valueE11, valueE12,
445 valueE21, valueE22,
446 xs11, xs12,
447 xs21, xs22,
448 valueT1, valueT2,
449 k, integrDiff );
450
451 return theta;
452}
G4double QuadInterpolator(G4double e11, G4double e12, G4double e21, G4double e22, G4double x11, G4double x12, G4double x21, G4double x22, G4double t1, G4double t2, G4double t, G4double e)
QuadInterpolator.

References diffCrossSectionData, G4Electron::ElectronDefinition(), eValuesVect, G4ParticleDefinition::GetParticleName(), QuadInterpolator(), and tValuesVec.

Referenced by RandomizeCosTheta().

Field Documentation

◆ diffCrossSectionData

TriDimensionMap G4DNAPTBElasticModel::diffCrossSectionData
private

A map: [materialName][particleName]=DiffCrossSectionTable.

Definition at line 126 of file G4DNAPTBElasticModel.hh.

Referenced by ReadDiffCSFile(), and Theta().

◆ eValuesVect

VecMap G4DNAPTBElasticModel::eValuesVect
private

map with vectors containing all the output energy (E) of the differential file

Definition at line 129 of file G4DNAPTBElasticModel.hh.

Referenced by ReadDiffCSFile(), and Theta().

◆ 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().

◆ fKillBelowEnergy

G4double G4DNAPTBElasticModel::fKillBelowEnergy
private

energy kill limit

Definition at line 123 of file G4DNAPTBElasticModel.hh.

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

◆ 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().

◆ killBelowEnergyTable

std::map<G4String, double > G4DNAPTBElasticModel::killBelowEnergyTable
private

map to save the different energy kill limits for the materials

Definition at line 122 of file G4DNAPTBElasticModel.hh.

◆ tValuesVec

std::map<G4String, std::map<G4String, std::vector<double> > > G4DNAPTBElasticModel::tValuesVec
private

map with vectors containing all the incident (T) energy of the differential file

Definition at line 130 of file G4DNAPTBElasticModel.hh.

Referenced by ReadDiffCSFile(), and Theta().

◆ verboseLevel

G4int G4DNAPTBElasticModel::verboseLevel
private

verbose level

Definition at line 121 of file G4DNAPTBElasticModel.hh.

Referenced by CrossSectionPerVolume(), G4DNAPTBElasticModel(), Initialise(), and SampleSecondaries().


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