Geant4  10.02
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 92859 2015-09-18 07:58:30Z 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  G4int Z = 8;
496  {
498 
499  if (ionizationShell <5 && ionizationShell >1)
500  {
501  as = G4AtomicShellEnumerator(4-ionizationShell);
502  }
503  else if (ionizationShell <2)
504  {
505  as = G4AtomicShellEnumerator(3);
506  }
507 
508  // DEBUG
509  // if (ionizationShell == 4) {
510  //
511  // G4cout << "Z: " << Z << " as: " << as
512  // << " ionizationShell: " << ionizationShell << " bindingEnergy: "<< bindingEnergy/eV << G4endl;
513  // G4cout << "Press <Enter> key to continue..." << G4endl;
514  // G4cin.ignore();
515  // }
516 
517  const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
518  secNumberInit = fvect->size();
519  fAtomDeexcitation->GenerateParticles(fvect, shell, Z, 0, 0);
520  secNumberFinal = fvect->size();
521  }
522 
523  G4double secondaryKinetic = RandomizeEjectedElectronEnergy(definition,k,ionizationShell);
524 
525  G4ThreeVector deltaDirection =
526  GetAngularDistribution()->SampleDirectionForShell(particle, secondaryKinetic,
527  Z, ionizationShell,
528  couple->GetMaterial());
529 
530  // Ignored for ions on electrons
531  /*
532  G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
533 
534  G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
535  G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
536  G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
537  G4double finalMomentum = std::sqrt(finalPx*finalPx+finalPy*finalPy+finalPz*finalPz);
538  finalPx /= finalMomentum;
539  finalPy /= finalMomentum;
540  finalPz /= finalMomentum;
541 
542  G4ThreeVector direction;
543  direction.set(finalPx,finalPy,finalPz);
544 
545  fParticleChangeForGamma->ProposeMomentumDirection(direction.unit()) ;
546  */
548 
549  G4double scatteredEnergy = k-bindingEnergy-secondaryKinetic;
550  G4double deexSecEnergy = 0;
551  for (G4int j=secNumberInit; j < secNumberFinal; j++)
552  {
553 
554  deexSecEnergy = deexSecEnergy + (*fvect)[j]->GetKineticEnergy();
555 
556  }
557 
558  if (!statCode)
559  {
561  fParticleChangeForGamma->ProposeLocalEnergyDeposit(k-scatteredEnergy-secondaryKinetic-deexSecEnergy);
562  }
563  else
564  {
567  }
568 
569  // debug
570  // k-scatteredEnergy-secondaryKinetic-deexSecEnergy = k-(k-bindingEnergy-secondaryKinetic)-secondaryKinetic-deexSecEnergy =
571  // = k-k+bindingEnergy+secondaryKinetic-secondaryKinetic-deexSecEnergy=
572  // = bindingEnergy-deexSecEnergy
573  // SO deexSecEnergy=0 => LocalEnergyDeposit = bindingEnergy
574 
575  G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic);
576  fvect->push_back(dp);
577 
578  const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
580  ionizationShell,
581  theIncomingTrack);
582  }
583 
584  // SI - not useful since low energy of model is 0 eV
585 
586  if (k < lowLim)
587  {
591  }
592 
593 }
594 
595 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
596 
598  G4double k,
599  G4int shell)
600 {
601  G4double maximumKineticEnergyTransfer = 0.;
602 
605 
606  if (particleDefinition == G4Proton::ProtonDefinition()
607  || particleDefinition == instance->GetIon("hydrogen"))
608  {
609  maximumKineticEnergyTransfer = 4. * (electron_mass_c2 / proton_mass_c2) * k;
610  }
611 
612  else if (particleDefinition == instance->GetIon("helium")
613  || particleDefinition == instance->GetIon("alpha+")
614  || particleDefinition == instance->GetIon("alpha++"))
615  {
616  maximumKineticEnergyTransfer = 4. * (0.511 / 3728) * k;
617  }
618 
619  G4double crossSectionMaximum = 0.;
620 
621  for (G4double value = waterStructure.IonisationEnergy(shell);
622  value <= 5. * waterStructure.IonisationEnergy(shell) && k >= value;
623  value += 0.1 * eV)
624  {
625  G4double differentialCrossSection =
626  DifferentialCrossSection(particleDefinition, k, value, shell);
627  if (differentialCrossSection >= crossSectionMaximum)
628  crossSectionMaximum = differentialCrossSection;
629  }
630 
631  G4double secElecKinetic = 0.;
632 
633  do
634  {
635  secElecKinetic = G4UniformRand()* maximumKineticEnergyTransfer;
636  }while(G4UniformRand()*crossSectionMaximum > DifferentialCrossSection(particleDefinition,
637  k,
638  secElecKinetic+waterStructure.IonisationEnergy(shell),
639  shell));
640 
641  return (secElecKinetic);
642 }
643 
644 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
645 
646 // The following section is not used anymore but is kept for memory
647 // GetAngularDistribution()->SampleDirectionForShell is used instead
648 
649 /*
650  void G4DNARuddIonisationModel::RandomizeEjectedElectronDirection(G4ParticleDefinition* particleDefinition,
651  G4double k,
652  G4double secKinetic,
653  G4double & cosTheta,
654  G4double & phi )
655  {
656  G4DNAGenericIonsManager *instance;
657  instance = G4DNAGenericIonsManager::Instance();
658 
659  G4double maxSecKinetic = 0.;
660 
661  if (particleDefinition == G4Proton::ProtonDefinition()
662  || particleDefinition == instance->GetIon("hydrogen"))
663  {
664  maxSecKinetic = 4.* (electron_mass_c2 / proton_mass_c2) * k;
665  }
666 
667  else if (particleDefinition == instance->GetIon("helium")
668  || particleDefinition == instance->GetIon("alpha+")
669  || particleDefinition == instance->GetIon("alpha++"))
670  {
671  maxSecKinetic = 4.* (0.511 / 3728) * k;
672  }
673 
674  phi = twopi * G4UniformRand();
675 
676  // cosTheta = std::sqrt(secKinetic / maxSecKinetic);
677 
678  // Restriction below 100 eV from Emfietzoglou (2000)
679 
680  if (secKinetic>100*eV) cosTheta = std::sqrt(secKinetic / maxSecKinetic);
681  else cosTheta = (2.*G4UniformRand())-1.;
682 
683  }
684  */
685 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
687  G4double k,
688  G4double energyTransfer,
689  G4int ionizationLevelIndex)
690 {
691  // Shells ids are 0 1 2 3 4 (4 is k shell)
692  // !!Attention, "energyTransfer" here is the energy transfered to the electron which means
693  // that the secondary kinetic energy is w = energyTransfer - bindingEnergy
694  //
695  // ds S F1(nu) + w * F2(nu)
696  // ---- = G(k) * ---- -------------------------------------------
697  // dw Bj (1+w)^3 * [1 + exp{alpha * (w - wc) / nu}]
698  //
699  // w is the secondary electron kinetic Energy in eV
700  //
701  // All the other parameters can be found in Rudd's Papers
702  //
703  // M.Eugene Rudd, 1988, User-Friendly model for the energy distribution of
704  // electrons from protons or electron collisions. Nucl. Tracks Rad. Meas.Vol 16 N0 2/3 pp 219-218
705  //
706 
707  const G4int j = ionizationLevelIndex;
708 
709  G4double A1;
710  G4double B1;
711  G4double C1;
712  G4double D1;
713  G4double E1;
714  G4double A2;
715  G4double B2;
716  G4double C2;
717  G4double D2;
718  G4double alphaConst;
719 
720  // const G4double Bj[5] = {12.61*eV, 14.73*eV, 18.55*eV, 32.20*eV, 539.7*eV};
721  // The following values are provided by M. dingfelder (priv. comm)
722  const G4double Bj[5] = { 12.60 * eV, 14.70 * eV, 18.40 * eV, 32.20 * eV, 540
723  * eV };
724 
725  if (j == 4)
726  {
727  //Data For Liquid Water K SHELL from Dingfelder (Protons in Water)
728  A1 = 1.25;
729  B1 = 0.5;
730  C1 = 1.00;
731  D1 = 1.00;
732  E1 = 3.00;
733  A2 = 1.10;
734  B2 = 1.30;
735  C2 = 1.00;
736  D2 = 0.00;
737  alphaConst = 0.66;
738  } else
739  {
740  //Data For Liquid Water from Dingfelder (Protons in Water)
741  A1 = 1.02;
742  B1 = 82.0;
743  C1 = 0.45;
744  D1 = -0.80;
745  E1 = 0.38;
746  A2 = 1.07;
747  // Value provided by M. Dingfelder (priv. comm)
748  B2 = 11.6;
749  //
750  C2 = 0.60;
751  D2 = 0.04;
752  alphaConst = 0.64;
753  }
754 
755  const G4double n = 2.;
756  const G4double Gj[5] = { 0.99, 1.11, 1.11, 0.52, 1. };
757 
760 
761  G4double wBig = (energyTransfer
762  - waterStructure.IonisationEnergy(ionizationLevelIndex));
763  if (wBig < 0)
764  return 0.;
765 
766  G4double w = wBig / Bj[ionizationLevelIndex];
767  // Note that the following (j==4) cases are provided by M. Dingfelder (priv. comm)
768  if (j == 4)
769  w = wBig / waterStructure.IonisationEnergy(ionizationLevelIndex);
770 
771  G4double Ry = 13.6 * eV;
772 
773  G4double tau = 0.;
774 
775  G4bool isProtonOrHydrogen = false;
776  G4bool isHelium = false;
777 
778  if (particleDefinition == G4Proton::ProtonDefinition()
779  || particleDefinition == instance->GetIon("hydrogen"))
780  {
781  isProtonOrHydrogen = true;
782  tau = (electron_mass_c2 / proton_mass_c2) * k;
783  }
784 
785  else if (particleDefinition == instance->GetIon("helium")
786  || particleDefinition == instance->GetIon("alpha+")
787  || particleDefinition == instance->GetIon("alpha++"))
788  {
789  isHelium = true;
790  tau = (0.511 / 3728.) * k;
791  }
792 
793  G4double S = 4. * pi * Bohr_radius * Bohr_radius * n
794  * std::pow((Ry / Bj[ionizationLevelIndex]), 2);
795  if (j == 4)
796  S = 4. * pi * Bohr_radius * Bohr_radius * n
797  * std::pow((Ry / waterStructure.IonisationEnergy(ionizationLevelIndex)),
798  2);
799 
800  G4double v2 = tau / Bj[ionizationLevelIndex];
801  if (j == 4)
802  v2 = tau / waterStructure.IonisationEnergy(ionizationLevelIndex);
803 
804  G4double v = std::sqrt(v2);
805  G4double wc = 4. * v2 - 2. * v - (Ry / (4. * Bj[ionizationLevelIndex]));
806  if (j == 4)
807  wc = 4. * v2 - 2. * v
808  - (Ry / (4. * waterStructure.IonisationEnergy(ionizationLevelIndex)));
809 
810  G4double L1 = (C1 * std::pow(v, (D1))) / (1. + E1 * std::pow(v, (D1 + 4.)));
811  G4double L2 = C2 * std::pow(v, (D2));
812  G4double H1 = (A1 * std::log(1. + v2)) / (v2 + (B1 / v2));
813  G4double H2 = (A2 / v2) + (B2 / (v2 * v2));
814 
815  G4double F1 = L1 + H1;
816  G4double F2 = (L2 * H2) / (L2 + H2);
817 
818  G4double sigma =
819  CorrectionFactor(particleDefinition, k) * Gj[j]
820  * (S / Bj[ionizationLevelIndex])
821  * ((F1 + w * F2)
822  / (std::pow((1. + w), 3)
823  * (1. + std::exp(alphaConst * (w - wc) / v))));
824 
825  if (j == 4)
826  sigma = CorrectionFactor(particleDefinition, k) * Gj[j]
827  * (S / waterStructure.IonisationEnergy(ionizationLevelIndex))
828  * ((F1 + w * F2)
829  / (std::pow((1. + w), 3)
830  * (1. + std::exp(alphaConst * (w - wc) / v))));
831 
832  if ((particleDefinition == instance->GetIon("hydrogen"))
833  && (ionizationLevelIndex == 4))
834  {
835  // sigma = Gj[j] * (S/Bj[ionizationLevelIndex])
836  sigma = Gj[j] * (S / waterStructure.IonisationEnergy(ionizationLevelIndex))
837  * ((F1 + w * F2)
838  / (std::pow((1. + w), 3)
839  * (1. + std::exp(alphaConst * (w - wc) / v))));
840  }
841 // if ( particleDefinition == G4Proton::ProtonDefinition()
842 // || particleDefinition == instance->GetIon("hydrogen")
843 // )
844  if (isProtonOrHydrogen)
845  {
846  return (sigma);
847  }
848 
849  if (particleDefinition == instance->GetIon("alpha++"))
850  {
851  slaterEffectiveCharge[0] = 0.;
852  slaterEffectiveCharge[1] = 0.;
853  slaterEffectiveCharge[2] = 0.;
854  sCoefficient[0] = 0.;
855  sCoefficient[1] = 0.;
856  sCoefficient[2] = 0.;
857  }
858 
859  else if (particleDefinition == instance->GetIon("alpha+"))
860  {
861  slaterEffectiveCharge[0] = 2.0;
862  // The following values are provided by M. Dingfelder (priv. comm)
863  slaterEffectiveCharge[1] = 2.0;
864  slaterEffectiveCharge[2] = 2.0;
865  //
866  sCoefficient[0] = 0.7;
867  sCoefficient[1] = 0.15;
868  sCoefficient[2] = 0.15;
869  }
870 
871  else if (particleDefinition == instance->GetIon("helium"))
872  {
873  slaterEffectiveCharge[0] = 1.7;
874  slaterEffectiveCharge[1] = 1.15;
875  slaterEffectiveCharge[2] = 1.15;
876  sCoefficient[0] = 0.5;
877  sCoefficient[1] = 0.25;
878  sCoefficient[2] = 0.25;
879  }
880 
881 // if ( particleDefinition == instance->GetIon("helium")
882 // || particleDefinition == instance->GetIon("alpha+")
883 // || particleDefinition == instance->GetIon("alpha++")
884 // )
885  if (isHelium)
886  {
887  sigma = Gj[j] * (S / Bj[ionizationLevelIndex])
888  * ((F1 + w * F2)
889  / (std::pow((1. + w), 3)
890  * (1. + std::exp(alphaConst * (w - wc) / v))));
891 
892  if (j == 4)
893  sigma = Gj[j]
894  * (S / waterStructure.IonisationEnergy(ionizationLevelIndex))
895  * ((F1 + w * F2)
896  / (std::pow((1. + w), 3)
897  * (1. + std::exp(alphaConst * (w - wc) / v))));
898 
899  G4double zEff = particleDefinition->GetPDGCharge() / eplus
900  + particleDefinition->GetLeptonNumber();
901 
902  zEff -= (sCoefficient[0]
903  * S_1s(k, energyTransfer, slaterEffectiveCharge[0], 1.)
904  + sCoefficient[1]
905  * S_2s(k, energyTransfer, slaterEffectiveCharge[1], 2.)
906  + sCoefficient[2]
907  * S_2p(k, energyTransfer, slaterEffectiveCharge[2], 2.));
908 
909  return zEff * zEff * sigma;
910  }
911 
912  return 0;
913 }
914 
915 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
916 
918  G4double energyTransferred,
919  G4double slaterEffectiveChg,
920  G4double shellNumber)
921 {
922  // 1 - e^(-2r) * ( 1 + 2 r + 2 r^2)
923  // Dingfelder, in Chattanooga 2005 proceedings, formula (7)
924 
925  G4double r = R(t, energyTransferred, slaterEffectiveChg, shellNumber);
926  G4double value = 1. - std::exp(-2 * r) * ((2. * r + 2.) * r + 1.);
927 
928  return value;
929 }
930 
931 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
932 
934  G4double energyTransferred,
935  G4double slaterEffectiveChg,
936  G4double shellNumber)
937 {
938  // 1 - e^(-2 r) * ( 1 + 2 r + 2 r^2 + 2 r^4)
939  // Dingfelder, in Chattanooga 2005 proceedings, formula (8)
940 
941  G4double r = R(t, energyTransferred, slaterEffectiveChg, shellNumber);
942  G4double value = 1.
943  - std::exp(-2 * r) * (((2. * r * r + 2.) * r + 2.) * r + 1.);
944 
945  return value;
946 
947 }
948 
949 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
950 
952  G4double energyTransferred,
953  G4double slaterEffectiveChg,
954  G4double shellNumber)
955 {
956  // 1 - e^(-2 r) * ( 1 + 2 r + 2 r^2 + 4/3 r^3 + 2/3 r^4)
957  // Dingfelder, in Chattanooga 2005 proceedings, formula (9)
958 
959  G4double r = R(t, energyTransferred, slaterEffectiveChg, shellNumber);
960  G4double value = 1.
961  - std::exp(-2 * r)
962  * ((((2. / 3. * r + 4. / 3.) * r + 2.) * r + 2.) * r + 1.);
963 
964  return value;
965 }
966 
967 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
968 
970  G4double energyTransferred,
971  G4double slaterEffectiveChg,
972  G4double shellNumber)
973 {
974  // tElectron = m_electron / m_alpha * t
975  // Dingfelder, in Chattanooga 2005 proceedings, p 4
976 
977  G4double tElectron = 0.511 / 3728. * t;
978  // The following values are provided by M. Dingfelder (priv. comm)
979  G4double H = 2. * 13.60569172 * eV;
980  G4double value = std::sqrt(2. * tElectron / H) / (energyTransferred / H)
981  * (slaterEffectiveChg / shellNumber);
982 
983  return value;
984 }
985 
986 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
987 
989  G4double k)
990 {
993 
994  if (particleDefinition == G4Proton::Proton())
995  {
996  return (1.);
997  } else if (particleDefinition == instance->GetIon("hydrogen"))
998  {
999  G4double value = (std::log10(k / eV) - 4.2) / 0.5;
1000  // The following values are provided by M. Dingfelder (priv. comm)
1001  return ((0.6 / (1 + std::exp(value))) + 0.9);
1002  } else
1003  {
1004  return (1.);
1005  }
1006 }
1007 
1008 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1009 
1011  const G4String& particle)
1012 {
1013 
1014  // BEGIN PART 1/2 OF ELECTRON CORRECTION
1015 
1016  // add ONE or TWO electron-water ionisation for alpha+ and helium
1017 
1018  G4int level = 0;
1019 
1020  // Retrieve data table corresponding to the current particle type
1021 
1022  std::map<G4String, G4DNACrossSectionDataSet*, std::less<G4String> >::iterator pos;
1023  pos = tableData.find(particle);
1024 
1025  if (pos != tableData.end())
1026  {
1027  G4DNACrossSectionDataSet* table = pos->second;
1028 
1029  if (table != 0)
1030  {
1031  G4double* valuesBuffer = new G4double[table->NumberOfComponents()];
1032 
1033  const size_t n(table->NumberOfComponents());
1034  size_t i(n);
1035  G4double value = 0.;
1036 
1037  while (i > 0)
1038  {
1039  i--;
1040  valuesBuffer[i] = table->GetComponent(i)->FindValue(k);
1041  value += valuesBuffer[i];
1042  }
1043 
1044  value *= G4UniformRand();
1045 
1046  i = n;
1047 
1048  while (i > 0)
1049  {
1050  i--;
1051 
1052  if (valuesBuffer[i] > value)
1053  {
1054  delete[] valuesBuffer;
1055  return i;
1056  }
1057  value -= valuesBuffer[i];
1058  }
1059 
1060  if (valuesBuffer)
1061  delete[] valuesBuffer;
1062 
1063  }
1064  } else
1065  {
1066  G4Exception("G4DNARuddIonisationModel::RandomSelect",
1067  "em0002",
1069  "Model not applicable to particle type.");
1070  }
1071 
1072  return level;
1073 }
1074 
1075 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1076 
1078 {
1079  G4double sigma = 0.;
1080 
1081  const G4DynamicParticle* particle = track.GetDynamicParticle();
1082  G4double k = particle->GetKineticEnergy();
1083 
1084  G4double lowLim = 0;
1085  G4double highLim = 0;
1086 
1087  const G4String& particleName = particle->GetDefinition()->GetParticleName();
1088 
1089  std::map<G4String, G4double, std::less<G4String> >::iterator pos1;
1090  pos1 = lowEnergyLimit.find(particleName);
1091 
1092  if (pos1 != lowEnergyLimit.end())
1093  {
1094  lowLim = pos1->second;
1095  }
1096 
1097  std::map<G4String, G4double, std::less<G4String> >::iterator pos2;
1098  pos2 = highEnergyLimit.find(particleName);
1099 
1100  if (pos2 != highEnergyLimit.end())
1101  {
1102  highLim = pos2->second;
1103  }
1104 
1105  if (k >= lowLim && k <= highLim)
1106  {
1107  std::map<G4String, G4DNACrossSectionDataSet*, std::less<G4String> >::iterator pos;
1108  pos = tableData.find(particleName);
1109 
1110  if (pos != tableData.end())
1111  {
1112  G4DNACrossSectionDataSet* table = pos->second;
1113  if (table != 0)
1114  {
1115  sigma = table->FindValue(k);
1116  }
1117  } else
1118  {
1119  G4Exception("G4DNARuddIonisationModel::PartialCrossSection",
1120  "em0002",
1122  "Model not applicable to particle type.");
1123  }
1124  }
1125 
1126  return sigma;
1127 }
1128 
1129 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1130 
1132  const G4String& /* particle */)
1133 {
1134  return 0;
1135 }
1136 
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)