Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNARuddIonisationDynamicModel.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// Created 11.02.2025 V.Ivanchenko & M. Vologzhin
27// on base of previous Rudd models
28//
29// Russian Goverment grant No 075-15-2024-667 23.08.2024
30//
31
34#include "G4SystemOfUnits.hh"
37#include "G4LossTableManager.hh"
38#include "G4EmCorrections.hh"
39#include "G4NistManager.hh"
42
44
46#include "G4MatUtils.hh"
48#include "G4EmParameters.hh"
49#include "G4NistManager.hh"
50
51#include "G4IonTable.hh"
52#include "G4DNARuddAngle.hh"
53#include "G4DeltaAngle.hh"
54#include "G4Exp.hh"
55#include "G4Log.hh"
56#include "G4Pow.hh"
57#include "G4Alpha.hh"
58#include "G4Proton.hh"
59#include "G4Electron.hh"
60
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
62
63G4ExtendedPhysicsVector* G4DNARuddIonisationDynamicModel::xsdata_alpha = nullptr;
64G4ExtendedPhysicsVector* G4DNARuddIonisationDynamicModel::xsdata_alphap = nullptr;
65G4ExtendedPhysicsVector* G4DNARuddIonisationDynamicModel::xsdata_hydrogen = nullptr;
66G4ExtendedPhysicsVector* G4DNARuddIonisationDynamicModel::xsdata_helium = nullptr;
67G4ExtendedPhysicsVector* G4DNARuddIonisationDynamicModel::xsdata_p = nullptr;
68const std::vector<G4double>* G4DNARuddIonisationDynamicModel::fpWaterDensity = nullptr;
69
70namespace
71{
72 const G4double scaleFactor = CLHEP::m*CLHEP::m;
73 const G4double tolerance = 1*CLHEP::eV;
74 const G4double Ry = 13.6*CLHEP::eV;
75
76 // Following values provided by M. Dingfelder (priv. comm)
77 const G4double Bj[5] = {12.60*CLHEP::eV, 14.70*CLHEP::eV, 18.40*CLHEP::eV,
78 32.20*CLHEP::eV, 539*CLHEP::eV};
79}
80
81//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
82
84 const G4String& nam)
85 : G4VEmModel(nam)
86{
87 fEmCorrections = G4LossTableManager::Instance()->EmCorrections();
88 fGpow = G4Pow::GetInstance();
89 fLowestEnergy = 100*CLHEP::eV;
90 fAbsorptionEnergy = 50*CLHEP::eV;
91
92 // Mark this model as "applicable" for atomic deexcitation
94
95 // Define default angular generator
97
98 if (nullptr == xsdata_p) {
99 isFirst = true;
100 LoadData();
101 }
102}
103
104//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
105
107{
108 if (isFirst) {
109 delete xsdata_alpha;
110 delete xsdata_alphap;
111 delete xsdata_p;
112 delete xsdata_hydrogen;
113 delete xsdata_helium;
114 }
115}
116
117//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
118
119void G4DNARuddIonisationDynamicModel::LoadData()
120{
121 // initialisation of static data once
122 const G4String& dirpath = G4EmParameters::Instance()->GetDirLEDATA();
123
124 G4String filename = "dna/sigma_ionisation_p_rudd.dat";
125 xsdata_p =
126 G4MatUtils::BuildExtendedVector(dirpath, filename, 5, 645, CLHEP::eV, scaleFactor);
127
128 filename = "dna/sigma_ionisation_alphaplusplus_rudd.dat";
129 xsdata_alpha =
130 G4MatUtils::BuildExtendedVector(dirpath, filename, 5, 540, CLHEP::eV, scaleFactor);
131
132 filename = "dna/sigma_ionisation_alphaplus_rudd.dat";
133 xsdata_alphap =
134 G4MatUtils::BuildExtendedVector(dirpath, filename, 5, 540, CLHEP::eV, scaleFactor);
135
136 filename = "dna/sigma_ionisation_h_rudd.dat";
137 xsdata_hydrogen =
138 G4MatUtils::BuildExtendedVector(dirpath, filename, 5, 600, CLHEP::eV, scaleFactor);
139
140 filename = "dna/sigma_ionisation_he_rudd.dat";
141 xsdata_helium =
142 G4MatUtils::BuildExtendedVector(dirpath, filename, 5, 540, CLHEP::eV, scaleFactor);
143
144 // to avoid possible threading problem fill this vector only once
145 auto water = G4NistManager::Instance()->FindMaterial("G4_WATER");
146 fpWaterDensity =
148}
149
150//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
151
153 const G4DataVector&)
154{
155 if (p != fParticle) { SetParticle(p); }
156
157 // particle change object may be externally set
158 if (nullptr == fParticleChangeForGamma) {
160 }
161 const G4String& pname = p->GetParticleName();
162
163 // the same definition of generic ion as in G4VEmProcess class
164 if (p->GetParticleType() == "nucleus" && p->GetParticleSubType() == "generic") {
165 if (pname != "deuteron" && pname != "triton" &&
166 pname != "He3" && pname != "alpha" && pname != "alpha+" &&
167 pname != "helium" && pname != "hydrogen") {
168 isIon = true;
169 }
170 }
171
172 // initialisation once in each thread
173 if (!isInitialised) {
174 isInitialised = true;
175 xsdata = xsdata_p;
176
177 if (pname == "helium") {
178 isHelium = true;
179 xsdata = xsdata_helium;
180 slaterEffectiveCharge[0]=1.7;
181 slaterEffectiveCharge[1]=1.15;
182 slaterEffectiveCharge[2]=1.15;
183 sCoefficient[0]=0.5;
184 sCoefficient[1]=0.25;
185 sCoefficient[2]=0.25;
186 fLowestEnergy = 1*CLHEP::keV;
187 } else if (pname == "alpha+") {
188 isHelium = true;
189 xsdata = xsdata_alphap;
190 // The following values are provided by M. Dingfelder (priv. comm)
191 slaterEffectiveCharge[0]=2.0;
192 slaterEffectiveCharge[1]=2.0;
193 slaterEffectiveCharge[2]=2.0;
194 sCoefficient[0]=0.7;
195 sCoefficient[1]=0.15;
196 sCoefficient[2]=0.15;
197 } else if (pname == "alpha") {
198 isHelium = true;
199 xsdata = xsdata_alpha;
200 } else if (pname == "hydrogen") {
201 xsdata = xsdata_hydrogen;
202 } else if (pname != "proton") {
203 isIon = true;
204 }
205 if (isHelium) { fLowestEnergy = 1*CLHEP::keV; }
206
207 // defined stationary mode
209
210 // initialise atomic de-excitation
211 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
212
213 // chemistry
215 if (chem->IsChemistryActivated()) {
216 fChemistry = chem;
217 }
218
219 InitialiseIntegrator(0.1, 0.25, 1.05, 4*CLHEP::eV, 0.2*CLHEP::eV, 10*CLHEP::keV);
220
221 if (verbose > 0) {
222 G4cout << "### G4DNARuddIonisationDynamicModel::Initialise(..) "
223 << fParticle->GetParticleName() << G4endl;
224 }
225 }
226}
227
228//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
229
230void G4DNARuddIonisationDynamicModel::SetParticle(const G4ParticleDefinition* p)
231{
232 fParticle = p;
233 fMass = p->GetPDGMass();
234 if (isIon) {
235 fMassRate = CLHEP::proton_mass_c2/fMass;
236 fMass = CLHEP::proton_mass_c2;
237 }
238}
239
240//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
241
243{
244 fTrack = track;
245}
246
247//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
248
251 const G4ParticleDefinition* part,
252 G4double kinE,
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 // check on kinetic energy (not scaled energy) to stop low-energy ion
261 if (kinE < fAbsorptionEnergy) { return DBL_MAX; }
262
263 // ion may be different
264 if (fParticle != part) { SetParticle(part); }
265
266 // cross section for scaled energy
267 G4double e = kinE*fMassRate;
268
269 G4double sigma = (e > fLowestEnergy) ? xsdata->LogLogValue(e, idx)
270 : xsdata->LogLogValue(fLowestEnergy, idx) * e / fLowestEnergy;
271
272 sigma *= density;
273 if (isIon) {
274 sigma *= fEmCorrections->EffectiveChargeSquareRatio(part, material, kinE);
275 }
276
277 if (verbose > 1) {
278 G4cout << "G4DNARuddIonisationDynamicModel for " << part->GetParticleName()
279 << " Ekin(keV)=" << kinE/CLHEP::keV
280 << " sigma(cm^2)=" << sigma/CLHEP::cm2 << G4endl;
281 }
282 return sigma;
283}
284
285//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
286
287void
288G4DNARuddIonisationDynamicModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
289 const G4MaterialCutsCouple* couple,
290 const G4DynamicParticle* dpart,
292{
293 const G4ParticleDefinition* pd = dpart->GetDefinition();
294 if (fParticle != pd) { SetParticle(pd); }
295
296 // stop ion with energy below low energy limit
297 G4double kinE = dpart->GetKineticEnergy();
298 // ion shoud be stopped - check on kinetic energy and not scaled energy
299 if (kinE <= fAbsorptionEnergy) {
300 fParticleChangeForGamma->SetProposedKineticEnergy(0.);
301 fParticleChangeForGamma->ProposeTrackStatus(fStopButAlive);
302 fParticleChangeForGamma->ProposeLocalEnergyDeposit(kinE);
303 return;
304 }
305
306 fScaledEnergy = kinE*fMassRate;
307 fSelectedShell = SelectShell();
308 G4double bindingEnergy = (useDNAWaterStructure)
309 ? waterStructure.IonisationEnergy(fSelectedShell) : Bj[fSelectedShell];
310
311 //Si: additional protection if tcs interpolation method is modified
312 if (kinE < bindingEnergy) { return; }
313
314 G4double esec = SampleElectronEnergy();
315 G4double esum = 0.0;
316
317 // sample deexcitation
318 // here we assume that H2O electronic levels are the same as Oxygen.
319 // this can be considered true with a rough 10% error in energy on K-shell,
320 G4int Z = 8;
321 G4ThreeVector deltaDir =
323 fSelectedShell,
324 couple->GetMaterial());
325
326 // SI: only atomic deexcitation from K shell is considered
327 if (fAtomDeexcitation != nullptr && fSelectedShell == 4) {
328 auto as = G4AtomicShellEnumerator(0);
329 auto ashell = fAtomDeexcitation->GetAtomicShell(Z, as);
330 fAtomDeexcitation->GenerateParticles(fvect, ashell, Z, 0, 0);
331
332 // compute energy sum from de-excitation
333 for (auto const & ptr : *fvect) {
334 esum += ptr->GetKineticEnergy();
335 }
336 }
337 // check energy balance
338 // remaining excitation energy of water molecule
339 G4double exc = std::max(bindingEnergy - esum, 0.0);
340
341 // remaining projectile energy
342 G4double scatteredEnergy = kinE - bindingEnergy - esec;
343 if(scatteredEnergy < -tolerance || exc < -tolerance) {
344 G4cout << "G4DNARuddIonisationDynamicModel::SampleSecondaries: "
345 << "negative final E(keV)=" << scatteredEnergy/CLHEP::keV << " Ein(keV)="
346 << kinE/CLHEP::keV << " " << pd->GetParticleName()
347 << " Edelta(keV)=" << esec/CLHEP::keV << " MeV, Exc(keV)=" << exc/CLHEP::keV
348 << G4endl;
349 }
350 scatteredEnergy = std::max(scatteredEnergy, 0.0);
351 /*
352 G4cout << "Eprim(keV)=" << kinE/CLHEP::keV << " Efin(keV)=" << scatteredEnergy/CLHEP::keV
353 << " Esec(keV)=" << esec/CLHEP::keV << " Exc(keV)=" << exc/CLHEP::keV
354 << " tolerance(keV)=" << tolerance/CLHEP::keV << G4endl;
355 */
356 // projectile
357 if (!statCode) {
358 fParticleChangeForGamma->SetProposedKineticEnergy(scatteredEnergy);
359 fParticleChangeForGamma->ProposeLocalEnergyDeposit(exc);
360 } else {
361 fParticleChangeForGamma->SetProposedKineticEnergy(kinE);
362 fParticleChangeForGamma->ProposeLocalEnergyDeposit(kinE - scatteredEnergy);
363 }
364
365 // delta-electron
366 auto dp = new G4DynamicParticle(G4Electron::Electron(), deltaDir, esec);
367 fvect->push_back(dp);
368
369 // create radical
370 if (nullptr != fChemistry) {
371 fChemistry->CreateWaterMolecule(eIonizedMolecule, fSelectedShell, fTrack);
372 }
373}
374
375//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
376
377G4int G4DNARuddIonisationDynamicModel::SelectShell()
378{
379 return xsdata->SampleReactionChannel(fScaledEnergy, G4UniformRand(), idx);
380}
381
382//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
383
385G4DNARuddIonisationDynamicModel::MaxEnergy()
386{
387 // kinematic limit
388 G4double tau = fScaledEnergy/fMass;
389 G4double gam = 1.0 + tau;
390 G4double emax = 2.0*CLHEP::electron_mass_c2*tau*(tau + 2.0);
391
392 // Initialisation of sampling
393 G4double A1, B1, C1, D1, E1, A2, B2, C2, D2;
394 if (fSelectedShell == 4) {
395 //Data For Liquid Water K SHELL from Dingfelder (Protons in Water)
396 A1 = 1.25;
397 B1 = 0.5;
398 C1 = 1.00;
399 D1 = 1.00;
400 E1 = 3.00;
401 A2 = 1.10;
402 B2 = 1.30;
403 C2 = 1.00;
404 D2 = 0.00;
405 alphaConst = 0.66;
406 } else {
407 //Data For Liquid Water from Dingfelder (Protons in Water)
408 A1 = 1.02;
409 B1 = 82.0;
410 C1 = 0.45;
411 D1 = -0.80;
412 E1 = 0.38;
413 A2 = 1.07;
414 // Value provided by M. Dingfelder (priv. comm)
415 B2 = 11.6;
416 C2 = 0.60;
417 D2 = 0.04;
418 alphaConst = 0.64;
419 }
420 bEnergy = Bj[fSelectedShell];
421 G4double v2 = 0.25*emax/(bEnergy*gam*gam);
422 v = std::sqrt(v2);
423 u = Ry/bEnergy;
424 wc = 4.*v2 - 2.*v - 0.25*u;
425
426 G4double L1 = (C1 * fGpow->powA(v, D1)) / (1. + E1 * fGpow->powA(v, (D1 + 4.)));
427 G4double L2 = C2 * fGpow->powA(v, D2);
428 G4double H1 = (A1 * G4Log(1. + v2)) / (v2 + (B1 / v2));
429 G4double H2 = (A2 / v2) + (B2 / (v2 * v2));
430
431 F1 = L1 + H1;
432 F2 = (L2 * H2) / (L2 + H2);
433 return emax;
434}
435
436//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
437
439G4DNARuddIonisationDynamicModel::SampleElectronEnergy()
440{
441 // sampling is performed for proton projectile
442 G4double emax = MaxEnergy();
443
444 ComputeIntegral(0.0, emax);
445 G4double e = SampleValue();
446 if (verbose > 1) {
447 G4cout << "G4DNARuddIonisationDynamicModel::SampleElectronEnergy: "
448 << fParticle->GetParticleName()
449 << " Escaled(keV)=" << fScaledEnergy/CLHEP::keV << " Ee(keV)=" << e/CLHEP::keV
450 << " Emax(keV)=" << emax/CLHEP::keV << " shell=" << fSelectedShell
451 << G4endl;
452 }
453 return e;
454}
455
456//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
457
459{
460 // Shells ids are 0 1 2 3 4 (4 is k shell)
461 // !!Attention, "energyTransfer" here is the energy transfered to the electron which means
462 // that the secondary kinetic energy is w = energyTransfer - bindingEnergy
463 //
464 // ds S F1(nu) + w * F2(nu)
465 // ---- = G(k) * ---- -------------------------------------------
466 // dw Bj (1+w)^3 * [1 + exp{alpha * (w - wc) / nu}]
467 //
468 // w is the secondary electron kinetic Energy in eV
469 //
470 // All the other parameters can be found in Rudd's Papers
471 //
472 // M.Eugene Rudd, 1988, User-Friendly model for the energy distribution of
473 // electrons from protons or electron collisions. Nucl. Tracks Rad. Meas.Vol 16 N0 2/3 pp 219-218
474 //
475 G4double w = e/bEnergy;
476 G4double x = alphaConst*(w - wc)/v;
477 G4double y = (x > -15.) ? 1.0 + G4Exp(x) : 1.0;
478
479 G4double res = CorrectionFactor() * (F1 + w*F2) /
480 (fGpow->powN((1. + w)/u, 3) * y);
481
482 if (isHelium) {
483 G4double energyTransfer = e + bEnergy;
484 G4double Zeff = 2.0 -
485 (sCoefficient[0] * S_1s(fScaledEnergy, energyTransfer, slaterEffectiveCharge[0], 1.) +
486 sCoefficient[1] * S_2s(fScaledEnergy, energyTransfer, slaterEffectiveCharge[1], 2.) +
487 sCoefficient[2] * S_2p(fScaledEnergy, energyTransfer, slaterEffectiveCharge[2], 2.) );
488
489 res *= Zeff * Zeff;
490 }
491 return std::max(res, 0.0);
492}
493
494//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
495
496G4double G4DNARuddIonisationDynamicModel::S_1s(G4double kine,
497 G4double energyTransfer,
498 G4double slaterEffCharge,
499 G4double shellNumber)
500{
501 // 1 - e^(-2r) * ( 1 + 2 r + 2 r^2)
502 // Dingfelder, in Chattanooga 2005 proceedings, formula (7)
503
504 G4double r = Rh(kine, energyTransfer, slaterEffCharge, shellNumber);
505 G4double value = 1. - G4Exp(-2 * r) * ( ( 2. * r + 2. ) * r + 1. );
506 return value;
507}
508
509//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
510
511G4double G4DNARuddIonisationDynamicModel::S_2s(G4double kine,
512 G4double energyTransfer,
513 G4double slaterEffCharge,
514 G4double shellNumber)
515{
516 // 1 - e^(-2 r) * ( 1 + 2 r + 2 r^2 + 2 r^4)
517 // Dingfelder, in Chattanooga 2005 proceedings, formula (8)
518
519 G4double r = Rh(kine, energyTransfer, slaterEffCharge, shellNumber);
520 G4double value =
521 1. - G4Exp(-2 * r) * (((2. * r * r + 2.) * r + 2.) * r + 1.);
522
523 return value;
524}
525
526//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
527
528G4double G4DNARuddIonisationDynamicModel::S_2p(G4double kine,
529 G4double energyTransfer,
530 G4double slaterEffCharge,
531 G4double shellNumber)
532{
533 // 1 - e^(-2 r) * ( 1 + 2 r + 2 r^2 + 4/3 r^3 + 2/3 r^4)
534 // Dingfelder, in Chattanooga 2005 proceedings, formula (9)
535
536 G4double r = Rh(kine, energyTransfer, slaterEffCharge, shellNumber);
537 G4double value =
538 1. - G4Exp(-2 * r) * (((( 2./3. * r + 4./3.) * r + 2.) * r + 2.) * r + 1.);
539
540 return value;
541}
542
543//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
544
545G4double G4DNARuddIonisationDynamicModel::Rh(G4double ekin, G4double etrans,
546 G4double q, G4double shell)
547{
548 // The following values are provided by M. Dingfelder (priv. comm)
549 // Dingfelder, in Chattanooga 2005 proceedings, p 4
550
551 G4double escaled = CLHEP::electron_mass_c2/fMass * ekin;
552 const G4double H = 13.60569172 * CLHEP::eV;
553 G4double value = 2.0*std::sqrt(escaled / H)*q*H /(etrans*shell);
554
555 return value;
556}
557
558//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
559
560G4double G4DNARuddIonisationDynamicModel::CorrectionFactor()
561{
562 // ZF Shortened
563 G4double res = 1.0;
564 if (fSelectedShell < 4) {
565 const G4double ln10 = fGpow->logZ(10);
566 G4double x = 2.0*((G4Log(fScaledEnergy/(fMassRate*CLHEP::eV))/ln10) - 4.2);
567 // The following values are provided by M. Dingfelder (priv. comm)
568 res = 0.6/(1.0 + G4Exp(x)) + 0.9;
569 }
570 return res;
571}
572
573//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ eIonizedMolecule
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:132
G4double G4Log(G4double x)
Definition G4Log.hh:169
CLHEP::Hep3Vector G4ThreeVector
@ fStopButAlive
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define C1
#define G4UniformRand()
Definition Randomize.hh:52
static G4DNAChemistryManager * 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()
G4DNARuddIonisationDynamicModel(const G4ParticleDefinition *p=nullptr, const G4String &nam="DNARuddIonisationDynamicModel")
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 Initialise(const G4ParticleDefinition *, const G4DataVector &) override
G4double ProbabilityDensityFunction(G4double ekin) override
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4EmParameters * Instance()
G4bool DNAStationary() const
const G4String & GetDirLEDATA() const
G4int SampleReactionChannel(const G4double energy, const G4double rand, std::size_t &lastidx) const
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
G4EmCorrections * EmCorrections()
static G4ExtendedPhysicsVector * BuildExtendedVector(const G4String &dirpath, const G4String &filename, const G4int N, const G4int length, const G4double unitE=1.0, const G4double unitS=1.0)
Definition G4MatUtils.cc:39
const G4Material * GetMaterial() const
std::size_t GetIndex() const
G4Material * FindMaterial(const G4String &name) const
static G4NistManager * Instance()
const G4String & GetParticleType() const
const G4String & GetParticleName() const
const G4String & GetParticleSubType() const
static G4Pow * GetInstance()
Definition G4Pow.cc:41
virtual G4ThreeVector & SampleDirectionForShell(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, G4int shellID, const G4Material *)
G4VEmAngularDistribution * GetAngularDistribution()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
void SetDeexcitationFlag(G4bool val)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67
void SetAngularDistribution(G4VEmAngularDistribution *)
G4double ComputeIntegral(const G4double emin, const G4double emax)
void InitialiseIntegrator(G4double accuracy, G4double fact1, G4double fact2, G4double de, G4double dmin, G4double dmax)
#define DBL_MAX
Definition templates.hh:62