Geant4  10.02.p02
G4DNARuddIonisationModel.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: G4DNARuddIonisationModel.cc 95457 2016-02-11 10:18:10Z gcosmo $
27 // GEANT4 tag $Name: $
28 //
29 
31 #include "G4PhysicalConstants.hh"
32 #include "G4SystemOfUnits.hh"
33 #include "G4UAtomicDeexcitation.hh"
34 #include "G4LossTableManager.hh"
35 #include "G4DNAChemistryManager.hh"
37 #include "G4DNARuddAngle.hh"
38 #include "G4DeltaAngle.hh"
39 
40 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
41 
42 using namespace std;
43 
44 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
45 
47  const G4String& nam) :
48  G4VEmModel(nam), isInitialised(false)
49 {
50  fpWaterDensity = 0;
51 
52  slaterEffectiveCharge[0] = 0.;
53  slaterEffectiveCharge[1] = 0.;
54  slaterEffectiveCharge[2] = 0.;
55  sCoefficient[0] = 0.;
56  sCoefficient[1] = 0.;
57  sCoefficient[2] = 0.;
58 
59  lowEnergyLimitForZ1 = 0 * eV;
60  lowEnergyLimitForZ2 = 0 * eV;
65 
66  verboseLevel = 0;
67  // Verbosity scale:
68  // 0 = nothing
69  // 1 = warning for energy non-conservation
70  // 2 = details of energy budget
71  // 3 = calculation of cross sections, file openings, sampling of atoms
72  // 4 = entering in methods
73 
74  if (verboseLevel > 0)
75  {
76  G4cout << "Rudd ionisation model is constructed " << G4endl;
77  }
78 
79  // define default angular generator
81 
82  //Mark this model as "applicable" for atomic deexcitation
83  SetDeexcitationFlag(true);
86 
87  // Selection of stationary mode
88 
89  statCode = false;
90 }
91 
92 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
93 
95 {
96  // Cross section
97 
98  std::map<G4String, G4DNACrossSectionDataSet*, std::less<G4String> >::iterator pos;
99  for (pos = tableData.begin(); pos != tableData.end(); ++pos)
100  {
101  G4DNACrossSectionDataSet* table = pos->second;
102  delete table;
103  }
104 
105  // The following removal is forbidden since G4VEnergyLossmodel takes care of deletion
106  // Coverity however will signal this as an error
107  //if (fAtomDeexcitation) {delete fAtomDeexcitation;}
108 
109 }
110 
111 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
112 
114  const G4DataVector& /*cuts*/)
115 {
116 
117  if (verboseLevel > 3)
118  {
119  G4cout << "Calling G4DNARuddIonisationModel::Initialise()" << G4endl;
120  }
121 
122  // Energy limits
123 
124  G4String fileProton("dna/sigma_ionisation_p_rudd");
125  G4String fileHydrogen("dna/sigma_ionisation_h_rudd");
126  G4String fileAlphaPlusPlus("dna/sigma_ionisation_alphaplusplus_rudd");
127  G4String fileAlphaPlus("dna/sigma_ionisation_alphaplus_rudd");
128  G4String fileHelium("dna/sigma_ionisation_he_rudd");
129 
133  G4ParticleDefinition* hydrogenDef = instance->GetIon("hydrogen");
134  G4ParticleDefinition* alphaPlusPlusDef = instance->GetIon("alpha++");
135  G4ParticleDefinition* alphaPlusDef = instance->GetIon("alpha+");
136  G4ParticleDefinition* heliumDef = instance->GetIon("helium");
137 
139  G4String hydrogen;
140  G4String alphaPlusPlus;
141  G4String alphaPlus;
142  G4String helium;
143 
144  G4double scaleFactor = 1 * m*m;
145 
146  // LIMITS AND DATA
147 
148  // ********************************************************
149 
150  proton = protonDef->GetParticleName();
151  tableFile[proton] = fileProton;
152 
154  highEnergyLimit[proton] = 500. * keV;
155 
156  // Cross section
157 
159  eV,
160  scaleFactor );
161  tableProton->LoadData(fileProton);
162  tableData[proton] = tableProton;
163 
164  // ********************************************************
165 
166  hydrogen = hydrogenDef->GetParticleName();
167  tableFile[hydrogen] = fileHydrogen;
168 
170  highEnergyLimit[hydrogen] = 100. * MeV;
171 
172  // Cross section
173 
175  eV,
176  scaleFactor );
177  tableHydrogen->LoadData(fileHydrogen);
178 
179  tableData[hydrogen] = tableHydrogen;
180 
181  // ********************************************************
182 
183  alphaPlusPlus = alphaPlusPlusDef->GetParticleName();
184  tableFile[alphaPlusPlus] = fileAlphaPlusPlus;
185 
186  lowEnergyLimit[alphaPlusPlus] = lowEnergyLimitForZ2;
187  highEnergyLimit[alphaPlusPlus] = 400. * MeV;
188 
189  // Cross section
190 
192  eV,
193  scaleFactor );
194  tableAlphaPlusPlus->LoadData(fileAlphaPlusPlus);
195 
196  tableData[alphaPlusPlus] = tableAlphaPlusPlus;
197 
198  // ********************************************************
199 
200  alphaPlus = alphaPlusDef->GetParticleName();
201  tableFile[alphaPlus] = fileAlphaPlus;
202 
203  lowEnergyLimit[alphaPlus] = lowEnergyLimitForZ2;
204  highEnergyLimit[alphaPlus] = 400. * MeV;
205 
206  // Cross section
207 
209  eV,
210  scaleFactor );
211  tableAlphaPlus->LoadData(fileAlphaPlus);
212  tableData[alphaPlus] = tableAlphaPlus;
213 
214  // ********************************************************
215 
216  helium = heliumDef->GetParticleName();
217  tableFile[helium] = fileHelium;
218 
220  highEnergyLimit[helium] = 400. * MeV;
221 
222  // Cross section
223 
225  eV,
226  scaleFactor );
227  tableHelium->LoadData(fileHelium);
228  tableData[helium] = tableHelium;
229 
230  //
231 
232  if (particle==protonDef)
233  {
236  }
237 
238  if (particle==hydrogenDef)
239  {
242  }
243 
244  if (particle==heliumDef)
245  {
248  }
249 
250  if (particle==alphaPlusDef)
251  {
252  SetLowEnergyLimit(lowEnergyLimit[alphaPlus]);
254  }
255 
256  if (particle==alphaPlusPlusDef)
257  {
258  SetLowEnergyLimit(lowEnergyLimit[alphaPlusPlus]);
259  SetHighEnergyLimit(highEnergyLimit[alphaPlusPlus]);
260  }
261 
262  if( verboseLevel>0 )
263  {
264  G4cout << "Rudd ionisation model is initialized " << G4endl
265  << "Energy range: "
266  << LowEnergyLimit() / eV << " eV - "
267  << HighEnergyLimit() / keV << " keV for "
268  << particle->GetParticleName()
269  << G4endl;
270  }
271 
272  // Initialize water density pointer
274 
275  //
276 
278 
279  if (isInitialised)
280  { return;}
282  isInitialised = true;
283 
284 }
285 
286 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
287 
289  const G4ParticleDefinition* particleDefinition,
290  G4double k,
291  G4double,
292  G4double)
293 {
294  if (verboseLevel > 3)
295  {
296  G4cout << "Calling CrossSectionPerVolume() of G4DNARuddIonisationModel"
297  << G4endl;
298  }
299 
300  // Calculate total cross section for model
301 
304 
305  if (
306  particleDefinition != G4Proton::ProtonDefinition()
307  &&
308  particleDefinition != instance->GetIon("hydrogen")
309  &&
310  particleDefinition != instance->GetIon("alpha++")
311  &&
312  particleDefinition != instance->GetIon("alpha+")
313  &&
314  particleDefinition != instance->GetIon("helium")
315  )
316 
317  return 0;
318 
319  G4double lowLim = 0;
320 
321  if ( particleDefinition == G4Proton::ProtonDefinition()
322  || particleDefinition == instance->GetIon("hydrogen")
323  )
324 
326 
327  if ( particleDefinition == instance->GetIon("alpha++")
328  || particleDefinition == instance->GetIon("alpha+")
329  || particleDefinition == instance->GetIon("helium")
330  )
331 
333 
334  G4double highLim = 0;
335  G4double sigma=0;
336 
337  G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
338 
339  if(waterDensity!= 0.0)
340  // if (material == nistwater || material->GetBaseMaterial() == nistwater)
341  {
342  const G4String& particleName = particleDefinition->GetParticleName();
343 
344  // SI - the following is useless since lowLim is already defined
345  /*
346  std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
347  pos1 = lowEnergyLimit.find(particleName);
348 
349  if (pos1 != lowEnergyLimit.end())
350  {
351  lowLim = pos1->second;
352  }
353  */
354 
355  std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
356  pos2 = highEnergyLimit.find(particleName);
357 
358  if (pos2 != highEnergyLimit.end())
359  {
360  highLim = pos2->second;
361  }
362 
363  if (k <= highLim)
364  {
365  //SI : XS must not be zero otherwise sampling of secondaries method ignored
366 
367  if (k < lowLim) k = lowLim;
368 
369  //
370 
371  std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
372  pos = tableData.find(particleName);
373 
374  if (pos != tableData.end())
375  {
376  G4DNACrossSectionDataSet* table = pos->second;
377  if (table != 0)
378  {
379  sigma = table->FindValue(k);
380  }
381  }
382  else
383  {
384  G4Exception("G4DNARuddIonisationModel::CrossSectionPerVolume","em0002",
385  FatalException,"Model not applicable to particle type.");
386  }
387 
388  } // if (k >= lowLim && k < highLim)
389 
390  if (verboseLevel > 2)
391  {
392  G4cout << "__________________________________" << G4endl;
393  G4cout << "G4DNARuddIonisationModel - XS INFO START" << G4endl;
394  G4cout << "Kinetic energy(eV)=" << k/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
395  G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
396  G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
397  // G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
398  G4cout << "G4DNARuddIonisationModel - XS INFO END" << G4endl;
399  }
400 
401  } // if (waterMaterial)
402  else
403  {
404  if (verboseLevel > 2)
405  {
406  G4cout << "Warning : RuddIonisationModel: WATER DENSITY IS NULL" << G4endl;
407  }
408  }
409 
410  return sigma*waterDensity;
411  // return sigma*material->GetAtomicNumDensityVector()[1];
412 
413 }
414 
415 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
416 
417 void G4DNARuddIonisationModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
418  const G4MaterialCutsCouple* couple,
419  const G4DynamicParticle* particle,
420  G4double,
421  G4double)
422 {
423  if (verboseLevel > 3)
424  {
425  G4cout << "Calling SampleSecondaries() of G4DNARuddIonisationModel"
426  << G4endl;
427  }
428 
429  G4double lowLim = 0;
430  G4double highLim = 0;
431 
434 
435  if ( particle->GetDefinition() == G4Proton::ProtonDefinition()
436  || particle->GetDefinition() == instance->GetIon("hydrogen")
437  )
438 
439  lowLim = killBelowEnergyForZ1;
440 
441  if ( particle->GetDefinition() == instance->GetIon("alpha++")
442  || particle->GetDefinition() == instance->GetIon("alpha+")
443  || particle->GetDefinition() == instance->GetIon("helium")
444  )
445 
446  lowLim = killBelowEnergyForZ2;
447 
448  G4double k = particle->GetKineticEnergy();
449 
450  const G4String& particleName = particle->GetDefinition()->GetParticleName();
451 
452  // SI - the following is useless since lowLim is already defined
453  /*
454  std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
455  pos1 = lowEnergyLimit.find(particleName);
456 
457  if (pos1 != lowEnergyLimit.end())
458  {
459  lowLim = pos1->second;
460  }
461  */
462 
463  std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
464  pos2 = highEnergyLimit.find(particleName);
465 
466  if (pos2 != highEnergyLimit.end())
467  {
468  highLim = pos2->second;
469  }
470 
471  if (k >= lowLim && k <= highLim)
472  {
473  G4ParticleDefinition* definition = particle->GetDefinition();
474  G4ParticleMomentum primaryDirection = particle->GetMomentumDirection();
475  /*
476  G4double particleMass = definition->GetPDGMass();
477  G4double totalEnergy = k + particleMass;
478  G4double pSquare = k*(totalEnergy+particleMass);
479  G4double totalMomentum = std::sqrt(pSquare);
480  */
481 
482  G4int ionizationShell = RandomSelect(k,particleName);
483 
484  // sample deexcitation
485  // here we assume that H_{2}O electronic levels are the same of Oxigen.
486  // this can be considered true with a rough 10% error in energy on K-shell,
487 
488  G4int secNumberInit = 0;// need to know at a certain point the enrgy of secondaries
489  G4int secNumberFinal = 0;// So I'll make the diference and then sum the energies
490 
492  bindingEnergy = waterStructure.IonisationEnergy(ionizationShell);
493 
494  //SI: additional protection if tcs interpolation method is modified
495  if (k<bindingEnergy) return;
496  //
497 
498  G4int Z = 8;
500  {
502 
503  if (ionizationShell <5 && ionizationShell >1)
504  {
505  as = G4AtomicShellEnumerator(4-ionizationShell);
506  }
507  else if (ionizationShell <2)
508  {
509  as = G4AtomicShellEnumerator(3);
510  }
511 
512  // DEBUG
513  // if (ionizationShell == 4) {
514  //
515  // G4cout << "Z: " << Z << " as: " << as
516  // << " ionizationShell: " << ionizationShell << " bindingEnergy: "<< bindingEnergy/eV << G4endl;
517  // G4cout << "Press <Enter> key to continue..." << G4endl;
518  // G4cin.ignore();
519  // }
520 
521  const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
522  secNumberInit = fvect->size();
523  fAtomDeexcitation->GenerateParticles(fvect, shell, Z, 0, 0);
524  secNumberFinal = fvect->size();
525  }
526 
527  G4double secondaryKinetic = RandomizeEjectedElectronEnergy(definition,k,ionizationShell);
528 
529  G4ThreeVector deltaDirection =
530  GetAngularDistribution()->SampleDirectionForShell(particle, secondaryKinetic,
531  Z, ionizationShell,
532  couple->GetMaterial());
533 
534  // Ignored for ions on electrons
535  /*
536  G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
537 
538  G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
539  G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
540  G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
541  G4double finalMomentum = std::sqrt(finalPx*finalPx+finalPy*finalPy+finalPz*finalPz);
542  finalPx /= finalMomentum;
543  finalPy /= finalMomentum;
544  finalPz /= finalMomentum;
545 
546  G4ThreeVector direction;
547  direction.set(finalPx,finalPy,finalPz);
548 
549  fParticleChangeForGamma->ProposeMomentumDirection(direction.unit()) ;
550  */
552 
553  G4double scatteredEnergy = k-bindingEnergy-secondaryKinetic;
554  G4double deexSecEnergy = 0;
555  for (G4int j=secNumberInit; j < secNumberFinal; j++)
556  {
557 
558  deexSecEnergy = deexSecEnergy + (*fvect)[j]->GetKineticEnergy();
559 
560  }
561 
562  if (!statCode)
563  {
565  fParticleChangeForGamma->ProposeLocalEnergyDeposit(k-scatteredEnergy-secondaryKinetic-deexSecEnergy);
566  }
567  else
568  {
571  }
572 
573  // debug
574  // k-scatteredEnergy-secondaryKinetic-deexSecEnergy = k-(k-bindingEnergy-secondaryKinetic)-secondaryKinetic-deexSecEnergy =
575  // = k-k+bindingEnergy+secondaryKinetic-secondaryKinetic-deexSecEnergy=
576  // = bindingEnergy-deexSecEnergy
577  // SO deexSecEnergy=0 => LocalEnergyDeposit = bindingEnergy
578 
579  G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic);
580  fvect->push_back(dp);
581 
582  const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
584  ionizationShell,
585  theIncomingTrack);
586  }
587 
588  // SI - not useful since low energy of model is 0 eV
589 
590  if (k < lowLim)
591  {
595  }
596 
597 }
598 
599 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
600 
602  G4double k,
603  G4int shell)
604 {
605  G4double maximumKineticEnergyTransfer = 0.;
606 
609 
610  if (particleDefinition == G4Proton::ProtonDefinition()
611  || particleDefinition == instance->GetIon("hydrogen"))
612  {
613  maximumKineticEnergyTransfer = 4. * (electron_mass_c2 / proton_mass_c2) * k;
614  }
615 
616  else if (particleDefinition == instance->GetIon("helium")
617  || particleDefinition == instance->GetIon("alpha+")
618  || particleDefinition == instance->GetIon("alpha++"))
619  {
620  maximumKineticEnergyTransfer = 4. * (0.511 / 3728) * k;
621  }
622 
623  G4double crossSectionMaximum = 0.;
624 
625  for (G4double value = waterStructure.IonisationEnergy(shell);
626  value <= 5. * waterStructure.IonisationEnergy(shell) && k >= value;
627  value += 0.1 * eV)
628  {
629  G4double differentialCrossSection =
630  DifferentialCrossSection(particleDefinition, k, value, shell);
631  if (differentialCrossSection >= crossSectionMaximum)
632  crossSectionMaximum = differentialCrossSection;
633  }
634 
635  G4double secElecKinetic = 0.;
636 
637  do
638  {
639  secElecKinetic = G4UniformRand()* maximumKineticEnergyTransfer;
640  }while(G4UniformRand()*crossSectionMaximum > DifferentialCrossSection(particleDefinition,
641  k,
642  secElecKinetic+waterStructure.IonisationEnergy(shell),
643  shell));
644 
645  return (secElecKinetic);
646 }
647 
648 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
649 
650 // The following section is not used anymore but is kept for memory
651 // GetAngularDistribution()->SampleDirectionForShell is used instead
652 
653 /*
654  void G4DNARuddIonisationModel::RandomizeEjectedElectronDirection(G4ParticleDefinition* particleDefinition,
655  G4double k,
656  G4double secKinetic,
657  G4double & cosTheta,
658  G4double & phi )
659  {
660  G4DNAGenericIonsManager *instance;
661  instance = G4DNAGenericIonsManager::Instance();
662 
663  G4double maxSecKinetic = 0.;
664 
665  if (particleDefinition == G4Proton::ProtonDefinition()
666  || particleDefinition == instance->GetIon("hydrogen"))
667  {
668  maxSecKinetic = 4.* (electron_mass_c2 / proton_mass_c2) * k;
669  }
670 
671  else if (particleDefinition == instance->GetIon("helium")
672  || particleDefinition == instance->GetIon("alpha+")
673  || particleDefinition == instance->GetIon("alpha++"))
674  {
675  maxSecKinetic = 4.* (0.511 / 3728) * k;
676  }
677 
678  phi = twopi * G4UniformRand();
679 
680  // cosTheta = std::sqrt(secKinetic / maxSecKinetic);
681 
682  // Restriction below 100 eV from Emfietzoglou (2000)
683 
684  if (secKinetic>100*eV) cosTheta = std::sqrt(secKinetic / maxSecKinetic);
685  else cosTheta = (2.*G4UniformRand())-1.;
686 
687  }
688  */
689 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
691  G4double k,
692  G4double energyTransfer,
693  G4int ionizationLevelIndex)
694 {
695  // Shells ids are 0 1 2 3 4 (4 is k shell)
696  // !!Attention, "energyTransfer" here is the energy transfered to the electron which means
697  // that the secondary kinetic energy is w = energyTransfer - bindingEnergy
698  //
699  // ds S F1(nu) + w * F2(nu)
700  // ---- = G(k) * ---- -------------------------------------------
701  // dw Bj (1+w)^3 * [1 + exp{alpha * (w - wc) / nu}]
702  //
703  // w is the secondary electron kinetic Energy in eV
704  //
705  // All the other parameters can be found in Rudd's Papers
706  //
707  // M.Eugene Rudd, 1988, User-Friendly model for the energy distribution of
708  // electrons from protons or electron collisions. Nucl. Tracks Rad. Meas.Vol 16 N0 2/3 pp 219-218
709  //
710 
711  const G4int j = ionizationLevelIndex;
712 
713  G4double A1;
714  G4double B1;
715  G4double C1;
716  G4double D1;
717  G4double E1;
718  G4double A2;
719  G4double B2;
720  G4double C2;
721  G4double D2;
722  G4double alphaConst;
723 
724  // const G4double Bj[5] = {12.61*eV, 14.73*eV, 18.55*eV, 32.20*eV, 539.7*eV};
725  // The following values are provided by M. dingfelder (priv. comm)
726  const G4double Bj[5] = { 12.60 * eV, 14.70 * eV, 18.40 * eV, 32.20 * eV, 540
727  * eV };
728 
729  if (j == 4)
730  {
731  //Data For Liquid Water K SHELL from Dingfelder (Protons in Water)
732  A1 = 1.25;
733  B1 = 0.5;
734  C1 = 1.00;
735  D1 = 1.00;
736  E1 = 3.00;
737  A2 = 1.10;
738  B2 = 1.30;
739  C2 = 1.00;
740  D2 = 0.00;
741  alphaConst = 0.66;
742  } else
743  {
744  //Data For Liquid Water from Dingfelder (Protons in Water)
745  A1 = 1.02;
746  B1 = 82.0;
747  C1 = 0.45;
748  D1 = -0.80;
749  E1 = 0.38;
750  A2 = 1.07;
751  // Value provided by M. Dingfelder (priv. comm)
752  B2 = 11.6;
753  //
754  C2 = 0.60;
755  D2 = 0.04;
756  alphaConst = 0.64;
757  }
758 
759  const G4double n = 2.;
760  const G4double Gj[5] = { 0.99, 1.11, 1.11, 0.52, 1. };
761 
764 
765  G4double wBig = (energyTransfer
766  - waterStructure.IonisationEnergy(ionizationLevelIndex));
767  if (wBig < 0)
768  return 0.;
769 
770  G4double w = wBig / Bj[ionizationLevelIndex];
771  // Note that the following (j==4) cases are provided by M. Dingfelder (priv. comm)
772  if (j == 4)
773  w = wBig / waterStructure.IonisationEnergy(ionizationLevelIndex);
774 
775  G4double Ry = 13.6 * eV;
776 
777  G4double tau = 0.;
778 
779  G4bool isProtonOrHydrogen = false;
780  G4bool isHelium = false;
781 
782  if (particleDefinition == G4Proton::ProtonDefinition()
783  || particleDefinition == instance->GetIon("hydrogen"))
784  {
785  isProtonOrHydrogen = true;
786  tau = (electron_mass_c2 / proton_mass_c2) * k;
787  }
788 
789  else if (particleDefinition == instance->GetIon("helium")
790  || particleDefinition == instance->GetIon("alpha+")
791  || particleDefinition == instance->GetIon("alpha++"))
792  {
793  isHelium = true;
794  tau = (0.511 / 3728.) * k;
795  }
796 
797  G4double S = 4. * pi * Bohr_radius * Bohr_radius * n
798  * std::pow((Ry / Bj[ionizationLevelIndex]), 2);
799  if (j == 4)
800  S = 4. * pi * Bohr_radius * Bohr_radius * n
801  * std::pow((Ry / waterStructure.IonisationEnergy(ionizationLevelIndex)),
802  2);
803 
804  G4double v2 = tau / Bj[ionizationLevelIndex];
805  if (j == 4)
806  v2 = tau / waterStructure.IonisationEnergy(ionizationLevelIndex);
807 
808  G4double v = std::sqrt(v2);
809  G4double wc = 4. * v2 - 2. * v - (Ry / (4. * Bj[ionizationLevelIndex]));
810  if (j == 4)
811  wc = 4. * v2 - 2. * v
812  - (Ry / (4. * waterStructure.IonisationEnergy(ionizationLevelIndex)));
813 
814  G4double L1 = (C1 * std::pow(v, (D1))) / (1. + E1 * std::pow(v, (D1 + 4.)));
815  G4double L2 = C2 * std::pow(v, (D2));
816  G4double H1 = (A1 * std::log(1. + v2)) / (v2 + (B1 / v2));
817  G4double H2 = (A2 / v2) + (B2 / (v2 * v2));
818 
819  G4double F1 = L1 + H1;
820  G4double F2 = (L2 * H2) / (L2 + H2);
821 
822  G4double sigma =
823  CorrectionFactor(particleDefinition, k) * Gj[j]
824  * (S / Bj[ionizationLevelIndex])
825  * ((F1 + w * F2)
826  / (std::pow((1. + w), 3)
827  * (1. + std::exp(alphaConst * (w - wc) / v))));
828 
829  if (j == 4)
830  sigma = CorrectionFactor(particleDefinition, k) * Gj[j]
831  * (S / waterStructure.IonisationEnergy(ionizationLevelIndex))
832  * ((F1 + w * F2)
833  / (std::pow((1. + w), 3)
834  * (1. + std::exp(alphaConst * (w - wc) / v))));
835 
836  if ((particleDefinition == instance->GetIon("hydrogen"))
837  && (ionizationLevelIndex == 4))
838  {
839  // sigma = Gj[j] * (S/Bj[ionizationLevelIndex])
840  sigma = Gj[j] * (S / waterStructure.IonisationEnergy(ionizationLevelIndex))
841  * ((F1 + w * F2)
842  / (std::pow((1. + w), 3)
843  * (1. + std::exp(alphaConst * (w - wc) / v))));
844  }
845 // if ( particleDefinition == G4Proton::ProtonDefinition()
846 // || particleDefinition == instance->GetIon("hydrogen")
847 // )
848  if (isProtonOrHydrogen)
849  {
850  return (sigma);
851  }
852 
853  if (particleDefinition == instance->GetIon("alpha++"))
854  {
855  slaterEffectiveCharge[0] = 0.;
856  slaterEffectiveCharge[1] = 0.;
857  slaterEffectiveCharge[2] = 0.;
858  sCoefficient[0] = 0.;
859  sCoefficient[1] = 0.;
860  sCoefficient[2] = 0.;
861  }
862 
863  else if (particleDefinition == instance->GetIon("alpha+"))
864  {
865  slaterEffectiveCharge[0] = 2.0;
866  // The following values are provided by M. Dingfelder (priv. comm)
867  slaterEffectiveCharge[1] = 2.0;
868  slaterEffectiveCharge[2] = 2.0;
869  //
870  sCoefficient[0] = 0.7;
871  sCoefficient[1] = 0.15;
872  sCoefficient[2] = 0.15;
873  }
874 
875  else if (particleDefinition == instance->GetIon("helium"))
876  {
877  slaterEffectiveCharge[0] = 1.7;
878  slaterEffectiveCharge[1] = 1.15;
879  slaterEffectiveCharge[2] = 1.15;
880  sCoefficient[0] = 0.5;
881  sCoefficient[1] = 0.25;
882  sCoefficient[2] = 0.25;
883  }
884 
885 // if ( particleDefinition == instance->GetIon("helium")
886 // || particleDefinition == instance->GetIon("alpha+")
887 // || particleDefinition == instance->GetIon("alpha++")
888 // )
889  if (isHelium)
890  {
891  sigma = Gj[j] * (S / Bj[ionizationLevelIndex])
892  * ((F1 + w * F2)
893  / (std::pow((1. + w), 3)
894  * (1. + std::exp(alphaConst * (w - wc) / v))));
895 
896  if (j == 4)
897  sigma = Gj[j]
898  * (S / waterStructure.IonisationEnergy(ionizationLevelIndex))
899  * ((F1 + w * F2)
900  / (std::pow((1. + w), 3)
901  * (1. + std::exp(alphaConst * (w - wc) / v))));
902 
903  G4double zEff = particleDefinition->GetPDGCharge() / eplus
904  + particleDefinition->GetLeptonNumber();
905 
906  zEff -= (sCoefficient[0]
907  * S_1s(k, energyTransfer, slaterEffectiveCharge[0], 1.)
908  + sCoefficient[1]
909  * S_2s(k, energyTransfer, slaterEffectiveCharge[1], 2.)
910  + sCoefficient[2]
911  * S_2p(k, energyTransfer, slaterEffectiveCharge[2], 2.));
912 
913  return zEff * zEff * sigma;
914  }
915 
916  return 0;
917 }
918 
919 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
920 
922  G4double energyTransferred,
923  G4double slaterEffectiveChg,
924  G4double shellNumber)
925 {
926  // 1 - e^(-2r) * ( 1 + 2 r + 2 r^2)
927  // Dingfelder, in Chattanooga 2005 proceedings, formula (7)
928 
929  G4double r = R(t, energyTransferred, slaterEffectiveChg, shellNumber);
930  G4double value = 1. - std::exp(-2 * r) * ((2. * r + 2.) * r + 1.);
931 
932  return value;
933 }
934 
935 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
936 
938  G4double energyTransferred,
939  G4double slaterEffectiveChg,
940  G4double shellNumber)
941 {
942  // 1 - e^(-2 r) * ( 1 + 2 r + 2 r^2 + 2 r^4)
943  // Dingfelder, in Chattanooga 2005 proceedings, formula (8)
944 
945  G4double r = R(t, energyTransferred, slaterEffectiveChg, shellNumber);
946  G4double value = 1.
947  - std::exp(-2 * r) * (((2. * r * r + 2.) * r + 2.) * r + 1.);
948 
949  return value;
950 
951 }
952 
953 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
954 
956  G4double energyTransferred,
957  G4double slaterEffectiveChg,
958  G4double shellNumber)
959 {
960  // 1 - e^(-2 r) * ( 1 + 2 r + 2 r^2 + 4/3 r^3 + 2/3 r^4)
961  // Dingfelder, in Chattanooga 2005 proceedings, formula (9)
962 
963  G4double r = R(t, energyTransferred, slaterEffectiveChg, shellNumber);
964  G4double value = 1.
965  - std::exp(-2 * r)
966  * ((((2. / 3. * r + 4. / 3.) * r + 2.) * r + 2.) * r + 1.);
967 
968  return value;
969 }
970 
971 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
972 
974  G4double energyTransferred,
975  G4double slaterEffectiveChg,
976  G4double shellNumber)
977 {
978  // tElectron = m_electron / m_alpha * t
979  // Dingfelder, in Chattanooga 2005 proceedings, p 4
980 
981  G4double tElectron = 0.511 / 3728. * t;
982  // The following values are provided by M. Dingfelder (priv. comm)
983  G4double H = 2. * 13.60569172 * eV;
984  G4double value = std::sqrt(2. * tElectron / H) / (energyTransferred / H)
985  * (slaterEffectiveChg / shellNumber);
986 
987  return value;
988 }
989 
990 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
991 
993  G4double k)
994 {
997 
998  if (particleDefinition == G4Proton::Proton())
999  {
1000  return (1.);
1001  } else if (particleDefinition == instance->GetIon("hydrogen"))
1002  {
1003  G4double value = (std::log10(k / eV) - 4.2) / 0.5;
1004  // The following values are provided by M. Dingfelder (priv. comm)
1005  return ((0.6 / (1 + std::exp(value))) + 0.9);
1006  } else
1007  {
1008  return (1.);
1009  }
1010 }
1011 
1012 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1013 
1015  const G4String& particle)
1016 {
1017 
1018  // BEGIN PART 1/2 OF ELECTRON CORRECTION
1019 
1020  // add ONE or TWO electron-water ionisation for alpha+ and helium
1021 
1022  G4int level = 0;
1023 
1024  // Retrieve data table corresponding to the current particle type
1025 
1026  std::map<G4String, G4DNACrossSectionDataSet*, std::less<G4String> >::iterator pos;
1027  pos = tableData.find(particle);
1028 
1029  if (pos != tableData.end())
1030  {
1031  G4DNACrossSectionDataSet* table = pos->second;
1032 
1033  if (table != 0)
1034  {
1035  G4double* valuesBuffer = new G4double[table->NumberOfComponents()];
1036 
1037  const size_t n(table->NumberOfComponents());
1038  size_t i(n);
1039  G4double value = 0.;
1040 
1041  while (i > 0)
1042  {
1043  i--;
1044  valuesBuffer[i] = table->GetComponent(i)->FindValue(k);
1045  value += valuesBuffer[i];
1046  }
1047 
1048  value *= G4UniformRand();
1049 
1050  i = n;
1051 
1052  while (i > 0)
1053  {
1054  i--;
1055 
1056  if (valuesBuffer[i] > value)
1057  {
1058  delete[] valuesBuffer;
1059  return i;
1060  }
1061  value -= valuesBuffer[i];
1062  }
1063 
1064  if (valuesBuffer)
1065  delete[] valuesBuffer;
1066 
1067  }
1068  } else
1069  {
1070  G4Exception("G4DNARuddIonisationModel::RandomSelect",
1071  "em0002",
1073  "Model not applicable to particle type.");
1074  }
1075 
1076  return level;
1077 }
1078 
1079 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1080 
1082 {
1083  G4double sigma = 0.;
1084 
1085  const G4DynamicParticle* particle = track.GetDynamicParticle();
1086  G4double k = particle->GetKineticEnergy();
1087 
1088  G4double lowLim = 0;
1089  G4double highLim = 0;
1090 
1091  const G4String& particleName = particle->GetDefinition()->GetParticleName();
1092 
1093  std::map<G4String, G4double, std::less<G4String> >::iterator pos1;
1094  pos1 = lowEnergyLimit.find(particleName);
1095 
1096  if (pos1 != lowEnergyLimit.end())
1097  {
1098  lowLim = pos1->second;
1099  }
1100 
1101  std::map<G4String, G4double, std::less<G4String> >::iterator pos2;
1102  pos2 = highEnergyLimit.find(particleName);
1103 
1104  if (pos2 != highEnergyLimit.end())
1105  {
1106  highLim = pos2->second;
1107  }
1108 
1109  if (k >= lowLim && k <= highLim)
1110  {
1111  std::map<G4String, G4DNACrossSectionDataSet*, std::less<G4String> >::iterator pos;
1112  pos = tableData.find(particleName);
1113 
1114  if (pos != tableData.end())
1115  {
1116  G4DNACrossSectionDataSet* table = pos->second;
1117  if (table != 0)
1118  {
1119  sigma = table->FindValue(k);
1120  }
1121  } else
1122  {
1123  G4Exception("G4DNARuddIonisationModel::PartialCrossSection",
1124  "em0002",
1126  "Model not applicable to particle type.");
1127  }
1128  }
1129 
1130  return sigma;
1131 }
1132 
1133 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1134 
1136  const G4String& /* particle */)
1137 {
1138  return 0;
1139 }
1140 
const double C2
static const double cm
Definition: G4SIunits.hh:118
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
static const double MeV
Definition: G4SIunits.hh:211
virtual G4double FindValue(G4double x, G4int componentId=0) const =0
static G4LossTableManager * Instance()
const double C1
double S(double temp)
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
const G4DynamicParticle * GetDynamicParticle() const
virtual const G4VEMDataSet * GetComponent(G4int componentId) const
G4ParticleChangeForGamma * fParticleChangeForGamma
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:604
size_t GetIndex() const
Definition: G4Material.hh:262
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
std::map< G4String, G4double, std::less< G4String > > lowEnergyLimit
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:617
const G4double w[NPOINTSGL]
G4double DifferentialCrossSection(G4ParticleDefinition *particleDefinition, G4double k, G4double energyTransfer, G4int shell)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
virtual G4bool LoadData(const G4String &argFileName)
G4ParticleDefinition * GetDefinition() const
G4double S_2s(G4double t, G4double energyTransferred, G4double slaterEffectiveChg, G4double shellNumber)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
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 G4ThreeVector & SampleDirectionForShell(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, G4int shellID, const G4Material *)
G4DNARuddIonisationModel(const G4ParticleDefinition *p=0, const G4String &nam="DNARuddIonisationModel")
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
bool G4bool
Definition: G4Types.hh:79
const G4ThreeVector & GetMomentumDirection() const
const std::vector< G4double > * fpWaterDensity
G4double R(G4double t, G4double energyTransferred, G4double slaterEffectiveChg, G4double shellNumber)
std::map< G4String, G4double, std::less< G4String > > highEnergyLimit
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4int RandomSelect(G4double energy, const G4String &particle)
static G4DNAGenericIonsManager * Instance(void)
const G4int n
virtual G4double FindValue(G4double e, G4int componentId=0) const
virtual size_t NumberOfComponents(void) const
G4double Sum(G4double energy, const G4String &particle)
G4double PartialCrossSection(const G4Track &track)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4DNAChemistryManager * Instance()
static const double eV
Definition: G4SIunits.hh:212
static G4DNAMolecularMaterial * Instance()
static const double pi
Definition: G4SIunits.hh:74
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
Method used by DNA physics model to create a water molecule.
G4double S_1s(G4double t, G4double energyTransferred, G4double slaterEffectiveChg, G4double shellNumber)
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:624
G4DNAWaterIonisationStructure waterStructure
G4double CorrectionFactor(G4ParticleDefinition *particleDefinition, G4double k)
const G4Track * GetCurrentTrack() const
G4double RandomizeEjectedElectronEnergy(G4ParticleDefinition *particleDefinition, G4double incomingParticleEnergy, G4int shell)
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
static const double m
Definition: G4SIunits.hh:128
static MCTruthManager * instance
G4VAtomDeexcitation * AtomDeexcitation()
static const double keV
Definition: G4SIunits.hh:213
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:781
static const double eplus
Definition: G4SIunits.hh:196
G4double GetPDGCharge() const
G4ThreeVector G4ParticleMomentum
G4double bindingEnergy(G4int A, G4int Z)
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
G4double S_2p(G4double t, G4double energyTransferred, G4double slaterEffectiveChg, G4double shellNumber)
G4VAtomDeexcitation * fAtomDeexcitation
static const G4double pos
G4AtomicShellEnumerator
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134
const G4Material * GetMaterial() const
G4ParticleDefinition * GetIon(const G4String &name)