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

#include <G4DNARuddIonisationModel.hh>

Inheritance diagram for G4DNARuddIonisationModel:

Public Member Functions

 G4DNARuddIonisationModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="DNARuddIonisationModel")
 ~G4DNARuddIonisationModel () override
G4DNARuddIonisationModeloperator= (const G4DNARuddIonisationModel &right)=delete
 G4DNARuddIonisationModel (const G4DNARuddIonisationModel &)=delete
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
void SelectStationary (G4bool input)
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 G4Material *, const G4ParticleDefinition *, const G4double kinEnergy, const G4double cutEnergy, 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)
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)
void SetMasterThread (G4bool)
G4VEmModeloperator= (const G4VEmModel &right)=delete
 G4VEmModel (const G4VEmModel &)=delete

Protected Attributes

G4ParticleChangeForGammafParticleChangeForGamma = nullptr
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

Additional Inherited Members

Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
G4ParticleChangeForGammaGetParticleChangeForGamma ()
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
const G4MaterialCutsCoupleCurrentCouple () const
void SetCurrentElement (const G4Element *)

Detailed Description

Definition at line 45 of file G4DNARuddIonisationModel.hh.

Constructor & Destructor Documentation

◆ G4DNARuddIonisationModel() [1/2]

G4DNARuddIonisationModel::G4DNARuddIonisationModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "DNARuddIonisationModel" )

Definition at line 49 of file G4DNARuddIonisationModel.cc.

50 :
51G4VEmModel(nam)
52{
53 slaterEffectiveCharge[0] = 0.;
54 slaterEffectiveCharge[1] = 0.;
55 slaterEffectiveCharge[2] = 0.;
56 sCoefficient[0] = 0.;
57 sCoefficient[1] = 0.;
58 sCoefficient[2] = 0.;
59
60 lowEnergyLimitForZ1 = 0 * eV;
61 lowEnergyLimitForZ2 = 0 * eV;
62 lowEnergyLimitOfModelForZ1 = 100 * eV;
63 lowEnergyLimitOfModelForZ2 = 1 * keV;
64 killBelowEnergyForZ1 = lowEnergyLimitOfModelForZ1;
65 killBelowEnergyForZ2 = lowEnergyLimitOfModelForZ2;
66
67 verboseLevel = 0;
68 // Verbosity scale:
69 // 0 = nothing
70 // 1 = warning for energy non-conservation
71 // 2 = details of energy budget
72 // 3 = calculation of cross sections, file openings, sampling of atoms
73 // 4 = entering in methods
74
75 if (verboseLevel > 0)
76 {
77 G4cout << "Rudd ionisation model is constructed " << G4endl;
78 }
79
80 // Define default angular generator
81 SetAngularDistribution(new G4DNARuddAngle());
82
83 // Mark this model as "applicable" for atomic deexcitation
85
86 // Selection of stationary mode
87
88 statCode = false;
89}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void SetDeexcitationFlag(G4bool val)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67
void SetAngularDistribution(G4VEmAngularDistribution *)

Referenced by G4DNARuddIonisationModel(), and operator=().

◆ ~G4DNARuddIonisationModel()

G4DNARuddIonisationModel::~G4DNARuddIonisationModel ( )
override

Definition at line 93 of file G4DNARuddIonisationModel.cc.

94{
95 // Cross section
96
97 std::map<G4String, G4DNACrossSectionDataSet*, std::less<G4String> >::iterator pos;
98 for (pos = tableData.begin(); pos != tableData.end(); ++pos)
99 {
100 G4DNACrossSectionDataSet* table = pos->second;
101 delete table;
102 }
103
104 // The following removal is forbidden since G4VEnergyLossmodel takes care of deletion
105 // Coverity however will signal this as an error
106 // if (fAtomDeexcitation) {delete fAtomDeexcitation;}
107
108}

◆ G4DNARuddIonisationModel() [2/2]

G4DNARuddIonisationModel::G4DNARuddIonisationModel ( const G4DNARuddIonisationModel & )
delete

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 289 of file G4DNARuddIonisationModel.cc.

294{
295 if (verboseLevel > 3)
296 {
297 G4cout << "Calling CrossSectionPerVolume() of G4DNARuddIonisationModel"
298 << G4endl;
299 }
300
301 // Calculate total cross section for model
302
303 if (
304 particleDefinition != protonDef
305 &&
306 particleDefinition != hydrogenDef
307 &&
308 particleDefinition != alphaPlusPlusDef
309 &&
310 particleDefinition != alphaPlusDef
311 &&
312 particleDefinition != heliumDef
313 )
314
315 return 0;
316
317 G4double lowLim = 0;
318
319 if ( particleDefinition == protonDef
320 || particleDefinition == hydrogenDef
321 )
322
323 lowLim = lowEnergyLimitOfModelForZ1;
324
325 if ( particleDefinition == alphaPlusPlusDef
326 || particleDefinition == alphaPlusDef
327 || particleDefinition == heliumDef
328 )
329
330 lowLim = lowEnergyLimitOfModelForZ2;
331
332 G4double highLim = 0;
333 G4double sigma=0;
334
335 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
336
337 const G4String& particleName = particleDefinition->GetParticleName();
338
339 // SI - the following is useless since lowLim is already defined
340 /*
341 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
342 pos1 = lowEnergyLimit.find(particleName);
343
344 if (pos1 != lowEnergyLimit.end())
345 {
346 lowLim = pos1->second;
347 }
348 */
349
350 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
351 pos2 = highEnergyLimit.find(particleName);
352
353 if (pos2 != highEnergyLimit.end())
354 {
355 highLim = pos2->second;
356 }
357
358 if (k <= highLim)
359 {
360 //SI : XS must not be zero otherwise sampling of secondaries method ignored
361
362 if (k < lowLim) k = lowLim;
363
364 //
365
366 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
367 pos = tableData.find(particleName);
368
369 if (pos != tableData.end())
370 {
371 G4DNACrossSectionDataSet* table = pos->second;
372 if (table != nullptr)
373 {
374 sigma = table->FindValue(k);
375 }
376 }
377 else
378 {
379 G4Exception("G4DNARuddIonisationModel::CrossSectionPerVolume","em0002",
380 FatalException,"Model not applicable to particle type.");
381 }
382
383 }
384
385 if (verboseLevel > 2)
386 {
387 G4cout << "__________________________________" << G4endl;
388 G4cout << "G4DNARuddIonisationModel - XS INFO START" << G4endl;
389 G4cout << "Kinetic energy(eV)=" << k/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
390 G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
391 G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
392 //G4cout << " - Cross section per water molecule (cm^-1)="
393 //<< sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
394 G4cout << "G4DNARuddIonisationModel - XS INFO END" << G4endl;
395 }
396
397 return sigma*waterDensity;
398
399}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
double G4double
Definition G4Types.hh:83
std::size_t GetIndex() const

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 112 of file G4DNARuddIonisationModel.cc.

114{
115
116 if (verboseLevel > 3)
117 {
118 G4cout << "Calling G4DNARuddIonisationModel::Initialise()" << G4endl;
119 }
120
121 // Energy limits
122
123 G4String fileProton("dna/sigma_ionisation_p_rudd");
124 G4String fileHydrogen("dna/sigma_ionisation_h_rudd");
125 G4String fileAlphaPlusPlus("dna/sigma_ionisation_alphaplusplus_rudd");
126 G4String fileAlphaPlus("dna/sigma_ionisation_alphaplus_rudd");
127 G4String fileHelium("dna/sigma_ionisation_he_rudd");
128
129 G4DNAGenericIonsManager *instance;
131 protonDef = G4Proton::ProtonDefinition();
132 hydrogenDef = instance->GetIon("hydrogen");
133 alphaPlusPlusDef = G4Alpha::Alpha();
134 alphaPlusDef = instance->GetIon("alpha+");
135 heliumDef = instance->GetIon("helium");
136
137 G4String proton;
138 G4String hydrogen;
139 G4String alphaPlusPlus;
140 G4String alphaPlus;
141 G4String helium;
142
143 G4double scaleFactor = 1 * m*m;
144
145 // LIMITS AND DATA
146
147 // ********************************************************
148
149 proton = protonDef->GetParticleName();
150 tableFile[proton] = fileProton;
151
152 lowEnergyLimit[proton] = lowEnergyLimitForZ1;
153 highEnergyLimit[proton] = 500. * keV;
154
155 // Cross section
156
157 auto tableProton = new G4DNACrossSectionDataSet(new G4LogLogInterpolation,
158 eV,
159 scaleFactor );
160 tableProton->LoadData(fileProton);
161 tableData[proton] = tableProton;
162
163 // ********************************************************
164
165 hydrogen = hydrogenDef->GetParticleName();
166 tableFile[hydrogen] = fileHydrogen;
167
168 lowEnergyLimit[hydrogen] = lowEnergyLimitForZ1;
169 highEnergyLimit[hydrogen] = 100. * MeV;
170
171 // Cross section
172
173 auto tableHydrogen = new G4DNACrossSectionDataSet(new G4LogLogInterpolation,
174 eV,
175 scaleFactor );
176 tableHydrogen->LoadData(fileHydrogen);
177
178 tableData[hydrogen] = tableHydrogen;
179
180 // ********************************************************
181
182 alphaPlusPlus = alphaPlusPlusDef->GetParticleName();
183 tableFile[alphaPlusPlus] = fileAlphaPlusPlus;
184
185 lowEnergyLimit[alphaPlusPlus] = lowEnergyLimitForZ2;
186 highEnergyLimit[alphaPlusPlus] = 400. * MeV;
187
188 // Cross section
189
190 auto tableAlphaPlusPlus = new G4DNACrossSectionDataSet(new G4LogLogInterpolation,
191 eV,
192 scaleFactor );
193 tableAlphaPlusPlus->LoadData(fileAlphaPlusPlus);
194
195 tableData[alphaPlusPlus] = tableAlphaPlusPlus;
196
197 // ********************************************************
198
199 alphaPlus = alphaPlusDef->GetParticleName();
200 tableFile[alphaPlus] = fileAlphaPlus;
201
202 lowEnergyLimit[alphaPlus] = lowEnergyLimitForZ2;
203 highEnergyLimit[alphaPlus] = 400. * MeV;
204
205 // Cross section
206
207 auto tableAlphaPlus = new G4DNACrossSectionDataSet(new G4LogLogInterpolation,
208 eV,
209 scaleFactor );
210 tableAlphaPlus->LoadData(fileAlphaPlus);
211 tableData[alphaPlus] = tableAlphaPlus;
212
213 // ********************************************************
214
215 helium = heliumDef->GetParticleName();
216 tableFile[helium] = fileHelium;
217
218 lowEnergyLimit[helium] = lowEnergyLimitForZ2;
219 highEnergyLimit[helium] = 400. * MeV;
220
221 // Cross section
222
223 auto tableHelium = new G4DNACrossSectionDataSet(new G4LogLogInterpolation,
224 eV,
225 scaleFactor );
226 tableHelium->LoadData(fileHelium);
227 tableData[helium] = tableHelium;
228
229 //
230
231 if (particle==protonDef)
232 {
233 SetLowEnergyLimit(lowEnergyLimit[proton]);
234 SetHighEnergyLimit(highEnergyLimit[proton]);
235 }
236
237 if (particle==hydrogenDef)
238 {
239 SetLowEnergyLimit(lowEnergyLimit[hydrogen]);
240 SetHighEnergyLimit(highEnergyLimit[hydrogen]);
241 }
242
243 if (particle==heliumDef)
244 {
245 SetLowEnergyLimit(lowEnergyLimit[helium]);
246 SetHighEnergyLimit(highEnergyLimit[helium]);
247 }
248
249 if (particle==alphaPlusDef)
250 {
251 SetLowEnergyLimit(lowEnergyLimit[alphaPlus]);
252 SetHighEnergyLimit(highEnergyLimit[alphaPlus]);
253 }
254
255 if (particle==alphaPlusPlusDef)
256 {
257 SetLowEnergyLimit(lowEnergyLimit[alphaPlusPlus]);
258 SetHighEnergyLimit(highEnergyLimit[alphaPlusPlus]);
259 }
260
261 if (isInitialised) { return; }
262
263 // defined stationary mode
265
266 // Initialize water density pointer
268
269 // atomic de-excitation
270 if (!statCode)
271 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
272
274 isInitialised = true;
275
276 if (verboseLevel > 0)
277 {
278 G4cout << "Rudd ionisation model is initialized " << G4endl
279 << "Energy range: "
280 << LowEnergyLimit() / eV << " eV - "
281 << HighEnergyLimit() / keV << " keV for "
282 << particle->GetParticleName()
283 << G4endl;
284 }
285}
G4TemplateRNGHelper< G4long > * G4TemplateRNGHelper< G4long >::instance
static G4Alpha * Alpha()
Definition G4Alpha.cc:83
static G4DNAGenericIonsManager * Instance()
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()
G4ParticleChangeForGamma * fParticleChangeForGamma
static G4EmParameters * Instance()
G4bool DNAStationary() const
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
const G4String & GetParticleName() const
static G4Proton * ProtonDefinition()
Definition G4Proton.cc:85
void SetHighEnergyLimit(G4double)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4double LowEnergyLimit() const
G4double HighEnergyLimit() const
void SetLowEnergyLimit(G4double)

◆ operator=()

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

◆ SampleSecondaries()

void G4DNARuddIonisationModel::SampleSecondaries ( std::vector< G4DynamicParticle * > * fvect,
const G4MaterialCutsCouple * couple,
const G4DynamicParticle * particle,
G4double tmin,
G4double maxEnergy )
overridevirtual

Implements G4VEmModel.

Definition at line 403 of file G4DNARuddIonisationModel.cc.

408{
409 if (verboseLevel > 3)
410 {
411 G4cout << "Calling SampleSecondaries() of G4DNARuddIonisationModel"
412 << G4endl;
413 }
414
415 G4double lowLim = 0;
416 G4double highLim = 0;
417
418 if ( particle->GetDefinition() == protonDef
419 || particle->GetDefinition() == hydrogenDef
420 )
421
422 lowLim = killBelowEnergyForZ1;
423
424 if ( particle->GetDefinition() == alphaPlusPlusDef
425 || particle->GetDefinition() == alphaPlusDef
426 || particle->GetDefinition() == heliumDef
427 )
428
429 lowLim = killBelowEnergyForZ2;
430
431 G4double k = particle->GetKineticEnergy();
432
433 const G4String& particleName = particle->GetDefinition()->GetParticleName();
434
435 // SI - the following is useless since lowLim is already defined
436 /*
437 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
438 pos1 = lowEnergyLimit.find(particleName);
439
440 if (pos1 != lowEnergyLimit.end())
441 {
442 lowLim = pos1->second;
443 }
444 */
445
446 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
447 pos2 = highEnergyLimit.find(particleName);
448
449 if (pos2 != highEnergyLimit.end())
450 {
451 highLim = pos2->second;
452 }
453
454 if (k >= lowLim && k <= highLim)
455 {
456 G4ParticleDefinition* definition = particle->GetDefinition();
457 G4ParticleMomentum primaryDirection = particle->GetMomentumDirection();
458 /*
459 G4double particleMass = definition->GetPDGMass();
460 G4double totalEnergy = k + particleMass;
461 G4double pSquare = k*(totalEnergy+particleMass);
462 G4double totalMomentum = std::sqrt(pSquare);
463 */
464
465 G4int ionizationShell = RandomSelect(k,particleName);
466
468 bindingEnergy = waterStructure.IonisationEnergy(ionizationShell);
469
470 //SI: additional protection if tcs interpolation method is modified
471 if (k<bindingEnergy) return;
472 //
473
474 // SI - For atom. deexc. tagging - 23/05/2017
475 G4int Z = 8;
476 //
477
478 G4double secondaryKinetic = RandomizeEjectedElectronEnergy(definition,k,ionizationShell);
479
480 G4ThreeVector deltaDirection =
481 GetAngularDistribution()->SampleDirectionForShell(particle, secondaryKinetic,
482 Z, ionizationShell,
483 couple->GetMaterial());
484
485 auto dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic);
486 fvect->push_back(dp);
487
488 // Ignored for ions on electrons
489 /*
490 G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
491
492 G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
493 G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
494 G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
495 G4double finalMomentum = std::sqrt(finalPx*finalPx+finalPy*finalPy+finalPz*finalPz);
496 finalPx /= finalMomentum;
497 finalPy /= finalMomentum;
498 finalPz /= finalMomentum;
499
500 G4ThreeVector direction;
501 direction.set(finalPx,finalPy,finalPz);
502
503 fParticleChangeForGamma->ProposeMomentumDirection(direction.unit()) ;
504 */
505
506 fParticleChangeForGamma->ProposeMomentumDirection(primaryDirection);
507
508 // sample deexcitation
509 // here we assume that H_{2}O electronic levels are the same of Oxigen.
510 // this can be considered true with a rough 10% error in energy on K-shell,
511
512 size_t secNumberInit = 0;// need to know at a certain point the energy of secondaries
513 size_t secNumberFinal = 0;// So I'll make the diference and then sum the energies
514
515 G4double scatteredEnergy = k-bindingEnergy-secondaryKinetic;
516
517 // SI: only atomic deexcitation from K shell is considered
518 if((fAtomDeexcitation != nullptr) && ionizationShell == 4)
519 {
520 const G4AtomicShell* shell
521 = fAtomDeexcitation->GetAtomicShell(Z, G4AtomicShellEnumerator(0));
522 secNumberInit = fvect->size();
523 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, 0, 0);
524 secNumberFinal = fvect->size();
525
526 if(secNumberFinal > secNumberInit)
527 {
528 for (size_t i=secNumberInit; i<secNumberFinal; ++i)
529 {
530 //Check if there is enough residual energy
531 if (bindingEnergy >= ((*fvect)[i])->GetKineticEnergy())
532 {
533 //Ok, this is a valid secondary: keep it
534 bindingEnergy -= ((*fvect)[i])->GetKineticEnergy();
535 }
536 else
537 {
538 //Invalid secondary: not enough energy to create it!
539 //Keep its energy in the local deposit
540 delete (*fvect)[i];
541 (*fvect)[i]=nullptr;
542 }
543 }
544 }
545
546 }
547
548 //This should never happen
549 if(bindingEnergy < 0.0)
550 G4Exception("G4DNAEmfietzoglouIonisatioModel1::SampleSecondaries()",
551 "em2050",FatalException,"Negative local energy deposit");
552
553 //bindingEnergy has been decreased
554 //by the amount of energy taken away by deexc. products
555 if (!statCode)
556 {
557 fParticleChangeForGamma->SetProposedKineticEnergy(scatteredEnergy);
558 fParticleChangeForGamma->ProposeLocalEnergyDeposit(bindingEnergy);
559 }
560 else
561 {
562 fParticleChangeForGamma->SetProposedKineticEnergy(k);
563 fParticleChangeForGamma->ProposeLocalEnergyDeposit(k-scatteredEnergy);
564 }
565
566 // debug
567 // k-scatteredEnergy-secondaryKinetic-deexSecEnergy = k-(k-bindingEnergy-secondaryKinetic)-secondaryKinetic-deexSecEnergy =
568 // = k-k+bindingEnergy+secondaryKinetic-secondaryKinetic-deexSecEnergy=
569 // = bindingEnergy-deexSecEnergy
570 // SO deexSecEnergy=0 => LocalEnergyDeposit = bindingEnergy
571
572 // TEST //////////////////////////
573 // if (secondaryKinetic<0) abort();
574 // if (scatteredEnergy<0) abort();
575 // if (k-scatteredEnergy-secondaryKinetic-deexSecEnergy<0) abort();
576 // if (k-scatteredEnergy<0) abort();
577 /////////////////////////////////
578
579 const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
581 ionizationShell,
582 theIncomingTrack);
583 }
584
585 // SI - not useful since low energy of model is 0 eV
586
587 if (k < lowLim)
588 {
589 fParticleChangeForGamma->SetProposedKineticEnergy(0.);
590 fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill);
591 fParticleChangeForGamma->ProposeLocalEnergyDeposit(k);
592 }
593
594}
@ eIonizedMolecule
G4ThreeVector G4ParticleMomentum
CLHEP::Hep3Vector G4ThreeVector
@ fStopAndKill
int G4int
Definition G4Types.hh:85
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition G4Electron.cc:91
const G4Material * GetMaterial() const
virtual G4ThreeVector & SampleDirectionForShell(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, G4int shellID, const G4Material *)
G4VEmAngularDistribution * GetAngularDistribution()
G4double bindingEnergy(G4int A, G4int Z)

◆ SelectStationary()

void G4DNARuddIonisationModel::SelectStationary ( G4bool input)
inline

Definition at line 167 of file G4DNARuddIonisationModel.hh.

168{
169 statCode = input;
170}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNARuddIonisationModel::fParticleChangeForGamma = nullptr
protected

Definition at line 76 of file G4DNARuddIonisationModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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