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

#include <G4DNARuddIonisationExtendedModel.hh>

+ Inheritance diagram for G4DNARuddIonisationExtendedModel:

Public Member Functions

 G4DNARuddIonisationExtendedModel (const G4ParticleDefinition *p=0, const G4String &nam="DNARuddIonisationExtendedModel")
 
virtual ~G4DNARuddIonisationExtendedModel ()
 
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 *, 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 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 *)
 
G4bool IsLocked () const
 
void SetLocked (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
 
size_t currentCoupleIndex = 0
 
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 G4DNARuddIonisationExtendedModel.hh.

Constructor & Destructor Documentation

◆ G4DNARuddIonisationExtendedModel()

G4DNARuddIonisationExtendedModel::G4DNARuddIonisationExtendedModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNARuddIonisationExtendedModel" 
)

Definition at line 55 of file G4DNARuddIonisationExtendedModel.cc.

57:G4VEmModel(nam),isInitialised(false)
58{
59 slaterEffectiveCharge[0]=0.;
60 slaterEffectiveCharge[1]=0.;
61 slaterEffectiveCharge[2]=0.;
62 sCoefficient[0]=0.;
63 sCoefficient[1]=0.;
64 sCoefficient[2]=0.;
65
66 lowEnergyLimitForA[1] = 0 * eV;
67 lowEnergyLimitForA[2] = 0 * eV;
68 lowEnergyLimitForA[3] = 0 * eV;
69 lowEnergyLimitOfModelForA[1] = 100 * eV;
70 lowEnergyLimitOfModelForA[4] = 1 * keV;
71 lowEnergyLimitOfModelForA[5] = 0.5 * MeV; // For A = 3 or above, limit is MeV/uma
72 killBelowEnergyForA[1] = lowEnergyLimitOfModelForA[1];
73 killBelowEnergyForA[4] = lowEnergyLimitOfModelForA[4];
74 killBelowEnergyForA[5] = lowEnergyLimitOfModelForA[5];
75
76 verboseLevel= 0;
77 // Verbosity scale:
78 // 0 = nothing
79 // 1 = warning for energy non-conservation
80 // 2 = details of energy budget
81 // 3 = calculation of cross sections, file openings, sampling of atoms
82 // 4 = entering in methods
83
84 if( verboseLevel>0 )
85 {
86 G4cout << "Rudd ionisation model is constructed " << G4endl;
87 }
88
89 // Define default angular generator
91
92 // Mark this model as "applicable" for atomic deexcitation
94
95 // Selection of stationary mode
96
97 statCode = false;
98}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:802
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:607

◆ ~G4DNARuddIonisationExtendedModel()

G4DNARuddIonisationExtendedModel::~G4DNARuddIonisationExtendedModel ( )
virtual

Definition at line 102 of file G4DNARuddIonisationExtendedModel.cc.

103{
104 if(isIon) {
105 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
106 for (pos = tableData.begin(); pos != tableData.end(); ++pos)
107 {
108 G4DNACrossSectionDataSet* table = pos->second;
109 delete table;
110 }
111 } else {
112 delete mainTable;
113 }
114}

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 297 of file G4DNARuddIonisationExtendedModel.cc.

302{
303 if (verboseLevel > 3)
304 G4cout << "Calling CrossSectionPerVolume() of G4DNARuddIonisationExtendedModel" << G4endl;
305
306 currParticle = GetDNAIonParticleDefinition(partDef);
307 currentScaledEnergy = k;
308 G4double e = k;
309 G4double q2 = 1.0;
310 currentTable = mainTable;
311
312 if (isIon){
313 if (currParticle == nullptr) {//not DNA particle
314 currentScaledEnergy *= massC12/partDef->GetPDGMass();
315 G4double q = partDef->GetPDGCharge()/(eplus*6);
316 q2 *= q*q;
317 e = currentScaledEnergy;
318 currParticle = carbonDef;
319 }
320 G4String pname = currParticle->GetParticleName();
321 auto goodTable = tableData.find(pname);
322 currentTable = goodTable->second;
323 }
324 // below low the ion should be stopped
325 if (currentScaledEnergy < localMinEnergy) { return DBL_MAX; }
326
327 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
328 G4double sigma = 0.0;
329 if (nullptr != currentTable) {
330 sigma = currentTable->FindValue(e)*q2;
331 } else {
332 G4cout << "G4DNARuddIonisationExtendedModel - no data table for "
333 << partDef->GetParticleName() << G4endl;
334 G4Exception("G4DNARuddIonisationExtendedModel::CrossSectionPerVolume(...)","em0002",
335 FatalException,"Data table is not available for the model.");
336 }
337 if (verboseLevel > 2)
338 {
339 G4cout << "__________________________________" << G4endl;
340 G4cout << "G4DNARuddIonisationExtendedModel - XS INFO START" << G4endl;
341 G4cout << "Kinetic energy(eV)=" << k/eV << " particle : " << partDef->GetParticleName() << G4endl;
342 G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
343 G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
344 G4cout << "G4DNARuddIonisationExtendedModel - XS INFO END" << G4endl;
345 }
346 return sigma*waterDensity;
347}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
double G4double
Definition: G4Types.hh:83
virtual G4double FindValue(G4double e, G4int componentId=0) const
size_t GetIndex() const
Definition: G4Material.hh:255
const G4String & GetParticleName() const
#define DBL_MAX
Definition: templates.hh:62

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 118 of file G4DNARuddIonisationExtendedModel.cc.

120{
121 if (isInitialised) { return; }
122 if (verboseLevel > 3)
123 G4cout << "Calling G4DNARuddIonisationExtendedModel::Initialise()" << G4endl;
124
125 // Energy limits
126 G4String fileProton("dna/sigma_ionisation_p_rudd");
127 G4String fileHydrogen("dna/sigma_ionisation_h_rudd");
128 G4String fileAlphaPlusPlus("dna/sigma_ionisation_alphaplusplus_rudd");
129 G4String fileAlphaPlus("dna/sigma_ionisation_alphaplus_rudd");
130 G4String fileHelium("dna/sigma_ionisation_he_rudd");
131 G4String fileCarbon("dna/sigma_ionisation_c_rudd");
132 G4String fileNitrogen("dna/sigma_ionisation_n_rudd");
133 G4String fileOxygen("dna/sigma_ionisation_o_rudd");
134 G4String fileSilicon("dna/sigma_ionisation_si_rudd");
135 G4String fileIron("dna/sigma_ionisation_fe_rudd");
136
137 G4String pname = particle->GetParticleName();
138
139 G4DNAGenericIonsManager *instance;
141 protonDef = G4Proton::ProtonDefinition();
142 hydrogenDef = instance->GetIon("hydrogen");
143 alphaPlusPlusDef = G4Alpha::Alpha();
144 alphaPlusDef = instance->GetIon("alpha+");
145 heliumDef = instance->GetIon("helium");
146
147 carbonDef = instance->GetIon("carbon");
148 nitrogenDef = instance->GetIon("nitrogen");
149 oxygenDef = instance->GetIon("oxygen");
150 siliconDef = instance->GetIon("silicon");
151 ironDef = instance->GetIon("iron");
152
153 G4String carbon;
154 G4String nitrogen;
155 G4String oxygen;
156 G4String silicon;
157 G4String iron;
158
159 G4double scaleFactor = 1 * m*m;
160 massC12 = carbonDef->GetPDGMass();
161
162 // LIMITS AND DATA
163
164 // **********************************************************************************************
165
166 if(pname == "proton") {
167 localMinEnergy = lowEnergyLimitForA[1];
168
169 // Cross section
170 mainTable = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV, scaleFactor );
171 mainTable->LoadData(fileProton);
172
173 // **********************************************************************************************
174
175 } else if(pname == "hydrogen") {
176
177 localMinEnergy = lowEnergyLimitForA[1];
178
179 // Cross section
180 mainTable = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV, scaleFactor );
181 mainTable->LoadData(fileHydrogen);
182
183 // **********************************************************************************************
184
185 } else if(pname == "alpha") {
186
187 localMinEnergy = lowEnergyLimitForA[4];
188
189 // Cross section
190 mainTable = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV, scaleFactor );
191 mainTable->LoadData(fileAlphaPlusPlus);
192
193 // **********************************************************************************************
194
195 } else if(pname == "alpha+") {
196
197 localMinEnergy = lowEnergyLimitForA[4];
198
199 // Cross section
200 mainTable = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV, scaleFactor );
201 mainTable->LoadData(fileAlphaPlus);
202
203 // **********************************************************************************************
204
205 } else if(pname == "helium") {
206
207 localMinEnergy = lowEnergyLimitForA[4];
208
209 // Cross section
210 mainTable = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV, scaleFactor );
211 mainTable->LoadData(fileHelium);
212
213 // **********************************************************************************************
214
215 } else if(pname == "GenericIon") {
216
217 isIon = true;
218 carbon = carbonDef->GetParticleName();
219 localMinEnergy = lowEnergyLimitForA[5]*massC12/proton_mass_c2;
220
221 // Cross section
222 mainTable = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV, scaleFactor );
223 mainTable->LoadData(fileCarbon);
224
225 tableData[carbon] = mainTable;
226
227 // **********************************************************************************************
228
229 oxygen = oxygenDef->GetParticleName();
230 tableFile[oxygen] = fileOxygen;
231
232 // Cross section
234 eV, scaleFactor );
235 tableOxygen->LoadData(fileOxygen);
236 tableData[oxygen] = tableOxygen;
237
238 // **********************************************************************************************
239
240 nitrogen = nitrogenDef->GetParticleName();
241 tableFile[nitrogen] = fileNitrogen;
242
243 // Cross section
245 eV, scaleFactor );
246 tableNitrogen->LoadData(fileNitrogen);
247 tableData[nitrogen] = tableNitrogen;
248
249 // **********************************************************************************************
250
251 silicon = siliconDef->GetParticleName();
252 tableFile[silicon] = fileSilicon;
253
254 // Cross section
256 eV, scaleFactor );
257 tableSilicon->LoadData(fileSilicon);
258 tableData[silicon] = tableSilicon;
259
260 // **********************************************************************************************
261
262 iron = ironDef->GetParticleName();
263 tableFile[iron] = fileIron;
264
265 // Cross section
266
268 eV, scaleFactor );
269 tableIron->LoadData(fileIron);
270 tableData[iron] = tableIron;
271 }
272 // **********************************************************************************************
273
274 if( verboseLevel>0 )
275 {
276 G4cout << "Rudd ionisation model is initialized " << G4endl
277 << "Energy range for model: "
278 << LowEnergyLimit() / keV << " keV - "
279 << HighEnergyLimit() / keV << " keV for "
280 << pname
281 << " internal low energy limit E(keV)=" << localMinEnergy / keV
282 << G4endl;
283 }
284
285 // Initialize water density pointer
287
288 //
289 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
290
292 isInitialised = true;
293}
static G4Alpha * Alpha()
Definition: G4Alpha.cc:88
virtual G4bool LoadData(const G4String &argFileName)
static G4DNAGenericIonsManager * Instance(void)
G4ParticleDefinition * GetIon(const G4String &name)
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:691
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:87
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:124
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 351 of file G4DNARuddIonisationExtendedModel.cc.

356{
357 if (verboseLevel > 3)
358 G4cout << "Calling SampleSecondaries() of G4DNARuddIonisationExtendedModel" << G4endl;
359
360 // stop ion with energy below low energy limit
361 G4double k = particle->GetKineticEnergy();
362 if (currentScaledEnergy < localMinEnergy) {
366 return;
367 }
368
369 // sampling of final state
370 G4ParticleDefinition* definition = particle->GetDefinition();
371 const G4ThreeVector& primaryDirection = particle->GetMomentumDirection();
372
373 G4int ionizationShell = RandomSelect(currentScaledEnergy);
374
375 // sample deexcitation
376 // here we assume that H_{2}O electronic levels are the same as Oxygen.
377 // this can be considered true with a rough 10% error in energy on K-shell,
378
379 G4double bindingEnergy = waterStructure.IonisationEnergy(ionizationShell);
380
381 //SI: additional protection if tcs interpolation method is modified
382 if (k < bindingEnergy) return;
383 //
384
385 G4double secondaryKinetic = RandomizeEjectedElectronEnergy(definition,k,ionizationShell);
386
387 // is ionisation possible?
388 G4double scatteredEnergy = k - bindingEnergy - secondaryKinetic;
389 if(scatteredEnergy < 0.0) { return; }
390
391 G4int Z = 8;
392
393 G4ThreeVector deltaDirection =
394 GetAngularDistribution()->SampleDirectionForShell(particle, secondaryKinetic,
395 Z, ionizationShell,
396 couple->GetMaterial());
397
398 G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic) ;
399 fvect->push_back(dp);
400
402
403 size_t secNumberInit = 0;// need to know at a certain point the energy of secondaries
404 size_t secNumberFinal = 0;// So I'll make the diference and then sum the energies
405
406 // SI: only atomic deexcitation from K shell is considered
407 if(fAtomDeexcitation != nullptr && ionizationShell == 4)
408 {
409 const G4AtomicShell* shell
410 = fAtomDeexcitation->GetAtomicShell(Z, G4AtomicShellEnumerator(0));
411 secNumberInit = fvect->size();
412 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, 0, 0);
413 secNumberFinal = fvect->size();
414
415 if(secNumberFinal > secNumberInit)
416 {
417 for (size_t i=secNumberInit; i<secNumberFinal; ++i)
418 {
419 //Check if there is enough residual energy
420 if (bindingEnergy >= ((*fvect)[i])->GetKineticEnergy())
421 {
422 //Ok, this is a valid secondary: keep it
423 bindingEnergy -= ((*fvect)[i])->GetKineticEnergy();
424 }
425 else
426 {
427 //Invalid secondary: not enough energy to create it!
428 //Keep its energy in the local deposit
429 delete (*fvect)[i];
430 (*fvect)[i]=0;
431 }
432 }
433 }
434 }
435
436 //bindingEnergy has been decreased
437 //by the amount of energy taken away by deexc. products
438 if (!statCode)
439 {
442 }
443 else
444 {
447 }
448
449 // TEST //////////////////////////
450 // if (secondaryKinetic<0) abort();
451 // if (scatteredEnergy<0) abort();
452 // if (k-scatteredEnergy-secondaryKinetic-deexSecEnergy<0) abort();
453 // if (k-scatteredEnergy<0) abort();
454 /////////////////////////////////
455
456 const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
458 ionizationShell,
459 theIncomingTrack);
460 // SI - not useful since low energy of model is 0 eV
461}
G4AtomicShellEnumerator
@ eIonizedMolecule
@ fStopButAlive
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
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
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
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:600
void ProposeTrackStatus(G4TrackStatus status)
const G4Track * GetCurrentTrack() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double bindingEnergy(G4int A, G4int Z)

◆ SelectStationary()

void G4DNARuddIonisationExtendedModel::SelectStationary ( G4bool  input)
inline

Definition at line 184 of file G4DNARuddIonisationExtendedModel.hh.

185{
186 statCode = input;
187}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNARuddIonisationExtendedModel::fParticleChangeForGamma = nullptr
protected

Definition at line 73 of file G4DNARuddIonisationExtendedModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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