Geant4 10.7.0
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 Member Functions

 G4DNAPTBIonisationModel (const G4String &applyToMaterial="all", const G4ParticleDefinition *p=0, const G4String &nam="DNAPTBIonisationModel", const G4bool isAuger=true)
 G4DNAPTBIonisationModel Constructor.
 
virtual ~G4DNAPTBIonisationModel ()
 ~G4DNAPTBIonisationModel Destructor
 
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.
 
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.
 
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.
 
- Public Member Functions inherited from G4VDNAModel
 G4VDNAModel (const G4String &nam, const G4String &applyToMaterial)
 G4VDNAModel Constructeur of the G4VDNAModel class.
 
virtual ~G4VDNAModel ()
 ~G4VDNAModel
 
virtual void Initialise (const G4ParticleDefinition *particle, const G4DataVector &cuts, G4ParticleChangeForGamma *fpChangeForGamme=nullptr)=0
 Initialise Each model must implement an Initialize method.
 
virtual G4double CrossSectionPerVolume (const G4Material *material, const G4String &materialName, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)=0
 CrossSectionPerVolume Every model must implement its own CrossSectionPerVolume method. It is used by the process to determine the step path and must return a cross section times a number of molecules per volume unit.
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4String &materialName, const G4DynamicParticle *, G4ParticleChangeForGamma *particleChangeForGamma, G4double tmin=0, G4double tmax=DBL_MAX)=0
 SampleSecondaries Each model must implement SampleSecondaries to decide if a particle will be created after the ModelInterface or if any charateristic of the incident particle will change.
 
G4bool IsMaterialDefine (const G4String &materialName)
 IsMaterialDefine Check if the given material is defined in the simulation.
 
G4bool IsMaterialExistingInModel (const G4String &materialName)
 IsMaterialExistingInModel Check if the given material is defined in the current model class.
 
G4bool IsParticleExistingInModelForMaterial (const G4String &particleName, const G4String &materialName)
 IsParticleExistingInModelForMaterial To check two things: 1- is the material existing in model ? 2- if yes, is the particle defined for that material ?
 
G4String GetName ()
 GetName.
 
G4double GetHighELimit (const G4String &material, const G4String &particle)
 GetHighEnergyLimit.
 
G4double GetLowELimit (const G4String &material, const G4String &particle)
 GetLowEnergyLimit.
 
void SetHighELimit (const G4String &material, const G4String &particle, G4double lim)
 SetHighEnergyLimit.
 
void SetLowELimit (const G4String &material, const G4String &particle, G4double lim)
 SetLowEnergyLimit.
 

Additional Inherited Members

- Protected Types inherited from G4VDNAModel
typedef std::map< G4String, std::map< G4String, G4DNACrossSectionDataSet *, std::less< G4String > > > TableMapData
 
typedef std::map< G4String, std::map< G4String, G4double > > RatioMapData
 
typedef std::map< G4String, G4double >::const_iterator ItCompoMapData
 
- Protected Member Functions inherited from G4VDNAModel
TableMapDataGetTableData ()
 GetTableData.
 
std::vector< G4StringBuildApplyToMatVect (const G4String &materials)
 BuildApplyToMatVect Build the material name vector which is used to know the materials the user want to include in the model.
 
void ReadAndSaveCSFile (const G4String &materialName, const G4String &particleName, const G4String &file, G4double scaleFactor)
 ReadAndSaveCSFile Read and save a "simple" cross section file : use of G4DNACrossSectionDataSet->loadData()
 
G4int RandomSelectShell (G4double k, const G4String &particle, const G4String &materialName)
 RandomSelectShell Method to randomely select a shell from the data table uploaded. The size of the table (number of columns) is used to determine the total number of possible shells.
 
void AddCrossSectionData (G4String materialName, G4String particleName, G4String fileCS, G4String fileDiffCS, G4double scaleFactor)
 AddCrossSectionData Method used during the initialization of the model class to add a new material. It adds a material to the model and fills vectors with informations.
 
void AddCrossSectionData (G4String materialName, G4String particleName, G4String fileCS, G4double scaleFactor)
 AddCrossSectionData Method used during the initialization of the model class to add a new material. It adds a material to the model and fills vectors with informations. Not every model needs differential cross sections.
 
void LoadCrossSectionData (const G4String &particleName)
 LoadCrossSectionData Method to loop on all the registered materials in the model and load the corresponding data.
 
virtual void ReadDiffCSFile (const G4String &materialName, const G4String &particleName, const G4String &path, const G4double scaleFactor)
 ReadDiffCSFile Virtual method that need to be implemented if one wish to use the differential cross sections. The read method for that kind of information is not standardized yet.
 
void EnableForMaterialAndParticle (const G4String &materialName, const G4String &particleName)
 EnableMaterialAndParticle.
 

Detailed Description

The G4DNAPTBIonisationModel class Implements the PTB ionisation model.

Definition at line 53 of file G4DNAPTBIonisationModel.hh.

Constructor & Destructor Documentation

◆ G4DNAPTBIonisationModel()

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
64 fDNAPTBAugerModel = 0;
65 }
66}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
The G4DNAPTBAugerModel class Implement the PTB Auger model.
The G4VDNAModel class.
Definition: G4VDNAModel.hh:50

◆ ~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
73 if(fDNAPTBAugerModel) delete fDNAPTBAugerModel;
74}

Member Function Documentation

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

◆ 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
245 if(fDNAPTBAugerModel) fDNAPTBAugerModel->Initialise();
246}
virtual void Initialise()
Initialise Set the verbose value.
static G4Electron * ElectronDefinition()
Definition: G4Electron.cc:88
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:87
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

◆ 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)
422 if(fDNAPTBAugerModel)
423 {
424 // run the PTB Auger model
425 if(materialName!="G4_WATER") fDNAPTBAugerModel->ComputeAugerEffect(fvect, materialName, bindingEnergy);
426 }
427 }
428}
int G4int
Definition: G4Types.hh:85
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 IonisationEnergy(G4int level, const G4String &materialName)
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)
G4double bindingEnergy(G4int A, G4int Z)

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