Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAPTBIonisationModel Class Reference

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

#include <G4DNAPTBIonisationModel.hh>

+ Inheritance diagram for G4DNAPTBIonisationModel:

Public Types

using TriDimensionMap
 
using VecMap
 
using VecMapWithShell
 

Public Member Functions

 G4DNAPTBIonisationModel (const G4String &applyToMaterial="all", const G4ParticleDefinition *p=nullptr, const G4String &nam="DNAPTBIonisationModel", const G4bool isAuger=true)
 G4DNAPTBIonisationModel Constructor.
 
 ~G4DNAPTBIonisationModel () override=default
 ~G4DNAPTBIonisationModel Destructor
 
 G4DNAPTBIonisationModel (const G4DNAPTBIonisationModel &)=delete
 
G4DNAPTBIonisationModeloperator= (const G4DNAPTBIonisationModel &right)=delete
 
void Initialise (const G4ParticleDefinition *particle, const G4DataVector &data) override
 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.
 
G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
 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.
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double tmax) override
 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.
 
- Public Member Functions inherited from G4VDNAModel
 G4VDNAModel (const G4String &nam, const G4String &applyToMaterial="")
 G4VDNAModel Constructeur of the G4VDNAModel class.
 
 ~G4VDNAModel () override
 ~G4VDNAModel
 
G4bool IsMaterialDefine (const size_t &materialID)
 IsMaterialDefine Check if the given material is defined in the simulation.
 
G4bool IsParticleExistingInModelForMaterial (const G4ParticleDefinition *particleName, const size_t &materialID)
 IsParticleExistingInModelForMaterial To check two things: 1- is the material existing in model ? 2- if yes, is the particle defined for that material ?
 
G4String GetName ()
 GetName.
 
G4double GetHighELimit (const size_t &materialID, const G4ParticleDefinition *particle)
 GetHighEnergyLimit.
 
G4double GetLowELimit (const size_t &materialID, const G4ParticleDefinition *particle)
 GetLowEnergyLimit.
 
void SetHighELimit (const size_t &materialID, const G4ParticleDefinition *particle, G4double lim)
 SetHighEnergyLimit.
 
void SetLowELimit (const size_t &materialID, const G4ParticleDefinition *particle, G4double lim)
 SetLowEnergyLimit.
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double &length, G4double &eloss)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
void SetLPMFlag (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Public Attributes

G4ParticleChangeForGammafParticleChangeForGamma = nullptr
 

Additional Inherited Members

- Protected Types inherited from G4VDNAModel
using MaterialParticleMapData
 
- Protected Member Functions inherited from G4VDNAModel
MaterialParticleMapDataGetData ()
 GetTableData.
 
std::vector< G4StringBuildApplyToMatVect (const G4String &materials)
 BuildApplyToMatVect Build the material name vector which is used to know the materials the user want to include in the model.
 
void ReadAndSaveCSFile (const size_t &materialID, const G4ParticleDefinition *p, const G4String &file, const G4double &scaleFactor)
 ReadAndSaveCSFile Read and save a "simple" cross section file : use of G4DNACrossSectionDataSet->loadData()
 
G4int RandomSelectShell (const G4double &k, const G4ParticleDefinition *particle, const size_t &materialName)
 RandomSelectShell Method to randomely select a shell from the data table uploaded. The size of the table (number of columns) is used to determine the total number of possible shells.
 
void AddCrossSectionData (const size_t &materialName, const G4ParticleDefinition *particleName, const G4String &fileCS, const G4String &fileDiffCS, const G4double &scaleFactor)
 AddCrossSectionData Method used during the initialization of the model class to add a new material. It adds a material to the model and fills vectors with informations.
 
void AddCrossSectionData (const size_t &materialName, const G4ParticleDefinition *particleName, const G4String &fileCS, const G4double &scaleFactor)
 AddCrossSectionData Method used during the initialization of the model class to add a new material. It adds a material to the model and fills vectors with informations. Not every model needs differential cross sections.
 
void LoadCrossSectionData (const G4ParticleDefinition *particleName)
 LoadCrossSectionData Method to loop on all the registered materials in the model and load the corresponding data.
 
virtual void ReadDiffCSFile (const size_t &materialName, const G4ParticleDefinition *particleName, const G4String &path, const G4double &scaleFactor)
 ReadDiffCSFile Virtual method that need to be implemented if one wish to use the differential cross sections. The read method for that kind of information is not standardized yet.
 
void EnableForMaterialAndParticle (const size_t &materialID, const G4ParticleDefinition *p)
 EnableMaterialAndParticle.
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 
- Protected Attributes inherited from G4VDNAModel
G4bool isInitialised = false
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 
const G4MaterialpBaseMaterial = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4double inveplus
 
G4double pFactor = 1.0
 
std::size_t currentCoupleIndex = 0
 
std::size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Detailed Description

The G4DNAPTBIonisationModel class Implements the PTB ionisation model.

Definition at line 50 of file G4DNAPTBIonisationModel.hh.

Member Typedef Documentation

◆ TriDimensionMap

Initial value:
std::map<std::size_t, std::map<const G4ParticleDefinition*,
std::map<G4double, std::map<G4double, std::map<G4double, G4double>>>>>

Definition at line 53 of file G4DNAPTBIonisationModel.hh.

◆ VecMap

Initial value:
std::map<std::size_t,
std::map<const G4ParticleDefinition*, std::map<G4double, std::vector<G4double>>>>

Definition at line 56 of file G4DNAPTBIonisationModel.hh.

◆ VecMapWithShell

Initial value:
std::map<std::size_t, std::map<const G4ParticleDefinition*,
std::map<G4double, std::map<G4double, std::vector<G4double>>>>>

Definition at line 58 of file G4DNAPTBIonisationModel.hh.

Constructor & Destructor Documentation

◆ G4DNAPTBIonisationModel() [1/2]

G4DNAPTBIonisationModel::G4DNAPTBIonisationModel ( const G4String & applyToMaterial = "all",
const G4ParticleDefinition * p = nullptr,
const G4String & nam = "DNAPTBIonisationModel",
const G4bool isAuger = true )
explicit

G4DNAPTBIonisationModel Constructor.

Parameters
applyToMaterial
p
nam
isAuger

Definition at line 39 of file G4DNAPTBIonisationModel.cc.

42 : G4VDNAModel(nam, applyToMaterial)
43{
44 if (isAuger) {
45 // create the PTB Auger model
46 fpDNAPTBAugerModel = std::make_unique<G4DNAPTBAugerModel>("e-_G4DNAPTBAugerModel");
47 }
48 fpTHF = G4Material::GetMaterial("THF", false);
49 fpPY = G4Material::GetMaterial("PY", false);
50 fpPU = G4Material::GetMaterial("PU", false);
51 fpTMP = G4Material::GetMaterial("TMP", false);
52 fpG4_WATER = G4Material::GetMaterial("G4_WATER", false);
53 fpBackbone_THF = G4Material::GetMaterial("backbone_THF", false);
54 fpCytosine_PY = G4Material::GetMaterial("cytosine_PY", false);
55 fpThymine_PY = G4Material::GetMaterial("thymine_PY", false);
56 fpAdenine_PU = G4Material::GetMaterial("adenine_PU", false);
57 fpBackbone_TMP = G4Material::GetMaterial("backbone_TMP", false);
58 fpGuanine_PU = G4Material::GetMaterial("guanine_PU", false);
59 fpN2 = G4Material::GetMaterial("N2", false);
60}
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
G4VDNAModel(const G4String &nam, const G4String &applyToMaterial="")
G4VDNAModel Constructeur of the G4VDNAModel class.

◆ ~G4DNAPTBIonisationModel()

G4DNAPTBIonisationModel::~G4DNAPTBIonisationModel ( )
overridedefault

~G4DNAPTBIonisationModel Destructor

◆ G4DNAPTBIonisationModel() [2/2]

G4DNAPTBIonisationModel::G4DNAPTBIonisationModel ( const G4DNAPTBIonisationModel & )
delete

Member Function Documentation

◆ CrossSectionPerVolume()

G4double G4DNAPTBIonisationModel::CrossSectionPerVolume ( const G4Material * material,
const G4ParticleDefinition * p,
G4double ekin,
G4double emin,
G4double emax )
overridevirtual

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 255 of file G4DNAPTBIonisationModel.cc.

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 const std::size_t& matID = material->GetIndex();
266
267 // Set the low and high energy limits
268 G4double lowLim = fpModelData->GetLowELimit(matID, p);
269 G4double highLim = fpModelData->GetHighELimit(matID, p);
270
271 // Check that we are in the correct energy range
272 if (ekin >= lowLim && ekin < highLim) {
273 // Get the map with all the model data tables
274 auto tableData = fpModelData->GetData();
275 if ((*tableData)[matID][p] == nullptr) {
276 G4Exception("G4DNAPTBIonisationModel::CrossSectionPerVolume", "em00236", FatalException,
277 "No model is registered");
278 }
279 // Retrieve the cross section value for the current material, particle and energy values
280 sigma = (*tableData)[matID][p]->FindValue(ekin);
281
282 if (verboseLevel > 2) {
283 G4cout << "__________________________________" << G4endl;
284 G4cout << "°°° G4DNAPTBIonisationModel - XS INFO START" << G4endl;
285 G4cout << "°°° Kinetic energy(eV)=" << ekin / eV << " particle : " << particleName << G4endl;
286 G4cout << "°°° Cross section per " << matID << " index molecule (cm^2)=" << sigma / cm / cm
287 << G4endl;
288 G4cout << "°°° G4DNAPTBIonisationModel - XS INFO END" << G4endl;
289 }
290 }
291
292 // Return the cross section value
293 auto MolDensity =
295 return sigma * MolDensity;
296}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
double G4double
Definition G4Types.hh:83
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
const std::vector< G4double > * GetNumMolPerVolTableFor(const G4Material *) const
Retrieve a table of molecular densities (number of molecules per unit volume) in the G4 unit system f...
static G4DNAMolecularMaterial * Instance()
std::size_t GetIndex() const
const G4String & GetParticleName() const
G4double GetHighELimit(const size_t &materialID, const G4ParticleDefinition *particle)
GetHighEnergyLimit.
G4double GetLowELimit(const size_t &materialID, const G4ParticleDefinition *particle)
GetLowEnergyLimit.
MaterialParticleMapData * GetData()
GetTableData.

◆ Initialise()

void G4DNAPTBIonisationModel::Initialise ( const G4ParticleDefinition * particle,
const G4DataVector & data )
overridevirtual

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 64 of file G4DNAPTBIonisationModel.cc.

66{
67 if (isInitialised) {
68 return;
69 }
70 if (verboseLevel > 3) {
71 G4cout << "Calling G4DNAPTBIonisationModel::Initialise()" << G4endl;
72 }
73
74 G4double scaleFactor = 1e-16 * cm * cm;
75 G4double scaleFactorBorn = (1.e-22 / 3.343) * m * m;
76
79
80 //*******************************************************
81 // Cross section data
82 //*******************************************************
83 std::size_t index;
84 if (particle == electronDef) {
85 // Raw materials
86 //
87 // MPietrzak
88 if (fpN2 != nullptr) {
89 index = fpN2->GetIndex();
90 AddCrossSectionData(index, particle, "dna/sigma_ionisation_e-_PTB_N2",
91 "dna/sigmadiff_cumulated_ionisation_e-_PTB_N2", scaleFactor);
92 SetLowELimit(index, particle, 15.5 * eV);
93 SetHighELimit(index, particle, 1.02 * MeV);
94 }
95
96 // MPietrzak
97
98 if (fpTHF != nullptr) {
99 index = fpTHF->GetIndex();
100 AddCrossSectionData(index, particle, "dna/sigma_ionisation_e-_PTB_THF",
101 "dna/sigmadiff_cumulated_ionisation_e-_PTB_THF", scaleFactor);
102 SetLowELimit(index, particle, 12. * eV);
103 SetHighELimit(index, particle, 1. * keV);
104 }
105
106 if (fpPY != nullptr) {
107 index = fpPY->GetIndex();
108 AddCrossSectionData(index, particle, "dna/sigma_ionisation_e-_PTB_PY",
109 "dna/sigmadiff_cumulated_ionisation_e-_PTB_PY", scaleFactor);
110 SetLowELimit(index, particle, 12. * eV);
111 SetHighELimit(index, particle, 1. * keV);
112 }
113
114 if (fpPU != nullptr) {
115 index = fpPU->GetIndex();
116 AddCrossSectionData(index, particle, "dna/sigma_ionisation_e-_PTB_PU",
117 "dna/sigmadiff_cumulated_ionisation_e-_PTB_PU", scaleFactor);
118 SetLowELimit(index, particle, 12. * eV);
119 SetHighELimit(index, particle, 1. * keV);
120 }
121 if (fpTMP != nullptr) {
122 index = fpTMP->GetIndex();
123 AddCrossSectionData(index, particle, "dna/sigma_ionisation_e-_PTB_TMP",
124 "dna/sigmadiff_cumulated_ionisation_e-_PTB_TMP", scaleFactor);
125 SetLowELimit(index, particle, 12. * eV);
126 SetHighELimit(index, particle, 1. * keV);
127 }
128
129 if (fpG4_WATER != nullptr) {
130 index = fpG4_WATER->GetIndex();
131 AddCrossSectionData(index, particle, "dna/sigma_ionisation_e_born",
132 "dna/sigmadiff_ionisation_e_born", scaleFactorBorn);
133 SetLowELimit(index, particle, 12. * eV);
134 SetHighELimit(index, particle, 1. * keV);
135 }
136 // DNA materials
137 //
138 if (fpBackbone_THF != nullptr) {
139 index = fpBackbone_THF->GetIndex();
140 AddCrossSectionData(index, particle, "dna/sigma_ionisation_e-_PTB_THF",
141 "dna/sigmadiff_cumulated_ionisation_e-_PTB_THF", scaleFactor * 33. / 30);
142 SetLowELimit(index, particle, 12. * eV);
143 SetHighELimit(index, particle, 1. * keV);
144 }
145
146 if (fpCytosine_PY != nullptr) {
147 index = fpCytosine_PY->GetIndex();
148 AddCrossSectionData(index, particle, "dna/sigma_ionisation_e-_PTB_PY",
149 "dna/sigmadiff_cumulated_ionisation_e-_PTB_PY", scaleFactor * 42. / 30);
150 SetLowELimit(index, particle, 12. * eV);
151 SetHighELimit(index, particle, 1. * keV);
152 }
153
154 if (fpThymine_PY != nullptr) {
155 index = fpThymine_PY->GetIndex();
156 AddCrossSectionData(index, particle, "dna/sigma_ionisation_e-_PTB_PY",
157 "dna/sigmadiff_cumulated_ionisation_e-_PTB_PY", scaleFactor * 48. / 30);
158 SetLowELimit(index, particle, 12. * eV);
159 SetHighELimit(index, particle, 1. * keV);
160 }
161
162 if (fpAdenine_PU != nullptr) {
163 index = fpAdenine_PU->GetIndex();
164 AddCrossSectionData(index, particle, "dna/sigma_ionisation_e-_PTB_PU",
165 "dna/sigmadiff_cumulated_ionisation_e-_PTB_PU", scaleFactor * 50. / 44);
166 SetLowELimit(index, particle, 12. * eV);
167 SetHighELimit(index, particle, 1. * keV);
168 }
169
170 if (fpGuanine_PU != nullptr) {
171 index = fpGuanine_PU->GetIndex();
172 AddCrossSectionData(index, particle, "dna/sigma_ionisation_e-_PTB_PU",
173 "dna/sigmadiff_cumulated_ionisation_e-_PTB_PU", scaleFactor * 56. / 44);
174 SetLowELimit(index, particle, 12. * eV);
175 SetHighELimit(index, particle, 1. * keV);
176 }
177
178 if (fpBackbone_TMP != nullptr) {
179 index = fpBackbone_TMP->GetIndex();
180 AddCrossSectionData(index, particle, "dna/sigma_ionisation_e-_PTB_TMP",
181 "dna/sigmadiff_cumulated_ionisation_e-_PTB_TMP", scaleFactor * 33. / 50);
182 SetLowELimit(index, particle, 12. * eV);
183 SetHighELimit(index, particle, 1. * keV);
184 }
185 }
186
187 else if (particle == protonDef) {
188 G4String particleName = particle->GetParticleName();
189
190 // Raw materials
191 //
192 if (fpTHF != nullptr) {
193 index = fpTHF->GetIndex();
194 AddCrossSectionData(index, particle, "dna/sigma_ionisation_p_HKS_THF",
195 "dna/sigmadiff_cumulated_ionisation_p_PTB_THF", scaleFactor);
196 SetLowELimit(index, particle, 70. * keV);
197 SetHighELimit(index, particle, 10. * MeV);
198 }
199
200 if (fpPY != nullptr) {
201 index = fpPY->GetIndex();
202 AddCrossSectionData(index, particle, "dna/sigma_ionisation_p_HKS_PY",
203 "dna/sigmadiff_cumulated_ionisation_p_PTB_PY", scaleFactor);
204 SetLowELimit(index, particle, 70. * keV);
205 SetHighELimit(index, particle, 10. * MeV);
206 }
207 /*
208 AddCrossSectionData("PU",
209 particleName,
210 "dna/sigma_ionisation_e-_PTB_PU",
211 "dna/sigmadiff_cumulated_ionisation_e-_PTB_PU",
212 scaleFactor);
213 SetLowELimit("PU", particleName2, 70.*keV);
214 SetHighELimit("PU", particleName2, 10.*keV);
215*/
216
217 if (fpTMP != nullptr) {
218 index = fpTMP->GetIndex();
219 AddCrossSectionData(index, particle, "dna/sigma_ionisation_p_HKS_TMP",
220 "dna/sigmadiff_cumulated_ionisation_p_PTB_TMP", scaleFactor);
221 SetLowELimit(index, particle, 70. * keV);
222 SetHighELimit(index, particle, 10. * MeV);
223 }
224 }
225 // *******************************************************
226 // deal with composite materials
227 // *******************************************************
229 LoadCrossSectionData(particle);
231 fpModelData = this;
232 }
233 else {
234 auto dataModel = dynamic_cast<G4DNAPTBIonisationModel*>(
236 if (dataModel == nullptr) {
237 G4cout << "G4DNAPTBIonisationModel::Initialise:: not good modelData" << G4endl;
238 G4Exception("G4DNAPTBIonisationModel::Initialise", "PTB0004", FatalException,
239 "not good modelData");
240 }
241 else {
242 fpModelData = dataModel;
243 }
244 }
245 // initialise DNAPTBAugerModel
246 if (fpDNAPTBAugerModel) {
247 fpDNAPTBAugerModel->Initialise();
248 }
250 isInitialised = true;
251}
static G4DNAMaterialManager * Instance()
void SetMasterDataModel(const DNAModelType &t, G4VEmModel *m)
G4VEmModel * GetModel(const DNAModelType &t)
The G4DNAPTBIonisationModel class Implements the PTB ionisation model.
G4ParticleChangeForGamma * fParticleChangeForGamma
static G4Electron * ElectronDefinition()
Definition G4Electron.cc:86
static G4Proton * ProtonDefinition()
Definition G4Proton.cc:85
void LoadCrossSectionData(const G4ParticleDefinition *particleName)
LoadCrossSectionData Method to loop on all the registered materials in the model and load the corresp...
void SetLowELimit(const size_t &materialID, const G4ParticleDefinition *particle, G4double lim)
SetLowEnergyLimit.
void SetHighELimit(const size_t &materialID, const G4ParticleDefinition *particle, G4double lim)
SetHighEnergyLimit.
void AddCrossSectionData(const size_t &materialName, const G4ParticleDefinition *particleName, const G4String &fileCS, const G4String &fileDiffCS, const G4double &scaleFactor)
AddCrossSectionData Method used during the initialization of the model class to add a new material....
G4bool isInitialised
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4bool IsLocked() const

◆ operator=()

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

◆ SampleSecondaries()

void G4DNAPTBIonisationModel::SampleSecondaries ( std::vector< G4DynamicParticle * > * fvect,
const G4MaterialCutsCouple * pCouple,
const G4DynamicParticle * aDynamicParticle,
G4double tmin,
G4double tmax )
overridevirtual

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

param particleChangeForGamma

Parameters
tmin

param tmax

Implements G4VDNAModel.

Definition at line 300 of file G4DNAPTBIonisationModel.cc.

304{
305 // Get the current particle energy
306 G4double k = aDynamicParticle->GetKineticEnergy();
307 const std::size_t& materialID = pCouple->GetMaterial()->GetIndex();
308
309 // Get the current particle name
310 const auto& p = aDynamicParticle->GetDefinition();
311 auto materialName = pCouple->GetMaterial()->GetName();
312 // Get the energy limits
313 G4double lowLim = fpModelData->GetLowELimit(materialID, p);
314 G4double highLim = fpModelData->GetHighELimit(materialID, p);
315
316 // Check if we are in the correct energy range
317 if (k >= lowLim && k < highLim) {
318 G4ParticleMomentum primaryDirection = aDynamicParticle->GetMomentumDirection();
319 G4double particleMass = aDynamicParticle->GetDefinition()->GetPDGMass();
320 G4double totalEnergy = k + particleMass;
321 G4double pSquare = k * (totalEnergy + particleMass);
322 G4double totalMomentum = std::sqrt(pSquare);
323
324 // Get the ionisation shell from a random sampling
325 G4int ionizationShell = fpModelData->RandomSelectShell(k, p, materialID);
326
327 // Get the binding energy from the ptbStructure class
328 G4double bindingEnergy = ptbStructure.IonisationEnergy(ionizationShell, materialID);
329
330 // Initialize the secondary kinetic energy to a negative value.
331 G4double secondaryKinetic(-1000 * eV);
332
333 if (fpG4_WATER == nullptr || materialID != fpG4_WATER->GetIndex()) {
334 // Get the energy of the secondary particle
335 secondaryKinetic = fpModelData->RandomizeEjectedElectronEnergyFromCumulated(
336 aDynamicParticle->GetDefinition(), k / eV, ionizationShell, materialID);
337 }
338 else {
339 secondaryKinetic = fpModelData->RandomizeEjectedElectronEnergy(
340 aDynamicParticle->GetDefinition(), k, ionizationShell, materialID);
341 }
342
343 if (secondaryKinetic <= 0) {
344 G4cout << "Fatal error *************************************** " << secondaryKinetic / eV
345 << G4endl;
346 G4cout << "secondaryKinetic: " << secondaryKinetic / eV << G4endl;
347 G4cout << "k: " << k / eV << G4endl;
348 G4cout << "shell: " << ionizationShell << G4endl;
349 G4cout << "material:" << materialName << G4endl;
350 G4Exception("G4DNAPTBIonisationModel::SampleSecondaries", "em0026", FatalException,
351 "Fatal error:: scatteredEnergy <= 0");
352 }
353
354 G4double cosTheta = 0.;
355 G4double phi = 0.;
356 RandomizeEjectedElectronDirection(aDynamicParticle->GetDefinition(), k, secondaryKinetic,
357 cosTheta, phi);
358
359 G4double sinTheta = std::sqrt(1. - cosTheta * cosTheta);
360 G4double dirX = sinTheta * std::cos(phi);
361 G4double dirY = sinTheta * std::sin(phi);
362 G4double dirZ = cosTheta;
363 G4ThreeVector deltaDirection(dirX, dirY, dirZ);
364 deltaDirection.rotateUz(primaryDirection);
365
366 // The model is written only for electron and thus we want the change the direction of the
367 // incident electron after each ionization. However, if other particle are going to be
368 // introduced within this model the following should be added:
369 //
370 // Check if the particle is an electron
371 if (aDynamicParticle->GetDefinition() == G4Electron::ElectronDefinition()) {
372 // If yes do the following code until next commented "else" statement
373
374 G4double deltaTotalMomentum =
375 std::sqrt(secondaryKinetic * (secondaryKinetic + 2. * electron_mass_c2));
376 G4double finalPx =
377 totalMomentum * primaryDirection.x() - deltaTotalMomentum * deltaDirection.x();
378 G4double finalPy =
379 totalMomentum * primaryDirection.y() - deltaTotalMomentum * deltaDirection.y();
380 G4double finalPz =
381 totalMomentum * primaryDirection.z() - deltaTotalMomentum * deltaDirection.z();
382 G4double finalMomentum = std::sqrt(finalPx * finalPx + finalPy * finalPy + finalPz * finalPz);
383 finalPx /= finalMomentum;
384 finalPy /= finalMomentum;
385 finalPz /= finalMomentum;
386
387 G4ThreeVector direction(finalPx, finalPy, finalPz);
388 if (direction.unit().getX() > 1 || direction.unit().getY() > 1 || direction.unit().getZ() > 1)
389 {
390 G4cout << "Fatal error ****************************" << G4endl;
391 G4cout << "direction problem " << direction.unit() << G4endl;
392 G4Exception("G4DNAPTBIonisationModel::SampleSecondaries", "em0017", FatalException,
393 "Fatal error:: direction problem");
394 }
395
396 // Give the new direction to the particle
398 }
399 // If the particle is not an electron
400 else {
402 }
403
404 // note that secondaryKinetic is the energy of the delta ray, not of all secondaries.
405 G4double scatteredEnergy = k - bindingEnergy - secondaryKinetic;
406
407 if (scatteredEnergy <= 0) {
408 G4cout << "Fatal error ****************************" << G4endl;
409 G4cout << "k: " << k / eV << G4endl;
410 G4cout << "secondaryKinetic: " << secondaryKinetic / eV << G4endl;
411 G4cout << "shell: " << ionizationShell << G4endl;
412 G4cout << "bindingEnergy: " << bindingEnergy / eV << G4endl;
413 G4cout << "scatteredEnergy: " << scatteredEnergy / eV << G4endl;
414 G4cout << "material: " << materialName << G4endl;
415 G4Exception("G4DNAPTBIonisationModel::SampleSecondaries", "em0016", FatalException,
416 "Fatal error:: scatteredEnergy <= 0");
417 }
418
419 // Set the new energy of the particle
421
422 // Set the energy deposited by the ionization
423 fParticleChangeForGamma->ProposeLocalEnergyDeposit(k - scatteredEnergy - secondaryKinetic);
424
425 // Create the new particle with its characteristics
426 auto dp = new G4DynamicParticle(G4Electron::Electron(), deltaDirection, secondaryKinetic);
427 fvect->push_back(dp);
428
429 // Check if the auger model is activated (ie instanciated)
430 if (fpDNAPTBAugerModel) {
431 // run the PTB Auger model
432 if (materialName != "G4_WATER") {
433 fpDNAPTBAugerModel->ComputeAugerEffect(fvect, materialName, bindingEnergy);
434 }
435 }
436 }
437}
int G4int
Definition G4Types.hh:85
double z() const
double x() const
double y() const
G4double IonisationEnergy(G4int level, const size_t &materialName)
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition G4Electron.cc:91
const G4Material * GetMaterial() const
const G4String & GetName() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
G4int RandomSelectShell(const G4double &k, const G4ParticleDefinition *particle, const size_t &materialName)
RandomSelectShell Method to randomely select a shell from the data table uploaded....
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double bindingEnergy(G4int A, G4int Z)

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNAPTBIonisationModel::fParticleChangeForGamma = nullptr

Definition at line 118 of file G4DNAPTBIonisationModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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