Geant4 11.2.2
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=nullptr, const G4String &nam="DNARuddIonisationExtendedModel")
 
 ~G4DNARuddIonisationExtendedModel () override
 
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 val)
 
G4double ComputeProbabilityFunction (const G4ParticleDefinition *, G4double kine, G4double deltae, G4int shell)
 
G4DNARuddIonisationExtendedModeloperator= (const G4DNARuddIonisationExtendedModel &right)=delete
 
 G4DNARuddIonisationExtendedModel (const G4DNARuddIonisationExtendedModel &)=delete
 
- 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 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 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)
 
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)
 
void SetLPMFlag (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 51 of file G4DNARuddIonisationExtendedModel.hh.

Constructor & Destructor Documentation

◆ G4DNARuddIonisationExtendedModel() [1/2]

G4DNARuddIonisationExtendedModel::G4DNARuddIonisationExtendedModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "DNARuddIonisationExtendedModel" )
explicit

Definition at line 75 of file G4DNARuddIonisationExtendedModel.cc.

77 : G4VEmModel(nam)
78{
79 fEmCorrections = G4LossTableManager::Instance()->EmCorrections();
80 fGpow = G4Pow::GetInstance();
81 fLowestEnergy = 100*CLHEP::eV;
82 fLimitEnergy = 1*CLHEP::keV;
83
84 // Mark this model as "applicable" for atomic deexcitation
86
87 // Define default angular generator
89}
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
static G4Pow * GetInstance()
Definition G4Pow.cc:41
void SetDeexcitationFlag(G4bool val)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67
void SetAngularDistribution(G4VEmAngularDistribution *)

◆ ~G4DNARuddIonisationExtendedModel()

G4DNARuddIonisationExtendedModel::~G4DNARuddIonisationExtendedModel ( )
override

Definition at line 93 of file G4DNARuddIonisationExtendedModel.cc.

94{
95 if(isFirst) {
96 for(auto & i : xsdata) { delete i; }
97 }
98}

◆ G4DNARuddIonisationExtendedModel() [2/2]

G4DNARuddIonisationExtendedModel::G4DNARuddIonisationExtendedModel ( const G4DNARuddIonisationExtendedModel & )
delete

Member Function Documentation

◆ ComputeProbabilityFunction()

G4double G4DNARuddIonisationExtendedModel::ComputeProbabilityFunction ( const G4ParticleDefinition * p,
G4double kine,
G4double deltae,
G4int shell )

Definition at line 578 of file G4DNARuddIonisationExtendedModel.cc.

580{
581 if (fParticle != p) { SetParticle(p); }
582 G4double bEnergy = (useDNAWaterStructure)
583 ? waterStructure.IonisationEnergy(shell) : Bj[shell];
584 return ProbabilityFunction(e, deltae, bEnergy, shell);
585}
double G4double
Definition G4Types.hh:83

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 250 of file G4DNARuddIonisationExtendedModel.cc.

254{
255 // check if model is applicable for given material
256 G4double density = (material->GetIndex() < fpWaterDensity->size())
257 ? (*fpWaterDensity)[material->GetIndex()] : 0.0;
258 if (0.0 == density) { return 0.0; }
259
260 // ion may be different
261 if (fParticle != part) { SetParticle(part); }
262
263 // initilise mass rate
264 fMassRate = 1.0;
265
266 // ion shoud be stopped - check on kinetic energy and not scaled energy
267 if (kinE < fLowestEnergy) { return DBL_MAX; }
268
269 G4double sigma = 0.;
270
271 // use ion table if available for given energy
272 // for proton, hydrogen, alpha, alpha+, and helium no scaling to proton x-section
273 if (idx == 0 || idx == 1) {
274 sigma = (kinE > fElow) ? xscurrent->FindValue(kinE)
275 : xscurrent->FindValue(fElow)*kinE/fElow;
276
277 // for ions with data above limit energy
278 } else if (idx > 1) {
279 sigma = (kinE > fElow) ? xsdata[idx]->FindValue(kinE)
280 : xsdata[idx]->FindValue(fElow)*kinE/fElow;
281
282 // scaling from proton
283 } else {
284 fMassRate = CLHEP::proton_mass_c2/fMass;
285 G4double e = kinE*fMassRate;
286 sigma = (e > fLowestEnergy) ? xsdata[1]->FindValue(e)
287 : xsdata[1]->FindValue(fLowestEnergy)*e/fLowestEnergy;
288 sigma *= fEmCorrections->EffectiveChargeSquareRatio(part, material, kinE);
289 }
290 sigma *= density;
291
292 if (verbose > 1) {
293 G4cout << "G4DNARuddIonisationExtendedModel for " << part->GetParticleName()
294 << " Ekin(keV)=" << kinE/CLHEP::keV
295 << " sigma(cm^2)=" << sigma/CLHEP::cm2 << G4endl;
296 }
297 return sigma;
298}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4double FindValue(G4double e, G4int componentId=0) const override
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, const G4double kineticEnergy)
std::size_t GetIndex() const
#define DBL_MAX
Definition templates.hh:62

◆ Initialise()

void G4DNARuddIonisationExtendedModel::Initialise ( const G4ParticleDefinition * p,
const G4DataVector &  )
overridevirtual

Implements G4VEmModel.

Definition at line 102 of file G4DNARuddIonisationExtendedModel.cc.

104{
105 if(p != fParticle) { SetParticle(p); }
106
107 // initialisation of static data once
108 if(nullptr == xsdata[0]) {
109 G4AutoLock l(&ionDNAMutex);
110 if(nullptr == xsdata[0]) {
111 isFirst = true;
112 G4String filename("dna/sigma_ionisation_h_rudd");
113 xsdata[0] = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, CLHEP::eV, scaleFactor);
114 xsdata[0]->LoadData(filename);
115
116 filename = "dna/sigma_ionisation_p_rudd";
117 xsdata[1] = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, CLHEP::eV, scaleFactor);
118 xsdata[1]->LoadData(filename);
119
120 filename = "dna/sigma_ionisation_alphaplusplus_rudd";
121 xsdata[2] = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, CLHEP::eV, scaleFactor);
122 xsdata[2]->LoadData(filename);
123
124 filename = "dna/sigma_ionisation_li_rudd";
125 xsdata[3] = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, CLHEP::eV, scaleFactor);
126 xsdata[3]->LoadData(filename);
127
128 filename = "dna/sigma_ionisation_be_rudd";
129 xsdata[4] = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, CLHEP::eV, scaleFactor);
130 xsdata[4]->LoadData(filename);
131
132 filename = "dna/sigma_ionisation_b_rudd";
133 xsdata[5] = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, CLHEP::eV, scaleFactor);
134 xsdata[5]->LoadData(filename);
135
136 filename = "dna/sigma_ionisation_c_rudd";
137 xsdata[6] = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, CLHEP::eV, scaleFactor);
138 xsdata[6]->LoadData(filename);
139
140 filename = "dna/sigma_ionisation_n_rudd";
141 xsdata[7] = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, CLHEP::eV, scaleFactor);
142 xsdata[7]->LoadData(filename);
143
144 filename = "dna/sigma_ionisation_o_rudd";
145 xsdata[8] = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, CLHEP::eV, scaleFactor);
146 xsdata[8]->LoadData(filename);
147
148 filename = "dna/sigma_ionisation_si_rudd";
149 xsdata[14] = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, CLHEP::eV, scaleFactor);
150 xsdata[14]->LoadData(filename);
151
152 filename = "dna/sigma_ionisation_fe_rudd";
153 xsdata[26] = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, CLHEP::eV, scaleFactor);
154 xsdata[26]->LoadData(filename);
155 filename = "dna/sigma_ionisation_alphaplus_rudd";
156 xsalphaplus = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, CLHEP::eV, scaleFactor);
157 xsalphaplus->LoadData(filename);
158
159 filename = "dna/sigma_ionisation_he_rudd";
160 xshelium = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, CLHEP::eV, scaleFactor);
161 xshelium->LoadData(filename);
162 }
163 // to avoid possible threading problem fill this vector only once
164 auto water = G4NistManager::Instance()->FindMaterial("G4_WATER");
165 fpWaterDensity =
167
168 l.unlock();
169 }
170
171 // initialisation once in each thread
172 if(nullptr == fParticleChangeForGamma) {
174 const G4String& pname = fParticle->GetParticleName();
175 if(pname == "proton") {
176 idx = 1;
177 xscurrent = xsdata[1];
178 fElow = fLowestEnergy;
179 } else if(pname == "hydrogen") {
180 idx = 0;
181 xscurrent = xsdata[0];
182 fElow = fLowestEnergy;
183 } else if(pname == "alpha") {
184 idx = 1;
185 xscurrent = xsdata[2];
186 isHelium = true;
187 fElow = fLimitEnergy;
188 } else if(pname == "alpha+") {
189 idx = 1;
190 isHelium = true;
191 xscurrent = xsalphaplus;
192 fElow = fLimitEnergy;
193 // The following values are provided by M. Dingfelder (priv. comm)
194 slaterEffectiveCharge[0]=2.0;
195 slaterEffectiveCharge[1]=2.0;
196 slaterEffectiveCharge[2]=2.0;
197 sCoefficient[0]=0.7;
198 sCoefficient[1]=0.15;
199 sCoefficient[2]=0.15;
200 } else if(pname == "helium") {
201 idx = 0;
202 isHelium = true;
203 fElow = fLimitEnergy;
204 xscurrent = xshelium;
205 slaterEffectiveCharge[0]=1.7;
206 slaterEffectiveCharge[1]=1.15;
207 slaterEffectiveCharge[2]=1.15;
208 sCoefficient[0]=0.5;
209 sCoefficient[1]=0.25;
210 sCoefficient[2]=0.25;
211 } else {
212 isIon = true;
213 }
214 // defined stationary mode
216
217 // initialise atomic de-excitation
218 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
219
220 if (verbose > 0) {
221 G4cout << "### G4DNARuddIonisationExtendedModel::Initialise(..) " << pname
222 << "/n idx=" << idx << " Amass=" << fAmass
223 << " isIon=" << isIon << " isHelium=" << isHelium << G4endl;
224 }
225 }
226}
G4bool LoadData(const G4String &argFileName) override
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 G4EmParameters * Instance()
G4bool DNAStationary() const
G4VAtomDeexcitation * AtomDeexcitation()
G4Material * FindMaterial(const G4String &name) const
static G4NistManager * Instance()
const G4String & GetParticleName() const
G4ParticleChangeForGamma * GetParticleChangeForGamma()

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 303 of file G4DNARuddIonisationExtendedModel.cc.

307{
308 const G4ParticleDefinition* pd = dpart->GetDefinition();
309 if (fParticle != pd) { SetParticle(pd); }
310
311 // stop ion with energy below low energy limit
312 G4double kinE = dpart->GetKineticEnergy();
313 // ion shoud be stopped - check on kinetic energy and not scaled energy
314 if (kinE <= fLowestEnergy) {
318 return;
319 }
320
321 G4int shell = SelectShell(kinE);
322 G4double bindingEnergy = (useDNAWaterStructure)
323 ? waterStructure.IonisationEnergy(shell) : Bj[shell];
324
325 //Si: additional protection if tcs interpolation method is modified
326 if (kinE < bindingEnergy) return;
327
328 G4double esec = SampleElectronEnergy(kinE, bindingEnergy, shell);
329 G4double esum = 0.0;
330
331 // sample deexcitation
332 // here we assume that H2O electronic levels are the same as Oxygen.
333 // this can be considered true with a rough 10% error in energy on K-shell,
334 G4int Z = 8;
335 G4ThreeVector deltaDir =
336 GetAngularDistribution()->SampleDirectionForShell(dpart, esec, Z, shell, couple->GetMaterial());
337
338 // SI: only atomic deexcitation from K shell is considered
339 if(fAtomDeexcitation != nullptr && shell == 4) {
340 auto as = G4AtomicShellEnumerator(0);
341 auto ashell = fAtomDeexcitation->GetAtomicShell(Z, as);
342 fAtomDeexcitation->GenerateParticles(fvect, ashell, Z, 0, 0);
343
344 // compute energy sum from de-excitation
345 std::size_t nn = fvect->size();
346 for (std::size_t i=0; i<nn; ++i) {
347 esum += (*fvect)[i]->GetKineticEnergy();
348 }
349 }
350 // check energy balance
351 // remaining excitation energy of water molecule
352 G4double exc = bindingEnergy - esum;
353
354 // remaining projectile energy
355 G4double scatteredEnergy = kinE - bindingEnergy - esec;
356 if(scatteredEnergy < -tolerance || exc < -tolerance) {
357 G4cout << "G4DNARuddIonisationExtendedModel::SampleSecondaries: "
358 << "negative final E(keV)=" << scatteredEnergy/CLHEP::keV << " Ein(keV)="
359 << kinE/CLHEP::keV << " " << pd->GetParticleName()
360 << " Edelta(keV)=" << esec/CLHEP::keV << " MeV, Exc(keV)=" << exc/CLHEP::keV
361 << G4endl;
362 }
363
364 // projectile
365 if (!statCode) {
368 } else {
370 fParticleChangeForGamma->ProposeLocalEnergyDeposit(kinE - scatteredEnergy);
371 }
372
373 // delta-electron
374 auto dp = new G4DynamicParticle(G4Electron::Electron(), deltaDir, esec);
375 fvect->push_back(dp);
376
377 // create radical
378 const G4Track* theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
380 theIncomingTrack);
381}
@ eIonizedMolecule
@ fStopButAlive
int G4int
Definition G4Types.hh:85
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition G4Electron.cc:91
const G4Material * GetMaterial() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
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()
void ProposeTrackStatus(G4TrackStatus status)
const G4Track * GetCurrentTrack() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double bindingEnergy(G4int A, G4int Z)

◆ SelectStationary()

void G4DNARuddIonisationExtendedModel::SelectStationary ( G4bool val)
inline

Definition at line 74 of file G4DNARuddIonisationExtendedModel.hh.

74{ statCode = val; };

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNARuddIonisationExtendedModel::fParticleChangeForGamma {nullptr}
protected

Definition at line 120 of file G4DNARuddIonisationExtendedModel.hh.

120{nullptr};

Referenced by Initialise(), and SampleSecondaries().


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