Geant4 11.1.1
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 260 of file G4DNAPTBIonisationModel.cc.

266{
267 if(verboseLevel > 3)
268 G4cout << "Calling CrossSectionPerVolume() of G4DNAPTBIonisationModel" << G4endl;
269
270 // initialise the cross section value (output value)
271 G4double sigma(0);
272
273 // Get the current particle name
274 const G4String& particleName = p->GetParticleName();
275
276 // Set the low and high energy limits
277 G4double lowLim = GetLowELimit(materialName, particleName);
278 G4double highLim = GetHighELimit(materialName, particleName);
279
280 // Check that we are in the correct energy range
281 if (ekin >= lowLim && ekin < highLim)
282 {
283 // Get the map with all the model data tables
284 TableMapData* tableData = GetTableData();
285
286 // Retrieve the cross section value for the current material, particle and energy values
287 sigma = (*tableData)[materialName][particleName]->FindValue(ekin);
288
289 if (verboseLevel > 2)
290 {
291 G4cout << "__________________________________" << G4endl;
292 G4cout << "°°° G4DNAPTBIonisationModel - XS INFO START" << G4endl;
293 G4cout << "°°° Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
294 G4cout << "°°° Cross section per "<< materialName <<" molecule (cm^2)=" << sigma/cm/cm << G4endl;
295 G4cout << "°°° G4DNAPTBIonisationModel - XS INFO END" << G4endl;
296 }
297 }
298
299 // Return the cross section value
300 return sigma;
301}
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 //
101 // MPietrzak
103 particleName,
104 "dna/sigma_ionisation_e-_PTB_N2",
105 "dna/sigmadiff_cumulated_ionisation_e-_PTB_N2",
106 scaleFactor);
107 SetLowELimit("N2", particleName, 15.5*eV);
108 SetHighELimit("N2", particleName, 1.02*MeV);
109 // MPietrzak
110
112 particleName,
113 "dna/sigma_ionisation_e-_PTB_THF",
114 "dna/sigmadiff_cumulated_ionisation_e-_PTB_THF",
115 scaleFactor);
116 SetLowELimit("THF", particleName, 12.*eV);
117 SetHighELimit("THF", particleName, 1.*keV);
118
120 particleName,
121 "dna/sigma_ionisation_e-_PTB_PY",
122 "dna/sigmadiff_cumulated_ionisation_e-_PTB_PY",
123 scaleFactor);
124 SetLowELimit("PY", particleName, 12.*eV);
125 SetHighELimit("PY", particleName, 1.*keV);
126
128 particleName,
129 "dna/sigma_ionisation_e-_PTB_PU",
130 "dna/sigmadiff_cumulated_ionisation_e-_PTB_PU",
131 scaleFactor);
132 SetLowELimit("PU", particleName, 12.*eV);
133 SetHighELimit("PU", particleName, 1.*keV);
134
136 particleName,
137 "dna/sigma_ionisation_e-_PTB_TMP",
138 "dna/sigmadiff_cumulated_ionisation_e-_PTB_TMP",
139 scaleFactor);
140 SetLowELimit("TMP", particleName, 12.*eV);
141 SetHighELimit("TMP", particleName, 1.*keV);
142
143 AddCrossSectionData("G4_WATER",
144 particleName,
145 "dna/sigma_ionisation_e_born",
146 "dna/sigmadiff_ionisation_e_born",
147 scaleFactorBorn);
148 SetLowELimit("G4_WATER", particleName, 12.*eV);
149 SetHighELimit("G4_WATER", particleName, 1.*keV);
150
151 // DNA materials
152 //
153 AddCrossSectionData("backbone_THF",
154 particleName,
155 "dna/sigma_ionisation_e-_PTB_THF",
156 "dna/sigmadiff_cumulated_ionisation_e-_PTB_THF",
157 scaleFactor*33./30);
158 SetLowELimit("backbone_THF", particleName, 12.*eV);
159 SetHighELimit("backbone_THF", particleName, 1.*keV);
160
161 AddCrossSectionData("cytosine_PY",
162 particleName,
163 "dna/sigma_ionisation_e-_PTB_PY",
164 "dna/sigmadiff_cumulated_ionisation_e-_PTB_PY",
165 scaleFactor*42./30);
166 SetLowELimit("cytosine_PY", particleName, 12.*eV);
167 SetHighELimit("cytosine_PY", particleName, 1.*keV);
168
169 AddCrossSectionData("thymine_PY",
170 particleName,
171 "dna/sigma_ionisation_e-_PTB_PY",
172 "dna/sigmadiff_cumulated_ionisation_e-_PTB_PY",
173 scaleFactor*48./30);
174 SetLowELimit("thymine_PY", particleName, 12.*eV);
175 SetHighELimit("thymine_PY", particleName, 1.*keV);
176
177 AddCrossSectionData("adenine_PU",
178 particleName,
179 "dna/sigma_ionisation_e-_PTB_PU",
180 "dna/sigmadiff_cumulated_ionisation_e-_PTB_PU",
181 scaleFactor*50./44);
182 SetLowELimit("adenine_PU", particleName, 12.*eV);
183 SetHighELimit("adenine_PU", particleName, 1.*keV);
184
185 AddCrossSectionData("guanine_PU",
186 particleName,
187 "dna/sigma_ionisation_e-_PTB_PU",
188 "dna/sigmadiff_cumulated_ionisation_e-_PTB_PU",
189 scaleFactor*56./44);
190 SetLowELimit("guanine_PU", particleName, 12.*eV);
191 SetHighELimit("guanine_PU", particleName, 1.*keV);
192
193 AddCrossSectionData("backbone_TMP",
194 particleName,
195 "dna/sigma_ionisation_e-_PTB_TMP",
196 "dna/sigmadiff_cumulated_ionisation_e-_PTB_TMP",
197 scaleFactor*33./50);
198 SetLowELimit("backbone_TMP", particleName, 12.*eV);
199 SetHighELimit("backbone_TMP", particleName, 1.*keV);
200 }
201
202 else if (particle == protonDef)
203 {
204 G4String particleName = particle->GetParticleName();
205
206 // Raw materials
207 //
209 particleName,
210 "dna/sigma_ionisation_p_HKS_THF",
211 "dna/sigmadiff_cumulated_ionisation_p_PTB_THF",
212 scaleFactor);
213 SetLowELimit("THF", particleName, 70.*keV);
214 SetHighELimit("THF", particleName, 10.*MeV);
215
216
218 particleName,
219 "dna/sigma_ionisation_p_HKS_PY",
220 "dna/sigmadiff_cumulated_ionisation_p_PTB_PY",
221 scaleFactor);
222 SetLowELimit("PY", particleName, 70.*keV);
223 SetHighELimit("PY", particleName, 10.*MeV);
224
225 /*
226 AddCrossSectionData("PU",
227 particleName,
228 "dna/sigma_ionisation_e-_PTB_PU",
229 "dna/sigmadiff_cumulated_ionisation_e-_PTB_PU",
230 scaleFactor);
231 SetLowELimit("PU", particleName2, 70.*keV);
232 SetHighELimit("PU", particleName2, 10.*keV);
233*/
234
236 particleName,
237 "dna/sigma_ionisation_p_HKS_TMP",
238 "dna/sigmadiff_cumulated_ionisation_p_PTB_TMP",
239 scaleFactor);
240 SetLowELimit("TMP", particleName, 70.*keV);
241 SetHighELimit("TMP", particleName, 10.*MeV);
242 }
243
244 // *******************************************************
245 // deal with composite materials
246 // *******************************************************
247
249
250 // *******************************************************
251 // Verbose
252 // *******************************************************
253
254 // initialise DNAPTBAugerModel
255 if(fDNAPTBAugerModel) fDNAPTBAugerModel->Initialise();
256}
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 305 of file G4DNAPTBIonisationModel.cc.

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