Geant4 10.7.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=0, const G4String &nam="DNARuddIonisationModel")
 
virtual ~G4DNARuddIonisationModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
void SelectStationary (G4bool input)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
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 CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=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 *, G4double &eloss, G4double &niel, G4double length)
 
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 ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual 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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
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 LPMFlag () 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 SetLPMFlag (G4bool val)
 
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 *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

G4ParticleChangeForGammafParticleChangeForGamma
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

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

G4DNARuddIonisationModel::G4DNARuddIonisationModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNARuddIonisationModel" 
)

Definition at line 45 of file G4DNARuddIonisationModel.cc.

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

◆ ~G4DNARuddIonisationModel()

G4DNARuddIonisationModel::~G4DNARuddIonisationModel ( )
virtual

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}

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 287 of file G4DNARuddIonisationModel.cc.

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

◆ Initialise()

void G4DNARuddIonisationModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector  
)
virtual

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;
132 G4ParticleDefinition* hydrogenDef = instance->GetIon("hydrogen");
133 G4ParticleDefinition* alphaPlusPlusDef = instance->GetIon("alpha++");
134 G4ParticleDefinition* alphaPlusDef = instance->GetIon("alpha+");
135 G4ParticleDefinition* heliumDef = instance->GetIon("helium");
136
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
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
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
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
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
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( verboseLevel>0 )
262 {
263 G4cout << "Rudd ionisation model is initialized " << G4endl
264 << "Energy range: "
265 << LowEnergyLimit() / eV << " eV - "
266 << HighEnergyLimit() / keV << " keV for "
267 << particle->GetParticleName()
268 << G4endl;
269 }
270
271 // Initialize water density pointer
273
274 //
275
276 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
277
278 if (isInitialised)
279 { return;}
281 isInitialised = true;
282
283}
virtual G4bool LoadData(const G4String &argFileName)
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()
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:651
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 404 of file G4DNARuddIonisationModel.cc.

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

◆ SelectStationary()

void G4DNARuddIonisationModel::SelectStationary ( G4bool  input)
inline

Definition at line 162 of file G4DNARuddIonisationModel.hh.

163{
164 statCode = input;
165}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNARuddIonisationModel::fParticleChangeForGamma
protected

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