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

The G4DNAPTBIonisationModel class Implements the PTB ionisation model. More...

#include <G4DNAPTBIonisationModel.hh>

Inheritance diagram for G4DNAPTBIonisationModel:
G4VDNAModel

Public Member Functions

virtual G4double CrossSectionPerVolume (const G4Material *material, const G4String &materialName, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
 CrossSectionPerVolume Mandatory for every model the CrossSectionPerVolume method is in charge of returning the cross section value corresponding to the material, particle and energy current values. More...
 
 G4DNAPTBIonisationModel (const G4String &applyToMaterial="all", const G4ParticleDefinition *p=0, const G4String &nam="DNAPTBIonisationModel", const G4bool isAuger=true)
 G4DNAPTBIonisationModel 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 Method called once at the beginning of the simulation. It is used to setup the list of the materials managed by the model and the energy limits. All the materials are setup but only a part of them can be activated by the user through the constructor. 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 SampleSecondaries will be called. The method sets the characteristics of the particles implied with the physical process after the ModelInterface (energy, momentum...). This method is mandatory for every model. 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 ~G4DNAPTBIonisationModel ()
 ~G4DNAPTBIonisationModel 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, std::map< double, double > > > > > TriDimensionMap
 
typedef std::map< G4String, std::map< G4String, std::map< double, std::vector< double > > > > VecMap
 
typedef std::map< G4String, std::map< G4String, std::map< double, std::map< double, std::vector< double > > > > > VecMapWithShell
 

Private Member Functions

double DifferentialCrossSection (G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell, const G4String &materialName)
 
 G4DNAPTBIonisationModel (const G4DNAPTBIonisationModel &)
 
G4double LogLogInterpolate (G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
 LogLogInterpolate. More...
 
G4DNAPTBIonisationModeloperator= (const G4DNAPTBIonisationModel &right)
 
G4double QuadInterpolator (G4double e11, G4double e12, G4double e21, G4double e22, G4double xs11, G4double xs12, G4double xs21, G4double xs22, G4double t1, G4double t2, G4double t, G4double e)
 QuadInterpolator. More...
 
void RandomizeEjectedElectronDirection (G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4double outgoingParticleEnergy, G4double &cosTheta, G4double &phi)
 RandomizeEjectedElectronDirection Method to calculate the ejected electron direction. More...
 
G4double RandomizeEjectedElectronEnergy (G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell, const G4String &materialName)
 
G4double RandomizeEjectedElectronEnergyFromCumulated (G4ParticleDefinition *particleDefinition, G4double k, G4int shell, const G4String &materialName)
 RandomizeEjectedElectronEnergyFromCumulated Uses the cumulated tables to find the energy of the ejected particle (electron) More...
 
void ReadDiffCSFile (const G4String &materialName, const G4String &particleName, const G4String &file, const G4double scaleFactor)
 ReadDiffCSFile Method to read the differential cross section files. More...
 

Private Attributes

TriDimensionMap diffCrossSectionData
 
G4DNAPTBAugerModelfDNAPTBAugerModel
 PTB Auger model instanciated in the constructor and deleted in the destructor of the class. More...
 
VecMap fEMapWithVector
 
TriDimensionMap fEnergySecondaryData
 
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...
 
VecMapWithShell fProbaShellMap
 
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, std::map< G4String, std::vector< double > > > fTMapWithVec
 
G4DNAPTBIonisationStructure ptbStructure
 
G4int verboseLevel
 verbose level More...
 

Detailed Description

The G4DNAPTBIonisationModel class Implements the PTB ionisation model.

Definition at line 53 of file G4DNAPTBIonisationModel.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, std::map<double, double> > > > > G4DNAPTBIonisationModel::TriDimensionMap
private

Definition at line 130 of file G4DNAPTBIonisationModel.hh.

◆ VecMap

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

Definition at line 134 of file G4DNAPTBIonisationModel.hh.

◆ VecMapWithShell

typedef std::map<G4String, std::map<G4String, std::map<double, std::map<double, std::vector<double> > > > > G4DNAPTBIonisationModel::VecMapWithShell
private

Definition at line 136 of file G4DNAPTBIonisationModel.hh.

Constructor & Destructor Documentation

◆ G4DNAPTBIonisationModel() [1/2]

G4DNAPTBIonisationModel::G4DNAPTBIonisationModel ( const G4String applyToMaterial = "all",
const G4ParticleDefinition p = 0,
const G4String nam = "DNAPTBIonisationModel",
const G4bool  isAuger = true 
)

G4DNAPTBIonisationModel Constructor.

Parameters
applyToMaterial
p
nam
isAuger

Definition at line 38 of file G4DNAPTBIonisationModel.cc.

41 : G4VDNAModel(nam, applyToMaterial)
42{
43 verboseLevel= 0;
44 // Verbosity scale:
45 // 0 = nothing
46 // 1 = warning for energy non-conservation
47 // 2 = details of energy budget
48 // 3 = calculation of cross sections, file openings, sampling of atoms
49 // 4 = entering in methods
50
51 if( verboseLevel>0 )
52 {
53 G4cout << "PTB ionisation model is constructed " << G4endl;
54 }
55
56 if(isAuger)
57 {
58 // create the PTB Auger model
59 fDNAPTBAugerModel = new G4DNAPTBAugerModel("e-_G4DNAPTBAugerModel");
60 }
61 else
62 {
63 // no PTB Auger model
65 }
66}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
The G4DNAPTBAugerModel class Implement the PTB Auger model.
G4DNAPTBAugerModel * fDNAPTBAugerModel
PTB Auger model instanciated in the constructor and deleted in the destructor of the class.
G4VDNAModel(const G4String &nam, const G4String &applyToMaterial)
G4VDNAModel Constructeur of the G4VDNAModel class.
Definition: G4VDNAModel.cc:35

References fDNAPTBAugerModel, G4cout, G4endl, and verboseLevel.

◆ ~G4DNAPTBIonisationModel()

G4DNAPTBIonisationModel::~G4DNAPTBIonisationModel ( )
virtual

~G4DNAPTBIonisationModel Destructor

Definition at line 70 of file G4DNAPTBIonisationModel.cc.

71{
72 // To delete the DNAPTBAugerModel created at initialisation of the ionisation class
74}

References fDNAPTBAugerModel.

◆ G4DNAPTBIonisationModel() [2/2]

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

CrossSectionPerVolume Mandatory for every model the CrossSectionPerVolume method is in charge of returning the cross section value corresponding to the material, particle and energy current values.

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

Implements G4VDNAModel.

Definition at line 250 of file G4DNAPTBIonisationModel.cc.

256{
257 if(verboseLevel > 3)
258 G4cout << "Calling CrossSectionPerVolume() of G4DNAPTBIonisationModel" << G4endl;
259
260 // initialise the cross section value (output value)
261 G4double sigma(0);
262
263 // Get the current particle name
264 const G4String& particleName = p->GetParticleName();
265
266 // Set the low and high energy limits
267 G4double lowLim = GetLowELimit(materialName, particleName);
268 G4double highLim = GetHighELimit(materialName, particleName);
269
270 // Check that we are in the correct energy range
271 if (ekin >= lowLim && ekin < highLim)
272 {
273 // Get the map with all the model data tables
274 TableMapData* tableData = GetTableData();
275
276 // Retrieve the cross section value for the current material, particle and energy values
277 sigma = (*tableData)[materialName][particleName]->FindValue(ekin);
278
279 if (verboseLevel > 2)
280 {
281 G4cout << "__________________________________" << G4endl;
282 G4cout << "°°° G4DNAPTBIonisationModel - XS INFO START" << G4endl;
283 G4cout << "°°° Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
284 G4cout << "°°° Cross section per "<< materialName <<" molecule (cm^2)=" << sigma/cm/cm << G4endl;
285 G4cout << "°°° G4DNAPTBIonisationModel - XS INFO END" << G4endl;
286 }
287 }
288
289 // Return the cross section value
290 return sigma;
291}
static constexpr double eV
Definition: G4SIunits.hh:201
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

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

◆ DifferentialCrossSection()

double G4DNAPTBIonisationModel::DifferentialCrossSection ( G4ParticleDefinition aParticleDefinition,
G4double  k,
G4double  energyTransfer,
G4int  shell,
const G4String materialName 
)
private

Definition at line 673 of file G4DNAPTBIonisationModel.cc.

678{
679 G4double sigma = 0.;
680 const G4String& particleName = particleDefinition->GetParticleName();
681
682 G4double shellEnergy (ptbStructure.IonisationEnergy(ionizationLevelIndex, materialName));
683 G4double kSE (energyTransfer-shellEnergy);
684
685 if (energyTransfer >= shellEnergy)
686 {
687 G4double valueT1 = 0;
688 G4double valueT2 = 0;
689 G4double valueE21 = 0;
690 G4double valueE22 = 0;
691 G4double valueE12 = 0;
692 G4double valueE11 = 0;
693
694 G4double xs11 = 0;
695 G4double xs12 = 0;
696 G4double xs21 = 0;
697 G4double xs22 = 0;
698
699 if (particleDefinition == G4Electron::ElectronDefinition())
700 {
701 // k should be in eV and energy transfer eV also
702 std::vector<double>::iterator t2 = std::upper_bound(fTMapWithVec[materialName][particleName].begin(),fTMapWithVec[materialName][particleName].end(), k);
703 std::vector<double>::iterator t1 = t2-1;
704
705 // SI : the following condition avoids situations where energyTransfer >last vector element
706 if (kSE <= fEMapWithVector[materialName][particleName][(*t1)].back() && kSE <= fEMapWithVector[materialName][particleName][(*t2)].back() )
707 {
708 std::vector<double>::iterator e12 = std::upper_bound(fEMapWithVector[materialName][particleName][(*t1)].begin(),fEMapWithVector[materialName][particleName][(*t1)].end(), kSE);
709 std::vector<double>::iterator e11 = e12-1;
710
711 std::vector<double>::iterator e22 = std::upper_bound(fEMapWithVector[materialName][particleName][(*t2)].begin(),fEMapWithVector[materialName][particleName][(*t2)].end(), kSE);
712 std::vector<double>::iterator e21 = e22-1;
713
714 valueT1 =*t1;
715 valueT2 =*t2;
716 valueE21 =*e21;
717 valueE22 =*e22;
718 valueE12 =*e12;
719 valueE11 =*e11;
720
721 xs11 = diffCrossSectionData[materialName][particleName][ionizationLevelIndex][valueT1][valueE11];
722 xs12 = diffCrossSectionData[materialName][particleName][ionizationLevelIndex][valueT1][valueE12];
723 xs21 = diffCrossSectionData[materialName][particleName][ionizationLevelIndex][valueT2][valueE21];
724 xs22 = diffCrossSectionData[materialName][particleName][ionizationLevelIndex][valueT2][valueE22];
725 }
726 }
727
728 if (particleDefinition == G4Proton::ProtonDefinition())
729 {
730 // k should be in eV and energy transfer eV also
731 std::vector<double>::iterator t2 = std::upper_bound(fTMapWithVec[materialName][particleName].begin(),fTMapWithVec[materialName][particleName].end(), k);
732 std::vector<double>::iterator t1 = t2-1;
733
734 std::vector<double>::iterator e12 = std::upper_bound(fEMapWithVector[materialName][particleName][(*t1)].begin(),fEMapWithVector[materialName][particleName][(*t1)].end(), kSE);
735 std::vector<double>::iterator e11 = e12-1;
736
737 std::vector<double>::iterator e22 = std::upper_bound(fEMapWithVector[materialName][particleName][(*t2)].begin(),fEMapWithVector[materialName][particleName][(*t2)].end(), kSE);
738 std::vector<double>::iterator e21 = e22-1;
739
740 valueT1 =*t1;
741 valueT2 =*t2;
742 valueE21 =*e21;
743 valueE22 =*e22;
744 valueE12 =*e12;
745 valueE11 =*e11;
746
747 xs11 = diffCrossSectionData[materialName][particleName][ionizationLevelIndex][valueT1][valueE11];
748 xs12 = diffCrossSectionData[materialName][particleName][ionizationLevelIndex][valueT1][valueE12];
749 xs21 = diffCrossSectionData[materialName][particleName][ionizationLevelIndex][valueT2][valueE21];
750 xs22 = diffCrossSectionData[materialName][particleName][ionizationLevelIndex][valueT2][valueE22];
751 }
752
753 G4double xsProduct = xs11 * xs12 * xs21 * xs22;
754
755 if (xsProduct != 0.)
756 {
757 sigma = QuadInterpolator(valueE11, valueE12,
758 valueE21, valueE22,
759 xs11, xs12,
760 xs21, xs22,
761 valueT1, valueT2,
762 k, kSE);
763 }
764 }
765
766
767 return sigma;
768}
std::map< G4String, std::map< G4String, std::vector< double > > > fTMapWithVec
G4double QuadInterpolator(G4double e11, G4double e12, G4double e21, G4double e22, G4double xs11, G4double xs12, G4double xs21, G4double xs22, G4double t1, G4double t2, G4double t, G4double e)
QuadInterpolator.
G4DNAPTBIonisationStructure ptbStructure
G4double IonisationEnergy(G4int level, const G4String &materialName)
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:88
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:87

References diffCrossSectionData, G4Electron::ElectronDefinition(), fEMapWithVector, fTMapWithVec, G4ParticleDefinition::GetParticleName(), G4DNAPTBIonisationStructure::IonisationEnergy(), G4Proton::ProtonDefinition(), ptbStructure, and QuadInterpolator().

Referenced by RandomizeEjectedElectronEnergy().

◆ 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(), G4DNAPTBExcitationModel::CrossSectionPerVolume(), CrossSectionPerVolume(), G4DNAPTBElasticModel::SampleSecondaries(), and 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(), G4DNAPTBExcitationModel::CrossSectionPerVolume(), CrossSectionPerVolume(), G4DNAPTBElasticModel::SampleSecondaries(), and 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(), G4DNAPTBExcitationModel::CrossSectionPerVolume(), CrossSectionPerVolume(), and G4VDNAModel::RandomSelectShell().

◆ Initialise()

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

Initialise Method called once at the beginning of the simulation. It is used to setup the list of the materials managed by the model and the energy limits. All the materials are setup but only a part of them can be activated by the user through the constructor.

Implements G4VDNAModel.

Definition at line 78 of file G4DNAPTBIonisationModel.cc.

80{
81
82 if (verboseLevel > 3)
83 G4cout << "Calling G4DNAPTBIonisationModel::Initialise()" << G4endl;
84
85 G4double scaleFactor = 1e-16 * cm*cm;
86 G4double scaleFactorBorn = (1.e-22 / 3.343) * m*m;
87
90
91 //*******************************************************
92 // Cross section data
93 //*******************************************************
94
95 if(particle == electronDef)
96 {
97 G4String particleName = particle->GetParticleName();
98
99 // Raw materials
100 //
102 particleName,
103 "dna/sigma_ionisation_e-_PTB_THF",
104 "dna/sigmadiff_cumulated_ionisation_e-_PTB_THF",
105 scaleFactor);
106 SetLowELimit("THF", particleName, 12.*eV);
107 SetHighELimit("THF", particleName, 1.*keV);
108
110 particleName,
111 "dna/sigma_ionisation_e-_PTB_PY",
112 "dna/sigmadiff_cumulated_ionisation_e-_PTB_PY",
113 scaleFactor);
114 SetLowELimit("PY", particleName, 12.*eV);
115 SetHighELimit("PY", particleName, 1.*keV);
116
118 particleName,
119 "dna/sigma_ionisation_e-_PTB_PU",
120 "dna/sigmadiff_cumulated_ionisation_e-_PTB_PU",
121 scaleFactor);
122 SetLowELimit("PU", particleName, 12.*eV);
123 SetHighELimit("PU", particleName, 1.*keV);
124
126 particleName,
127 "dna/sigma_ionisation_e-_PTB_TMP",
128 "dna/sigmadiff_cumulated_ionisation_e-_PTB_TMP",
129 scaleFactor);
130 SetLowELimit("TMP", particleName, 12.*eV);
131 SetHighELimit("TMP", particleName, 1.*keV);
132
133 AddCrossSectionData("G4_WATER",
134 particleName,
135 "dna/sigma_ionisation_e_born",
136 "dna/sigmadiff_ionisation_e_born",
137 scaleFactorBorn);
138 SetLowELimit("G4_WATER", particleName, 12.*eV);
139 SetHighELimit("G4_WATER", particleName, 1.*keV);
140
141 // DNA materials
142 //
143 AddCrossSectionData("backbone_THF",
144 particleName,
145 "dna/sigma_ionisation_e-_PTB_THF",
146 "dna/sigmadiff_cumulated_ionisation_e-_PTB_THF",
147 scaleFactor*33./30);
148 SetLowELimit("backbone_THF", particleName, 12.*eV);
149 SetHighELimit("backbone_THF", particleName, 1.*keV);
150
151 AddCrossSectionData("cytosine_PY",
152 particleName,
153 "dna/sigma_ionisation_e-_PTB_PY",
154 "dna/sigmadiff_cumulated_ionisation_e-_PTB_PY",
155 scaleFactor*42./30);
156 SetLowELimit("cytosine_PY", particleName, 12.*eV);
157 SetHighELimit("cytosine_PY", particleName, 1.*keV);
158
159 AddCrossSectionData("thymine_PY",
160 particleName,
161 "dna/sigma_ionisation_e-_PTB_PY",
162 "dna/sigmadiff_cumulated_ionisation_e-_PTB_PY",
163 scaleFactor*48./30);
164 SetLowELimit("thymine_PY", particleName, 12.*eV);
165 SetHighELimit("thymine_PY", particleName, 1.*keV);
166
167 AddCrossSectionData("adenine_PU",
168 particleName,
169 "dna/sigma_ionisation_e-_PTB_PU",
170 "dna/sigmadiff_cumulated_ionisation_e-_PTB_PU",
171 scaleFactor*50./44);
172 SetLowELimit("adenine_PU", particleName, 12.*eV);
173 SetHighELimit("adenine_PU", particleName, 1.*keV);
174
175 AddCrossSectionData("guanine_PU",
176 particleName,
177 "dna/sigma_ionisation_e-_PTB_PU",
178 "dna/sigmadiff_cumulated_ionisation_e-_PTB_PU",
179 scaleFactor*56./44);
180 SetLowELimit("guanine_PU", particleName, 12.*eV);
181 SetHighELimit("guanine_PU", particleName, 1.*keV);
182
183 AddCrossSectionData("backbone_TMP",
184 particleName,
185 "dna/sigma_ionisation_e-_PTB_TMP",
186 "dna/sigmadiff_cumulated_ionisation_e-_PTB_TMP",
187 scaleFactor*33./50);
188 SetLowELimit("backbone_TMP", particleName, 12.*eV);
189 SetHighELimit("backbone_TMP", particleName, 1.*keV);
190 }
191
192 else if (particle == protonDef)
193 {
194 G4String particleName = particle->GetParticleName();
195
196 // Raw materials
197 //
199 particleName,
200 "dna/sigma_ionisation_p_HKS_THF",
201 "dna/sigmadiff_cumulated_ionisation_p_PTB_THF",
202 scaleFactor);
203 SetLowELimit("THF", particleName, 70.*keV);
204 SetHighELimit("THF", particleName, 10.*MeV);
205
206
208 particleName,
209 "dna/sigma_ionisation_p_HKS_PY",
210 "dna/sigmadiff_cumulated_ionisation_p_PTB_PY",
211 scaleFactor);
212 SetLowELimit("PY", particleName, 70.*keV);
213 SetHighELimit("PY", particleName, 10.*MeV);
214
215 /*
216 AddCrossSectionData("PU",
217 particleName,
218 "dna/sigma_ionisation_e-_PTB_PU",
219 "dna/sigmadiff_cumulated_ionisation_e-_PTB_PU",
220 scaleFactor);
221 SetLowELimit("PU", particleName2, 70.*keV);
222 SetHighELimit("PU", particleName2, 10.*keV);
223*/
224
226 particleName,
227 "dna/sigma_ionisation_p_HKS_TMP",
228 "dna/sigmadiff_cumulated_ionisation_p_PTB_TMP",
229 scaleFactor);
230 SetLowELimit("TMP", particleName, 70.*keV);
231 SetHighELimit("TMP", particleName, 10.*MeV);
232 }
233
234 // *******************************************************
235 // deal with composite materials
236 // *******************************************************
237
239
240 // *******************************************************
241 // Verbose
242 // *******************************************************
243
244 // initialise DNAPTBAugerModel
246}
static constexpr double m
Definition: G4SIunits.hh:109
static constexpr double keV
Definition: G4SIunits.hh:202
static constexpr double MeV
Definition: G4SIunits.hh:200
virtual void Initialise()
Initialise Set the verbose value.
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, fDNAPTBAugerModel, G4cout, G4endl, G4ParticleDefinition::GetParticleName(), G4DNAPTBAugerModel::Initialise(), keV, G4VDNAModel::LoadCrossSectionData(), m, MeV, G4Proton::ProtonDefinition(), 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(), G4DNAPTBExcitationModel::Initialise(), and Initialise().

◆ LogLogInterpolate()

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

LogLogInterpolate.

Parameters
e1
e2
e
xs1
xs2
Returns
the interpolate value

Definition at line 939 of file G4DNAPTBIonisationModel.cc.

944{
945 G4double value (0);
946
947 // Switch to log-lin interpolation for faster code
948
949 if ((e2-e1)!=0 && xs1 !=0 && xs2 !=0)
950 {
951 G4double d1 = std::log10(xs1);
952 G4double d2 = std::log10(xs2);
953 value = std::pow(10.,(d1 + (d2 - d1)*(e - e1)/ (e2 - e1)) );
954 }
955
956 // Switch to lin-lin interpolation for faster code
957 // in case one of xs1 or xs2 (=cum proba) value is zero
958
959 if ((e2-e1)!=0 && (xs1 ==0 || xs2 ==0))
960 {
961 G4double d1 = xs1;
962 G4double d2 = xs2;
963 value = (d1 + (d2 - d1)*(e - e1)/ (e2 - e1));
964 }
965
966 return value;
967}
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().

◆ operator=()

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

◆ QuadInterpolator()

G4double G4DNAPTBIonisationModel::QuadInterpolator ( G4double  e11,
G4double  e12,
G4double  e21,
G4double  e22,
G4double  xs11,
G4double  xs12,
G4double  xs21,
G4double  xs22,
G4double  t1,
G4double  t2,
G4double  t,
G4double  e 
)
private

QuadInterpolator.

Parameters
e11
e12
e21
e22
xs11
xs12
xs21
xs22
t1
t2
t
e
Returns
the interpolated value

Definition at line 971 of file G4DNAPTBIonisationModel.cc.

977{
978 G4double interpolatedvalue1 (-1);
979 if(xs11!=xs12) interpolatedvalue1 = LogLogInterpolate(e11, e12, e, xs11, xs12);
980 else interpolatedvalue1 = xs11;
981
982 G4double interpolatedvalue2 (-1);
983 if(xs21!=xs22) interpolatedvalue2 = LogLogInterpolate(e21, e22, e, xs21, xs22);
984 else interpolatedvalue2 = xs21;
985
986 G4double value (-1);
987 if(interpolatedvalue1!=interpolatedvalue2) value = LogLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
988 else value = interpolatedvalue1;
989
990 return value;
991
992 // G4double interpolatedvalue1 = LogLogInterpolate(e11, e12, e, xs11, xs12);
993 // G4double interpolatedvalue2 = LogLogInterpolate(e21, e22, e, xs21, xs22);
994 // G4double value = LogLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
995 // return value;
996}
G4double LogLogInterpolate(G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2)
LogLogInterpolate.

References LogLogInterpolate().

Referenced by DifferentialCrossSection(), and RandomizeEjectedElectronEnergyFromCumulated().

◆ RandomizeEjectedElectronDirection()

void G4DNAPTBIonisationModel::RandomizeEjectedElectronDirection ( G4ParticleDefinition aParticleDefinition,
G4double  incomingParticleEnergy,
G4double  outgoingParticleEnergy,
G4double cosTheta,
G4double phi 
)
private

RandomizeEjectedElectronDirection Method to calculate the ejected electron direction.

Parameters
aParticleDefinition
incomingParticleEnergy
outgoingParticleEnergy
cosTheta
phi

Definition at line 634 of file G4DNAPTBIonisationModel.cc.

639{
640 if (particleDefinition == G4Electron::ElectronDefinition())
641 {
642 phi = twopi * G4UniformRand();
643 if (secKinetic < 50.*eV) cosTheta = (2.*G4UniformRand())-1.;
644 else if (secKinetic <= 200.*eV)
645 {
646 if (G4UniformRand() <= 0.1) cosTheta = (2.*G4UniformRand())-1.;
647 else cosTheta = G4UniformRand()*(std::sqrt(2.)/2);
648 }
649 else
650 {
651 G4double sin2O = (1.-secKinetic/k) / (1.+secKinetic/(2.*electron_mass_c2));
652 cosTheta = std::sqrt(1.-sin2O);
653 }
654 }
655
656 else if (particleDefinition == G4Proton::ProtonDefinition())
657 {
658 G4double maxSecKinetic = 4.* (electron_mass_c2 / proton_mass_c2) * k;
659 phi = twopi * G4UniformRand();
660
661 // cosTheta = std::sqrt(secKinetic / maxSecKinetic);
662
663 // Restriction below 100 eV from Emfietzoglou (2000)
664
665 if (secKinetic>100*eV) cosTheta = std::sqrt(secKinetic / maxSecKinetic);
666 else cosTheta = (2.*G4UniformRand())-1.;
667
668 }
669}
static constexpr double twopi
Definition: G4SIunits.hh:56
#define G4UniformRand()
Definition: Randomize.hh:52
float electron_mass_c2
Definition: hepunit.py:273
float proton_mass_c2
Definition: hepunit.py:274

References source.hepunit::electron_mass_c2, G4Electron::ElectronDefinition(), eV, G4UniformRand, source.hepunit::proton_mass_c2, G4Proton::ProtonDefinition(), and twopi.

Referenced by SampleSecondaries().

◆ RandomizeEjectedElectronEnergy()

G4double G4DNAPTBIonisationModel::RandomizeEjectedElectronEnergy ( G4ParticleDefinition aParticleDefinition,
G4double  incomingParticleEnergy,
G4int  shell,
const G4String materialName 
)
private

Definition at line 537 of file G4DNAPTBIonisationModel.cc.

539{
540 if (particleDefinition == G4Electron::ElectronDefinition())
541 {
542 //G4double Tcut=25.0E-6;
543 G4double maximumEnergyTransfer=0.;
544 if ((k+ptbStructure.IonisationEnergy(shell, materialName))/2. > k) maximumEnergyTransfer=k;
545 else maximumEnergyTransfer = (k+ptbStructure.IonisationEnergy(shell,materialName))/2.;
546
547 // SI : original method
548 /*
549 G4double crossSectionMaximum = 0.;
550 for(G4double value=waterStructure.IonisationEnergy(shell); value<=maximumEnergyTransfer; value+=0.1*eV)
551 {
552 G4double differentialCrossSection = DifferentialCrossSection(particleDefinition, k/eV, value/eV, shell);
553 if(differentialCrossSection >= crossSectionMaximum) crossSectionMaximum = differentialCrossSection;
554 }
555 */
556
557
558 // SI : alternative method
559
560 //if (k > Tcut)
561 //{
562 G4double crossSectionMaximum = 0.;
563
564 G4double minEnergy = ptbStructure.IonisationEnergy(shell, materialName);
565 G4double maxEnergy = maximumEnergyTransfer;
566 G4int nEnergySteps = 50;
567 G4double value(minEnergy);
568 G4double stpEnergy(std::pow(maxEnergy/value, 1./static_cast<G4double>(nEnergySteps-1)));
569 G4int step(nEnergySteps);
570 while (step>0)
571 {
572 step--;
573 G4double differentialCrossSection = DifferentialCrossSection(particleDefinition, k/eV, value/eV, shell, materialName);
574 if(differentialCrossSection >= crossSectionMaximum) crossSectionMaximum = differentialCrossSection;
575 value *= stpEnergy;
576
577 }
578 //
579
580
581 G4double secondaryElectronKineticEnergy=0.;
582
583 do
584 {
585 secondaryElectronKineticEnergy = G4UniformRand() * (maximumEnergyTransfer-ptbStructure.IonisationEnergy(shell, materialName));
586
587 } while(G4UniformRand()*crossSectionMaximum >
588 DifferentialCrossSection(particleDefinition, k/eV,(secondaryElectronKineticEnergy+ptbStructure.IonisationEnergy(shell, materialName))/eV,shell, materialName));
589
590 return secondaryElectronKineticEnergy;
591
592 // }
593
594 // else if (k < Tcut)
595 // {
596
597 // G4double bindingEnergy = ptbStructure.IonisationEnergy(shell, materialName);
598 // G4double maxEnergy = ((k-bindingEnergy)/2.);
599
600 // G4double secondaryElectronKineticEnergy = G4UniformRand()*maxEnergy;
601 // return secondaryElectronKineticEnergy;
602 // }
603 }
604
605
606 else if (particleDefinition == G4Proton::ProtonDefinition())
607 {
608 G4double maximumKineticEnergyTransfer = 4.* (electron_mass_c2 / proton_mass_c2) * k;
609
610 G4double crossSectionMaximum = 0.;
611 for (G4double value = ptbStructure.IonisationEnergy(shell, materialName);
612 value<=4.*ptbStructure.IonisationEnergy(shell, materialName) ;
613 value+=0.1*eV)
614 {
615 G4double differentialCrossSection = DifferentialCrossSection(particleDefinition, k/eV, value/eV, shell, materialName);
616 if (differentialCrossSection >= crossSectionMaximum) crossSectionMaximum = differentialCrossSection;
617 }
618
619 G4double secondaryElectronKineticEnergy = 0.;
620 do
621 {
622 secondaryElectronKineticEnergy = G4UniformRand() * maximumKineticEnergyTransfer;
623 } while(G4UniformRand()*crossSectionMaximum >=
624 DifferentialCrossSection(particleDefinition, k/eV,(secondaryElectronKineticEnergy+ptbStructure.IonisationEnergy(shell, materialName))/eV,shell, materialName));
625
626 return secondaryElectronKineticEnergy;
627 }
628
629 return 0;
630}
int G4int
Definition: G4Types.hh:85
double DifferentialCrossSection(G4ParticleDefinition *aParticleDefinition, G4double k, G4double energyTransfer, G4int shell, const G4String &materialName)

References DifferentialCrossSection(), source.hepunit::electron_mass_c2, G4Electron::ElectronDefinition(), eV, G4UniformRand, G4DNAPTBIonisationStructure::IonisationEnergy(), source.hepunit::proton_mass_c2, G4Proton::ProtonDefinition(), and ptbStructure.

Referenced by SampleSecondaries().

◆ RandomizeEjectedElectronEnergyFromCumulated()

G4double G4DNAPTBIonisationModel::RandomizeEjectedElectronEnergyFromCumulated ( G4ParticleDefinition particleDefinition,
G4double  k,
G4int  shell,
const G4String materialName 
)
private

RandomizeEjectedElectronEnergyFromCumulated Uses the cumulated tables to find the energy of the ejected particle (electron)

Parameters
particleDefinition
k
shell
materialName
Returns
the ejected electron energy

Definition at line 772 of file G4DNAPTBIonisationModel.cc.

773{
774 // k should be in eV
775
776 // Schematic explanation.
777 // We will do an interpolation to get a final E value (ejected electron energy).
778 // 1/ We choose a random number between 0 and 1 (ie we select a cumulated cross section).
779 // 2/ We look for T_lower and T_upper.
780 // 3/ We look for the cumulated corresponding cross sections and their associated E values.
781 //
782 // T_low | CS_low_1 -> E_low_1
783 // | CS_low_2 -> E_low_2
784 // T_up | CS_up_1 -> E_up_1
785 // | CS_up_2 -> E_up_2
786 //
787 // 4/ We interpolate to get our E value.
788 //
789 // T_low | CS_low_1 -> E_low_1 -----
790 // | |----> E_low --
791 // | CS_low_2 -> E_low_2 ----- |
792 // | ---> E_final
793 // T_up | CS_up_1 -> E_up_1 ------- |
794 // | |----> E_up ---
795 // | CS_up_2 -> E_up_2 -------
796
797 // Initialize some values
798 //
799 G4double ejectedElectronEnergy = 0.;
800 G4double valueK1 = 0;
801 G4double valueK2 = 0;
802 G4double valueCumulCS21 = 0;
803 G4double valueCumulCS22 = 0;
804 G4double valueCumulCS12 = 0;
805 G4double valueCumulCS11 = 0;
806 G4double secElecE11 = 0;
807 G4double secElecE12 = 0;
808 G4double secElecE21 = 0;
809 G4double secElecE22 = 0;
810 G4String particleName = particleDefinition->GetParticleName();
811
812 // ***************************************************************************
813 // Get a random number between 0 and 1 to compare with the cumulated CS
814 // ***************************************************************************
815 //
816 // It will allow us to choose an ejected electron energy with respect to the CS.
817 G4double random = G4UniformRand();
818
819 // **********************************************
820 // Take the input from the data tables
821 // **********************************************
822
823 // Cumulated tables are like this: T E cumulatedCS1 cumulatedCS2 cumulatedCS3
824 // We have two sets of loaded data: fTMapWithVec which contains data about T (incident particle energy)
825 // and fProbaShellMap which contains cumulated cross section data.
826 // Since we already have a specific T energy value which could not be explicitly in the table, we must interpolate all the values.
827
828 // First, we select the upper and lower T data values surrounding our T value (ie "k").
829 std::vector<double>::iterator k2 = std::upper_bound(fTMapWithVec[materialName][particleName].begin(),fTMapWithVec[materialName][particleName].end(), k);
830 std::vector<double>::iterator k1 = k2-1;
831
832 // Check if we have found a k2 value (0 if we did not found it).
833 // A missing k2 value can be caused by a energy to high for the data table,
834 // Ex : table done for 12*eV -> 1000*eV and k=2000*eV
835 // then k2 = 0 and k1 = max of the table.
836 // To detect this, we check that k1 is not superior to k2.
837 if(*k1 > *k2)
838 {
839 // Error
840 G4cerr<<"**************** Fatal error ******************"<<G4endl;
841 G4cerr<<"G4DNAPTBIonisationModel::RandomizeEjectedElectronEnergyFromCumulated"<<G4endl;
842 G4cerr<<"You have *k1 > *k2 with k1 "<<*k1<<" and k2 "<<*k2<<G4endl;
843 G4cerr<<"This may be because the energy of the incident particle is to high for the data table."<<G4endl;
844 G4cerr<<"Particle energy (eV): "<<k<<G4endl;
845 exit(EXIT_FAILURE);
846 }
847
848
849 // We have a random number and we select the cumulated cross section data values surrounding our random number.
850 // But we need to do that for each T value (ie two T values) previously selected.
851 //
852 // First one.
853 std::vector<double>::iterator cumulCS12 = std::upper_bound(fProbaShellMap[materialName][particleName][ionizationLevelIndex][(*k1)].begin(),
854 fProbaShellMap[materialName][particleName][ionizationLevelIndex][(*k1)].end(), random);
855 std::vector<double>::iterator cumulCS11 = cumulCS12-1;
856 // Second one.
857 std::vector<double>::iterator cumulCS22 = std::upper_bound(fProbaShellMap[materialName][particleName][ionizationLevelIndex][(*k2)].begin(),
858 fProbaShellMap[materialName][particleName][ionizationLevelIndex][(*k2)].end(), random);
859 std::vector<double>::iterator cumulCS21 = cumulCS22-1;
860
861 // Now that we have the "values" through pointers, we access them.
862 valueK1 = *k1;
863 valueK2 = *k2;
864 valueCumulCS11 = *cumulCS11;
865 valueCumulCS12 = *cumulCS12;
866 valueCumulCS21 = *cumulCS21;
867 valueCumulCS22 = *cumulCS22;
868
869 // *************************************************************
870 // Do the interpolation to get the ejected electron energy
871 // *************************************************************
872
873 // Here we will get four E values corresponding to our four cumulated cross section values previously selected.
874 // But we need to take into account a specific case: we have selected a shell by using the ionisation cross section table
875 // and, since we get two T values, we could have differential cross sections (or cumulated) equal to 0 for the lower T
876 // and not for the upper T. When looking for the cumulated cross section values which surround the selected random number (for the lower T),
877 // the upper_bound method will only found 0 values. Thus, the upper_bound method will return the last E value present in the table for the
878 // selected T. The last E value being the highest, we will later perform an interpolation between a high E value (for the lower T) and
879 // a small E value (for the upper T). This is inconsistent because if the cross section are equal to zero for the lower T then it
880 // means it is not possible to ionize and, thus, to have a secondary electron. But, in our situation, it is possible to ionize for the upper T
881 // AND for an interpolate T value between Tupper Tlower. That's why the final E value should be interpolate between 0 and the E value (upper T).
882 //
883 if(cumulCS12==fProbaShellMap[materialName][particleName][ionizationLevelIndex][(*k1)].end())
884 {
885 // Here we are in the special case and we force Elower1 and Elower2 to be equal at 0 for the interpolation.
886 secElecE11 = 0;
887 secElecE12 = 0;
888 secElecE21 = fEnergySecondaryData[materialName][particleName][ionizationLevelIndex][valueK2][valueCumulCS21];
889 secElecE22 = fEnergySecondaryData[materialName][particleName][ionizationLevelIndex][valueK2][valueCumulCS22];
890
891 valueCumulCS11 = 0;
892 valueCumulCS12 = 0;
893 }
894 else
895 {
896 // No special case, interpolation will happen as usual.
897 secElecE11 = fEnergySecondaryData[materialName][particleName][ionizationLevelIndex][valueK1][valueCumulCS11];
898 secElecE12 = fEnergySecondaryData[materialName][particleName][ionizationLevelIndex][valueK1][valueCumulCS12];
899 secElecE21 = fEnergySecondaryData[materialName][particleName][ionizationLevelIndex][valueK2][valueCumulCS21];
900 secElecE22 = fEnergySecondaryData[materialName][particleName][ionizationLevelIndex][valueK2][valueCumulCS22];
901 }
902
903 ejectedElectronEnergy = QuadInterpolator(valueCumulCS11, valueCumulCS12,
904 valueCumulCS21, valueCumulCS22,
905 secElecE11, secElecE12,
906 secElecE21, secElecE22,
907 valueK1, valueK2,
908 k, random);
909
910 // **********************************************
911 // Some tests for debugging
912 // **********************************************
913
914 G4double bindingEnergy (ptbStructure.IonisationEnergy(ionizationLevelIndex, materialName)/eV);
915 if(k-ejectedElectronEnergy-bindingEnergy<=0 || ejectedElectronEnergy<=0)
916 {
917 G4cout<<"k "<<k<<G4endl;
918 G4cout<<"material "<<materialName<<G4endl;
919 G4cout<<"secondaryKin "<<ejectedElectronEnergy<<G4endl;
920 G4cout<<"shell "<<ionizationLevelIndex<<G4endl;
921 G4cout<<"bindingEnergy "<<bindingEnergy<<G4endl;
922 G4cout<<"scatteredEnergy "<<k-ejectedElectronEnergy-bindingEnergy<<G4endl;
923 G4cout<<"rand "<<random<<G4endl;
924 G4cout<<"surrounding k values: valueK1 valueK2\n"<<valueK1<<" "<<valueK2<<G4endl;
925 G4cout<<"surrounding E values: secElecE11 secElecE12 secElecE21 secElecE22\n"
926 <<secElecE11<<" "<<secElecE12<<" "<<secElecE21<<" "<<secElecE22<<" "<<G4endl;
927 G4cout<<"surrounding cumulCS values: valueCumulCS11 valueCumulCS12 valueCumulCS21 valueCumulCS22\n"
928 <<valueCumulCS11<<" "<<valueCumulCS12<<" "<<valueCumulCS21<<" "<<valueCumulCS22<<" "<<G4endl;
929 G4cerr<<"*****************************"<<G4endl;
930 G4cerr<<"Fatal error, EXIT."<<G4endl;
931 exit(EXIT_FAILURE);
932 }
933
934 return ejectedElectronEnergy*eV;
935}
G4GLOB_DLL std::ostream G4cerr
G4double bindingEnergy(G4int A, G4int Z)
def exit()
Definition: g4zmq.py:99

References G4InuclSpecialFunctions::bindingEnergy(), eV, g4zmq::exit(), fEnergySecondaryData, fProbaShellMap, fTMapWithVec, G4cerr, G4cout, G4endl, G4UniformRand, G4ParticleDefinition::GetParticleName(), G4DNAPTBIonisationStructure::IonisationEnergy(), ptbStructure, and QuadInterpolator().

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

ReadDiffCSFile Method to read the differential cross section files.

Parameters
materialName
particleName
file
scaleFactor

Reimplemented from G4VDNAModel.

Definition at line 432 of file G4DNAPTBIonisationModel.cc.

436{
437 // To read and save the informations contained within the differential cross section files
438
439 // get the path of the G4LEDATA data folder
440 char *path = std::getenv("G4LEDATA");
441 // if it is not found then quit and print error message
442 if(!path)
443 {
444 G4Exception("G4DNAPTBIonisationModel::ReadAllDiffCSFiles","em0006",
445 FatalException,"G4LEDATA environment variable not set.");
446 return;
447 }
448
449 // build the fullFileName path of the data file
450 std::ostringstream fullFileName;
451 fullFileName << path <<"/"<< file<<".dat";
452
453 // open the data file
454 std::ifstream diffCrossSection (fullFileName.str().c_str());
455 // error if file is not there
456 std::stringstream endPath;
457 if (!diffCrossSection)
458 {
459 endPath << "Missing data file: "<<file;
460 G4Exception("G4DNAPTBIonisationModel::Initialise","em0003",
461 FatalException, endPath.str().c_str());
462 }
463
464 // load data from the file
465 fTMapWithVec[materialName][particleName].push_back(0.);
466
467 G4String line;
468
469 // read the file until we reach the end of file point
470 // fill fTMapWithVec, diffCrossSectionData, fEnergyTransferData, fProbaShellMap and fEMapWithVector
471 while(std::getline(diffCrossSection, line))
472 {
473 // check if the line is comment or empty
474 //
475 std::istringstream testIss(line);
477 testIss >> test;
478 // check first caracter to determine if following information is data or comments
479 if(test=="#")
480 {
481 // skip the line by beginning a new while loop.
482 continue;
483 }
484 // check if line is empty
485 else if(line.empty())
486 {
487 // skip the line by beginning a new while loop.
488 continue;
489 }
490 //
491 // end of the check
492
493 // transform the line into a iss
494 std::istringstream iss(line);
495
496 // Initialise the variables to be filled
497 double T;
498 double E;
499
500 // Filled T and E with the first two numbers of each file line
501 iss>>T>>E;
502
503 // Fill the fTMapWithVec container with all the different T values contained within the file.
504 // Duplicate must be avoided and this is the purpose of the if statement
505 if (T != fTMapWithVec[materialName][particleName].back()) fTMapWithVec[materialName][particleName].push_back(T);
506
507 // iterate on each shell of the corresponding material
508 for (int shell=0, eshell=ptbStructure.NumberOfLevels(materialName); shell<eshell; ++shell)
509 {
510 // map[material][particle][shell][T][E]=diffCrossSectionValue
511 // Fill the map with the informations of the input file
512 iss>>diffCrossSectionData[materialName][particleName][shell][T][E];
513
514 if(materialName!="G4_WATER")
515 {
516 // map[material][particle][shell][T][CS]=E
517 // Fill the map
518 fEnergySecondaryData[materialName][particleName][shell][T][diffCrossSectionData[materialName][particleName][shell][T][E] ]=E;
519
520 // map[material][particle][shell][T]=CS_vector
521 // Fill the vector within the map
522 fProbaShellMap[materialName][particleName][shell][T].push_back(diffCrossSectionData[materialName][particleName][shell][T][E]);
523 }
524 else
525 {
526 diffCrossSectionData[materialName][particleName][shell][T][E]*=scaleFactor;
527
528 fEMapWithVector[materialName][particleName][T].push_back(E);
529 }
530 }
531 }
532}
G4int NumberOfLevels(const G4String &materialName)
def test()
Definition: mcscore.py:117
Definition: test.py:1

References diffCrossSectionData, FatalException, fEMapWithVector, fEnergySecondaryData, geant4_check_module_cycles::file, fProbaShellMap, fTMapWithVec, G4Exception(), G4DNAPTBIonisationStructure::NumberOfLevels(), ptbStructure, and mcscore::test().

◆ SampleSecondaries()

void G4DNAPTBIonisationModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
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 SampleSecondaries will be called. The method sets the characteristics of the particles implied with the physical process after the ModelInterface (energy, momentum...). This method is mandatory for every model.

Parameters
materialName
particleChangeForGamma
tmin
tmax

Implements G4VDNAModel.

Definition at line 295 of file G4DNAPTBIonisationModel.cc.

302{
303 if (verboseLevel > 3)
304 G4cout << "Calling SampleSecondaries() of G4DNAPTBIonisationModel" << G4endl;
305
306 // Get the current particle energy
307 G4double k = aDynamicParticle->GetKineticEnergy();
308
309 // Get the current particle name
310 const G4String& particleName = aDynamicParticle->GetDefinition()->GetParticleName();
311
312 // Get the energy limits
313 G4double lowLim = GetLowELimit(materialName, particleName);
314 G4double highLim = GetHighELimit(materialName, particleName);
315
316 // Check if we are in the correct energy range
317 if (k >= lowLim && k < highLim)
318 {
319 G4ParticleMomentum primaryDirection = aDynamicParticle->GetMomentumDirection();
320 G4double particleMass = aDynamicParticle->GetDefinition()->GetPDGMass();
321 G4double totalEnergy = k + particleMass;
322 G4double pSquare = k * (totalEnergy + particleMass);
323 G4double totalMomentum = std::sqrt(pSquare);
324
325 // Get the ionisation shell from a random sampling
326 G4int ionizationShell = RandomSelectShell(k, particleName, materialName);
327
328 // Get the binding energy from the ptbStructure class
329 G4double bindingEnergy = ptbStructure.IonisationEnergy(ionizationShell, materialName);
330
331 // Initialize the secondary kinetic energy to a negative value.
332 G4double secondaryKinetic (-1000*eV);
333
334 if(materialName!="G4_WATER")
335 {
336 // Get the energy of the secondary particle
337 secondaryKinetic = RandomizeEjectedElectronEnergyFromCumulated(aDynamicParticle->GetDefinition(),k/eV,ionizationShell, materialName);
338 }
339 else
340 {
341 secondaryKinetic = RandomizeEjectedElectronEnergy(aDynamicParticle->GetDefinition(),k,ionizationShell, materialName);
342 }
343
344 if(secondaryKinetic<=0)
345 {
346 G4cout<<"Fatal error *************************************** "<<secondaryKinetic/eV<<G4endl;
347 G4cout<<"secondaryKinetic: "<<secondaryKinetic/eV<<G4endl;
348 G4cout<<"k: "<<k/eV<<G4endl;
349 G4cout<<"shell: "<<ionizationShell<<G4endl;
350 G4cout<<"material:"<<materialName<<G4endl;
351 exit(EXIT_FAILURE);
352 }
353
354 G4double cosTheta = 0.;
355 G4double phi = 0.;
356 RandomizeEjectedElectronDirection(aDynamicParticle->GetDefinition(), k, secondaryKinetic, cosTheta, phi);
357
358 G4double sinTheta = std::sqrt(1.-cosTheta*cosTheta);
359 G4double dirX = sinTheta*std::cos(phi);
360 G4double dirY = sinTheta*std::sin(phi);
361 G4double dirZ = cosTheta;
362 G4ThreeVector deltaDirection(dirX,dirY,dirZ);
363 deltaDirection.rotateUz(primaryDirection);
364
365 // The model is written only for electron and thus we want the change the direction of the incident electron
366 // after each ionization. However, if other particle are going to be introduced within this model the following should be added:
367 //
368 // Check if the particle is an electron
369 if(aDynamicParticle->GetDefinition() == G4Electron::ElectronDefinition() )
370 {
371 // If yes do the following code until next commented "else" statement
372
373 G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
374 G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
375 G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
376 G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
377 G4double finalMomentum = std::sqrt(finalPx*finalPx + finalPy*finalPy + finalPz*finalPz);
378 finalPx /= finalMomentum;
379 finalPy /= finalMomentum;
380 finalPz /= finalMomentum;
381
382 G4ThreeVector direction(finalPx,finalPy,finalPz);
383 if(direction.unit().getX()>1||direction.unit().getY()>1||direction.unit().getZ()>1)
384 {
385 G4cout<<"Fatal error ****************************"<<G4endl;
386 G4cout<<"direction problem "<<direction.unit()<<G4endl;
387 exit(EXIT_FAILURE);
388 }
389
390 // Give the new direction to the particle
391 particleChangeForGamma->ProposeMomentumDirection(direction.unit()) ;
392 }
393 // If the particle is not an electron
394 else particleChangeForGamma->ProposeMomentumDirection(primaryDirection) ;
395
396 // note that secondaryKinetic is the energy of the delta ray, not of all secondaries.
397 G4double scatteredEnergy = k-bindingEnergy-secondaryKinetic;
398
399 if(scatteredEnergy<=0)
400 {
401 G4cout<<"Fatal error ****************************"<<G4endl;
402 G4cout<<"k: "<<k/eV<<G4endl;
403 G4cout<<"secondaryKinetic: "<<secondaryKinetic/eV<<G4endl;
404 G4cout<<"shell: "<<ionizationShell<<G4endl;
405 G4cout<<"bindingEnergy: "<<bindingEnergy/eV<<G4endl;
406 G4cout<<"scatteredEnergy: "<<scatteredEnergy/eV<<G4endl;
407 G4cout<<"material: "<<materialName<<G4endl;
408 exit(EXIT_FAILURE);
409 }
410
411 // Set the new energy of the particle
412 particleChangeForGamma->SetProposedKineticEnergy(scatteredEnergy);
413
414 // Set the energy deposited by the ionization
415 particleChangeForGamma->ProposeLocalEnergyDeposit(k-scatteredEnergy-secondaryKinetic);
416
417 // Create the new particle with its characteristics
418 G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic) ;
419 fvect->push_back(dp);
420
421 // Check if the auger model is activated (ie instanciated)
423 {
424 // run the PTB Auger model
425 if(materialName!="G4_WATER") fDNAPTBAugerModel->ComputeAugerEffect(fvect, materialName, bindingEnergy);
426 }
427 }
428}
double z() const
double x() const
double y() const
void ComputeAugerEffect(std::vector< G4DynamicParticle * > *fvect, const G4String &materialNameIni, G4double bindingEnergy)
ComputeAugerEffect Main method to be called by the ionisation model.
G4double RandomizeEjectedElectronEnergy(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4int shell, const G4String &materialName)
G4double RandomizeEjectedElectronEnergyFromCumulated(G4ParticleDefinition *particleDefinition, G4double k, G4int shell, const G4String &materialName)
RandomizeEjectedElectronEnergyFromCumulated Uses the cumulated tables to find the energy of the eject...
void RandomizeEjectedElectronDirection(G4ParticleDefinition *aParticleDefinition, G4double incomingParticleEnergy, G4double outgoingParticleEnergy, G4double &cosTheta, G4double &phi)
RandomizeEjectedElectronDirection Method to calculate the ejected electron direction.
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:93
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 G4InuclSpecialFunctions::bindingEnergy(), G4DNAPTBAugerModel::ComputeAugerEffect(), G4Electron::Electron(), source.hepunit::electron_mass_c2, G4Electron::ElectronDefinition(), eV, g4zmq::exit(), fDNAPTBAugerModel, G4cout, G4endl, G4DynamicParticle::GetDefinition(), G4VDNAModel::GetHighELimit(), G4DynamicParticle::GetKineticEnergy(), G4VDNAModel::GetLowELimit(), G4DynamicParticle::GetMomentumDirection(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetPDGMass(), CLHEP::Hep3Vector::getX(), CLHEP::Hep3Vector::getY(), CLHEP::Hep3Vector::getZ(), G4DNAPTBIonisationStructure::IonisationEnergy(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChangeForGamma::ProposeMomentumDirection(), ptbStructure, RandomizeEjectedElectronDirection(), RandomizeEjectedElectronEnergy(), RandomizeEjectedElectronEnergyFromCumulated(), G4VDNAModel::RandomSelectShell(), CLHEP::Hep3Vector::rotateUz(), G4ParticleChangeForGamma::SetProposedKineticEnergy(), CLHEP::Hep3Vector::unit(), verboseLevel, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

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

Field Documentation

◆ diffCrossSectionData

TriDimensionMap G4DNAPTBIonisationModel::diffCrossSectionData
private

Definition at line 131 of file G4DNAPTBIonisationModel.hh.

Referenced by DifferentialCrossSection(), and ReadDiffCSFile().

◆ fDNAPTBAugerModel

G4DNAPTBAugerModel* G4DNAPTBIonisationModel::fDNAPTBAugerModel
private

PTB Auger model instanciated in the constructor and deleted in the destructor of the class.

Definition at line 124 of file G4DNAPTBIonisationModel.hh.

Referenced by G4DNAPTBIonisationModel(), Initialise(), SampleSecondaries(), and ~G4DNAPTBIonisationModel().

◆ fEMapWithVector

VecMap G4DNAPTBIonisationModel::fEMapWithVector
private

Definition at line 135 of file G4DNAPTBIonisationModel.hh.

Referenced by DifferentialCrossSection(), and ReadDiffCSFile().

◆ fEnergySecondaryData

TriDimensionMap G4DNAPTBIonisationModel::fEnergySecondaryData
private

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

◆ fProbaShellMap

VecMapWithShell G4DNAPTBIonisationModel::fProbaShellMap
private

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

◆ fTMapWithVec

std::map<G4String, std::map<G4String, std::vector<double> > > G4DNAPTBIonisationModel::fTMapWithVec
private

◆ ptbStructure

G4DNAPTBIonisationStructure G4DNAPTBIonisationModel::ptbStructure
private

ptbStructure class which contains the shell binding energies

Definition at line 128 of file G4DNAPTBIonisationModel.hh.

Referenced by DifferentialCrossSection(), RandomizeEjectedElectronEnergy(), RandomizeEjectedElectronEnergyFromCumulated(), ReadDiffCSFile(), and SampleSecondaries().

◆ verboseLevel

G4int G4DNAPTBIonisationModel::verboseLevel
private

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