Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4DNAMillerGreenExcitationModel.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// $Id$
27// GEANT4 tag $Name: $
28//
29
31#include "G4SystemOfUnits.hh"
34
35//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
36
37using namespace std;
38
39//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
40
42 const G4String& nam)
43 :G4VEmModel(nam),isInitialised(false)
44{
45 // nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
46 fpMolWaterDensity = 0;
47
48 nLevels=0;
49 kineticEnergyCorrection[0]=0.;
50 kineticEnergyCorrection[1]=0.;
51 kineticEnergyCorrection[2]=0.;
52 kineticEnergyCorrection[3]=0.;
53
54 verboseLevel= 0;
55 // Verbosity scale:
56 // 0 = nothing
57 // 1 = warning for energy non-conservation
58 // 2 = details of energy budget
59 // 3 = calculation of cross sections, file openings, sampling of atoms
60 // 4 = entering in methods
61
62 if( verboseLevel>0 )
63 {
64 G4cout << "Miller & Green excitation model is constructed " << G4endl;
65 }
67}
68
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
70
72{}
73
74//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
75
77 const G4DataVector& /*cuts*/)
78{
79
80 if (verboseLevel > 3)
81 G4cout << "Calling G4DNAMillerGreenExcitationModel::Initialise()" << G4endl;
82
83 // Energy limits
84
88 G4ParticleDefinition* hydrogenDef = instance->GetIon("hydrogen");
89 G4ParticleDefinition* alphaPlusPlusDef = instance->GetIon("alpha++");
90 G4ParticleDefinition* alphaPlusDef = instance->GetIon("alpha+");
91 G4ParticleDefinition* heliumDef = instance->GetIon("helium");
92
93 G4String proton;
94 G4String hydrogen;
95 G4String alphaPlusPlus;
96 G4String alphaPlus;
97 G4String helium;
98
99 // LIMITS AND CONSTANTS
100
101 proton = protonDef->GetParticleName();
102 lowEnergyLimit[proton] = 10. * eV;
103 highEnergyLimit[proton] = 500. * keV;
104
105 kineticEnergyCorrection[0] = 1.;
106 slaterEffectiveCharge[0][0] = 0.;
107 slaterEffectiveCharge[1][0] = 0.;
108 slaterEffectiveCharge[2][0] = 0.;
109 sCoefficient[0][0] = 0.;
110 sCoefficient[1][0] = 0.;
111 sCoefficient[2][0] = 0.;
112
113 hydrogen = hydrogenDef->GetParticleName();
114 lowEnergyLimit[hydrogen] = 10. * eV;
115 highEnergyLimit[hydrogen] = 500. * keV;
116
117 kineticEnergyCorrection[0] = 1.;
118 slaterEffectiveCharge[0][0] = 0.;
119 slaterEffectiveCharge[1][0] = 0.;
120 slaterEffectiveCharge[2][0] = 0.;
121 sCoefficient[0][0] = 0.;
122 sCoefficient[1][0] = 0.;
123 sCoefficient[2][0] = 0.;
124
125 alphaPlusPlus = alphaPlusPlusDef->GetParticleName();
126 lowEnergyLimit[alphaPlusPlus] = 1. * keV;
127 highEnergyLimit[alphaPlusPlus] = 400. * MeV;
128
129 kineticEnergyCorrection[1] = 0.9382723/3.727417;
130 slaterEffectiveCharge[0][1]=0.;
131 slaterEffectiveCharge[1][1]=0.;
132 slaterEffectiveCharge[2][1]=0.;
133 sCoefficient[0][1]=0.;
134 sCoefficient[1][1]=0.;
135 sCoefficient[2][1]=0.;
136
137 alphaPlus = alphaPlusDef->GetParticleName();
138 lowEnergyLimit[alphaPlus] = 1. * keV;
139 highEnergyLimit[alphaPlus] = 400. * MeV;
140
141 kineticEnergyCorrection[2] = 0.9382723/3.727417;
142 slaterEffectiveCharge[0][2]=2.0;
143
144 // Following values provided by M. Dingfelder
145 slaterEffectiveCharge[1][2]=2.00;
146 slaterEffectiveCharge[2][2]=2.00;
147 //
148 sCoefficient[0][2]=0.7;
149 sCoefficient[1][2]=0.15;
150 sCoefficient[2][2]=0.15;
151
152 helium = heliumDef->GetParticleName();
153 lowEnergyLimit[helium] = 1. * keV;
154 highEnergyLimit[helium] = 400. * MeV;
155
156 kineticEnergyCorrection[3] = 0.9382723/3.727417;
157 slaterEffectiveCharge[0][3]=1.7;
158 slaterEffectiveCharge[1][3]=1.15;
159 slaterEffectiveCharge[2][3]=1.15;
160 sCoefficient[0][3]=0.5;
161 sCoefficient[1][3]=0.25;
162 sCoefficient[2][3]=0.25;
163
164 //
165
166 if (particle==protonDef)
167 {
168 SetLowEnergyLimit(lowEnergyLimit[proton]);
169 SetHighEnergyLimit(highEnergyLimit[proton]);
170 }
171
172 if (particle==hydrogenDef)
173 {
174 SetLowEnergyLimit(lowEnergyLimit[hydrogen]);
175 SetHighEnergyLimit(highEnergyLimit[hydrogen]);
176 }
177
178 if (particle==alphaPlusPlusDef)
179 {
180 SetLowEnergyLimit(lowEnergyLimit[alphaPlusPlus]);
181 SetHighEnergyLimit(highEnergyLimit[alphaPlusPlus]);
182 }
183
184 if (particle==alphaPlusDef)
185 {
186 SetLowEnergyLimit(lowEnergyLimit[alphaPlus]);
187 SetHighEnergyLimit(highEnergyLimit[alphaPlus]);
188 }
189
190 if (particle==heliumDef)
191 {
192 SetLowEnergyLimit(lowEnergyLimit[helium]);
193 SetHighEnergyLimit(highEnergyLimit[helium]);
194 }
195
196 //
197
198 nLevels = waterExcitation.NumberOfLevels();
199
200 //
201 if( verboseLevel>0 )
202 {
203 G4cout << "Miller & Green excitation model is initialized " << G4endl
204 << "Energy range: "
205 << LowEnergyLimit() / eV << " eV - "
206 << HighEnergyLimit() / keV << " keV for "
207 << particle->GetParticleName()
208 << G4endl;
209 }
210
211 // Initialize water density pointer
213
214 if (isInitialised) { return; }
216 isInitialised = true;
217
218}
219
220//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
221
223 const G4ParticleDefinition* particleDefinition,
224 G4double k,
225 G4double,
226 G4double)
227{
228 if (verboseLevel > 3)
229 G4cout << "Calling CrossSectionPerVolume() of G4DNAMillerGreenExcitationModel" << G4endl;
230
231 // Calculate total cross section for model
232
233 G4DNAGenericIonsManager *instance;
235
236 if (
237 particleDefinition != G4Proton::ProtonDefinition()
238 &&
239 particleDefinition != instance->GetIon("hydrogen")
240 &&
241 particleDefinition != instance->GetIon("alpha++")
242 &&
243 particleDefinition != instance->GetIon("alpha+")
244 &&
245 particleDefinition != instance->GetIon("helium")
246 )
247
248 return 0;
249
250 G4double lowLim = 0;
251 G4double highLim = 0;
252 G4double crossSection = 0.;
253
254 G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
255
256 if(waterDensity!= 0.0)
257 // if (material == nistwater || material->GetBaseMaterial() == nistwater)
258 {
259// G4cout << "WATER DENSITY = " << waterDensity*G4Material::GetMaterial("G4_WATER")->GetMassOfMolecule()/(g/cm3)
260// << G4endl;
261 const G4String& particleName = particleDefinition->GetParticleName();
262
263 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
264 pos1 = lowEnergyLimit.find(particleName);
265
266 if (pos1 != lowEnergyLimit.end())
267 {
268 lowLim = pos1->second;
269 }
270
271 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
272 pos2 = highEnergyLimit.find(particleName);
273
274 if (pos2 != highEnergyLimit.end())
275 {
276 highLim = pos2->second;
277 }
278
279 if (k >= lowLim && k < highLim)
280 {
281 crossSection = Sum(k,particleDefinition);
282
283 // add ONE or TWO electron-water excitation for alpha+ and helium
284 /*
285 if ( particleDefinition == instance->GetIon("alpha+")
286 ||
287 particleDefinition == instance->GetIon("helium")
288 )
289 {
290
291 G4DNAEmfietzoglouExcitationModel * excitationXS = new G4DNAEmfietzoglouExcitationModel();
292 excitationXS->Initialise(G4Electron::ElectronDefinition());
293
294 G4double sigmaExcitation=0;
295 G4double tmp =0.;
296
297 if (k*0.511/3728 > 8.23*eV && k*0.511/3728 < 10*MeV ) sigmaExcitation =
298 excitationXS->CrossSectionPerVolume(material,G4Electron::ElectronDefinition(),k*0.511/3728,tmp,tmp)
299 /material->GetAtomicNumDensityVector()[1];
300
301 if ( particleDefinition == instance->GetIon("alpha+") )
302 crossSection = crossSection + sigmaExcitation ;
303
304 if ( particleDefinition == instance->GetIon("helium") )
305 crossSection = crossSection + 2*sigmaExcitation ;
306
307 delete excitationXS;
308
309 // Alternative excitation model
310
311 G4DNABornExcitationModel * excitationXS = new G4DNABornExcitationModel();
312 excitationXS->Initialise(G4Electron::ElectronDefinition());
313
314 G4double sigmaExcitation=0;
315 G4double tmp=0;
316
317 if (k*0.511/3728 > 9*eV && k*0.511/3728 < 1*MeV ) sigmaExcitation =
318 excitationXS->CrossSectionPerVolume(material,G4Electron::ElectronDefinition(),k*0.511/3728,tmp,tmp)
319 /material->GetAtomicNumDensityVector()[1];
320
321 if ( particleDefinition == instance->GetIon("alpha+") )
322 crossSection = crossSection + sigmaExcitation ;
323
324 if ( particleDefinition == instance->GetIon("helium") )
325 crossSection = crossSection + 2*sigmaExcitation ;
326
327 delete excitationXS;
328
329 }
330*/
331
332 }
333
334 if (verboseLevel > 2)
335 {
336 G4cout << "__________________________________" << G4endl;
337 G4cout << "°°° G4DNAMillerGreenExcitationModel - XS INFO START" << G4endl;
338 G4cout << "°°° Kinetic energy(eV)=" << k/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
339 G4cout << "°°° Cross section per water molecule (cm^2)=" << crossSection/cm/cm << G4endl;
340 G4cout << "°°° Cross section per water molecule (cm^-1)=" << crossSection*waterDensity/(1./cm) << G4endl;
341 // G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
342 G4cout << "°°° G4DNAMillerGreenExcitationModel - XS INFO END" << G4endl;
343 }
344 }
345 else
346 {
347 if (verboseLevel > 2)
348 {
349 G4cout << "MillerGreenExcitationModel : WARNING Water density is NULL" << G4endl;
350 }
351 }
352
353 return crossSection*waterDensity;
354 // return crossSection*material->GetAtomicNumDensityVector()[1];
355
356}
357
358//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
359
360void G4DNAMillerGreenExcitationModel::SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/,
361 const G4MaterialCutsCouple* /*couple*/,
362 const G4DynamicParticle* aDynamicParticle,
363 G4double,
364 G4double)
365{
366
367 if (verboseLevel > 3)
368 G4cout << "Calling SampleSecondaries() of G4DNAMillerGreenExcitationModel" << G4endl;
369
370 G4double particleEnergy0 = aDynamicParticle->GetKineticEnergy();
371
372 G4int level = RandomSelect(particleEnergy0,aDynamicParticle->GetDefinition());
373
374 // G4double excitationEnergy = waterExcitation.ExcitationEnergy(level);
375
376 // Dingfelder's excitation levels
377 const G4double excitation[]={ 8.17*eV, 10.13*eV, 11.31*eV, 12.91*eV, 14.50*eV};
378 G4double excitationEnergy = excitation[level];
379
380 G4double newEnergy = particleEnergy0 - excitationEnergy;
381
382 if (newEnergy>0)
383 {
387
388 const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
390 level,
391 theIncomingTrack);
392
393 }
394
395}
396
397//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
398
399G4double G4DNAMillerGreenExcitationModel::PartialCrossSection(G4double k, G4int excitationLevel,
400 const G4ParticleDefinition* particleDefinition)
401{
402 // ( ( z * aj ) ^ omegaj ) * ( t - ej ) ^ nu
403 // sigma(t) = zEff^2 * sigma0 * --------------------------------------------
404 // jj ^ ( omegaj + nu ) + t ^ ( omegaj + nu )
405 //
406 // where t is the kinetic energy corrected by Helium mass over proton mass for Helium ions
407 //
408 // zEff is:
409 // 1 for protons
410 // 2 for alpha++
411 // and 2 - c1 S_1s - c2 S_2s - c3 S_2p for alpha+ and He
412 //
413 // Dingfelder et al., RPC 59, 255-275, 2000 from Miller and Green (1973)
414 // Formula (34) and Table 2
415
416 const G4double sigma0(1.E+8 * barn);
417 const G4double nu(1.);
418 const G4double aj[]={876.*eV, 2084.* eV, 1373.*eV, 692.*eV, 900.*eV};
419 const G4double jj[]={19820.*eV, 23490.*eV, 27770.*eV, 30830.*eV, 33080.*eV};
420 const G4double omegaj[]={0.85, 0.88, 0.88, 0.78, 0.78};
421
422 // Dingfelder's excitation levels
423 const G4double Eliq[5]={ 8.17*eV, 10.13*eV, 11.31*eV, 12.91*eV, 14.50*eV};
424
425 G4int particleTypeIndex = 0;
426 G4DNAGenericIonsManager* instance;
428
429 if (particleDefinition == G4Proton::ProtonDefinition()) particleTypeIndex=0;
430 if (particleDefinition == instance->GetIon("hydrogen")) particleTypeIndex=0;
431 if (particleDefinition == instance->GetIon("alpha++")) particleTypeIndex=1;
432 if (particleDefinition == instance->GetIon("alpha+")) particleTypeIndex=2;
433 if (particleDefinition == instance->GetIon("helium")) particleTypeIndex=3;
434
435 G4double tCorrected;
436 tCorrected = k * kineticEnergyCorrection[particleTypeIndex];
437
438 // SI - added protection
439 if (tCorrected < Eliq[excitationLevel]) return 0;
440 //
441
442 G4int z = 10;
443
444 G4double numerator;
445 numerator = std::pow(z * aj[excitationLevel], omegaj[excitationLevel]) *
446 std::pow(tCorrected - Eliq[excitationLevel], nu);
447
448 // H case : see S. Uehara et al. IJRB 77, 2, 139-154 (2001) - section 3.3
449
450 if (particleDefinition == instance->GetIon("hydrogen"))
451 numerator = std::pow(z * 0.75*aj[excitationLevel], omegaj[excitationLevel]) *
452 std::pow(tCorrected - Eliq[excitationLevel], nu);
453
454
455 G4double power;
456 power = omegaj[excitationLevel] + nu;
457
458 G4double denominator;
459 denominator = std::pow(jj[excitationLevel], power) + std::pow(tCorrected, power);
460
461 G4double zEff = particleDefinition->GetPDGCharge() / eplus + particleDefinition->GetLeptonNumber();
462
463 zEff -= ( sCoefficient[0][particleTypeIndex] * S_1s(k, Eliq[excitationLevel], slaterEffectiveCharge[0][particleTypeIndex], 1.) +
464 sCoefficient[1][particleTypeIndex] * S_2s(k, Eliq[excitationLevel], slaterEffectiveCharge[1][particleTypeIndex], 2.) +
465 sCoefficient[2][particleTypeIndex] * S_2p(k, Eliq[excitationLevel], slaterEffectiveCharge[2][particleTypeIndex], 2.) );
466
467 if (particleDefinition == instance->GetIon("hydrogen")) zEff = 1.;
468
469 G4double cross = sigma0 * zEff * zEff * numerator / denominator;
470
471
472 return cross;
473}
474
475//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
476
477G4int G4DNAMillerGreenExcitationModel::RandomSelect(G4double k,const G4ParticleDefinition* particle)
478{
479 G4int i = nLevels;
480 G4double value = 0.;
481 std::deque<double> values;
482
483 G4DNAGenericIonsManager *instance;
485
486 if ( particle == instance->GetIon("alpha++") ||
487 particle == G4Proton::ProtonDefinition()||
488 particle == instance->GetIon("hydrogen") ||
489 particle == instance->GetIon("alpha+") ||
490 particle == instance->GetIon("helium")
491 )
492 {
493 while (i > 0)
494 {
495 i--;
496 G4double partial = PartialCrossSection(k,i,particle);
497 values.push_front(partial);
498 value += partial;
499 }
500
501 value *= G4UniformRand();
502
503 i = nLevels;
504
505 while (i > 0)
506 {
507 i--;
508 if (values[i] > value) return i;
509 value -= values[i];
510 }
511 }
512
513 /*
514 // add ONE or TWO electron-water excitation for alpha+ and helium
515
516 if ( particle == instance->GetIon("alpha+")
517 ||
518 particle == instance->GetIon("helium")
519 )
520 {
521 while (i>0)
522 {
523 i--;
524
525 G4DNAEmfietzoglouExcitationModel * excitationXS = new G4DNAEmfietzoglouExcitationModel();
526 excitationXS->Initialise(G4Electron::ElectronDefinition());
527
528 G4double sigmaExcitation=0;
529
530 if (k*0.511/3728 > 8.23*eV && k*0.511/3728 < 10*MeV ) sigmaExcitation = excitationXS->PartialCrossSection(k*0.511/3728,i);
531
532 G4double partial = PartialCrossSection(k,i,particle);
533
534 if (particle == instance->GetIon("alpha+")) partial = PartialCrossSection(k,i,particle) + sigmaExcitation;
535 if (particle == instance->GetIon("helium")) partial = PartialCrossSection(k,i,particle) + 2*sigmaExcitation;
536
537 values.push_front(partial);
538 value += partial;
539 delete excitationXS;
540 }
541
542 value*=G4UniformRand();
543
544 i=5;
545 while (i>0)
546 {
547 i--;
548
549 if (values[i]>value) return i;
550
551 value-=values[i];
552 }
553 }
554*/
555
556 return 0;
557}
558
559//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
560
561G4double G4DNAMillerGreenExcitationModel::Sum(G4double k, const G4ParticleDefinition* particle)
562{
563 G4double totalCrossSection = 0.;
564
565 for (G4int i=0; i<nLevels; i++)
566 {
567 totalCrossSection += PartialCrossSection(k,i,particle);
568 }
569 return totalCrossSection;
570}
571
572//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
573
574G4double G4DNAMillerGreenExcitationModel::S_1s(G4double t,
575 G4double energyTransferred,
576 G4double _slaterEffectiveCharge,
577 G4double shellNumber)
578{
579 // 1 - e^(-2r) * ( 1 + 2 r + 2 r^2)
580 // Dingfelder, in Chattanooga 2005 proceedings, formula (7)
581
582 G4double r = R(t, energyTransferred, _slaterEffectiveCharge, shellNumber);
583 G4double value = 1. - std::exp(-2 * r) * ( ( 2. * r + 2. ) * r + 1. );
584
585 return value;
586}
587
588
589//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
590
591G4double G4DNAMillerGreenExcitationModel::S_2s(G4double t,
592 G4double energyTransferred,
593 G4double _slaterEffectiveCharge,
594 G4double shellNumber)
595{
596 // 1 - e^(-2 r) * ( 1 + 2 r + 2 r^2 + 2 r^4)
597 // Dingfelder, in Chattanooga 2005 proceedings, formula (8)
598
599 G4double r = R(t, energyTransferred, _slaterEffectiveCharge, shellNumber);
600 G4double value = 1. - std::exp(-2 * r) * (((2. * r * r + 2.) * r + 2.) * r + 1.);
601
602 return value;
603
604}
605
606//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
607
608G4double G4DNAMillerGreenExcitationModel::S_2p(G4double t,
609 G4double energyTransferred,
610 G4double _slaterEffectiveCharge,
611 G4double shellNumber)
612{
613 // 1 - e^(-2 r) * ( 1 + 2 r + 2 r^2 + 4/3 r^3 + 2/3 r^4)
614 // Dingfelder, in Chattanooga 2005 proceedings, formula (9)
615
616 G4double r = R(t, energyTransferred, _slaterEffectiveCharge, shellNumber);
617 G4double value = 1. - std::exp(-2 * r) * (((( 2./3. * r + 4./3.) * r + 2.) * r + 2.) * r + 1.);
618
619 return value;
620}
621
622//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
623
624G4double G4DNAMillerGreenExcitationModel::R(G4double t,
625 G4double energyTransferred,
626 G4double _slaterEffectiveCharge,
627 G4double shellNumber)
628{
629 // tElectron = m_electron / m_alpha * t
630 // Dingfelder, in Chattanooga 2005 proceedings, p 4
631
632 G4double tElectron = 0.511/3728. * t;
633
634 // The following is provided by M. Dingfelder
635 G4double H = 2.*13.60569172 * eV;
636 G4double value = std::sqrt ( 2. * tElectron / H ) / ( energyTransferred / H ) * (_slaterEffectiveCharge/shellNumber);
637
638 return value;
639}
640
@ eExcitedMolecule
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
static G4DNAGenericIonsManager * Instance(void)
G4ParticleDefinition * GetIon(const G4String &name)
G4DNAMillerGreenExcitationModel(const G4ParticleDefinition *p=0, const G4String &nam="DNAMillerGreenExcitationModel")
G4ParticleChangeForGamma * fParticleChangeForGamma
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
static G4DNAMolecularMaterial * Instance()
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
size_t GetIndex() const
Definition: G4Material.hh:261
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:576
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4double GetPDGCharge() const
const G4String & GetParticleName() const
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592
void ProposeLocalEnergyDeposit(G4double anEnergyPart)