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