Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4DNADingfelderChargeIncreaseModel.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: G4DNADingfelderChargeIncreaseModel.cc 96060 2016-03-11 12:58:04Z gcosmo $
27 //
28 
30 #include "G4PhysicalConstants.hh"
31 #include "G4SystemOfUnits.hh"
33 
34 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
35 
36 using namespace std;
37 
38 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
39 
41  const G4String& nam) :
42  G4VEmModel(nam), isInitialised(false)
43 {
44  // nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
45  fpMolWaterDensity = 0;
46 
47  numberOfPartialCrossSections[0] = 0;
48  numberOfPartialCrossSections[1] = 0;
49 
50  verboseLevel = 0;
51  // Verbosity scale:
52  // 0 = nothing
53  // 1 = warning for energy non-conservation
54  // 2 = details of energy budget
55  // 3 = calculation of cross sections, file openings, sampling of atoms
56  // 4 = entering in methods
57 
58  if (verboseLevel > 0)
59  {
60  G4cout << "Dingfelder charge increase model is constructed " << G4endl;
61  }
63 
64  // Selection of stationary mode
65 
66  statCode = false;
67 }
68 
69 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
70 
72 {}
73 
74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
75 
77  const G4DataVector& /*cuts*/)
78 {
79 
80  if (verboseLevel > 3)
81  {
82  G4cout << "Calling G4DNADingfelderChargeIncreaseModel::Initialise()"
83  << G4endl;
84  }
85 
86  // Energy limits
87 
88  G4DNAGenericIonsManager *instance;
90  G4ParticleDefinition* hydrogenDef = instance->GetIon("hydrogen");
91  G4ParticleDefinition* alphaPlusDef = instance->GetIon("alpha+");
92  G4ParticleDefinition* heliumDef = instance->GetIon("helium");
93 
94  G4String hydrogen;
95  G4String alphaPlus;
96  G4String helium;
97 
98  // LIMITS
99 
100  hydrogen = hydrogenDef->GetParticleName();
101  lowEnergyLimit[hydrogen] = 100. * eV;
102  highEnergyLimit[hydrogen] = 100. * MeV;
103 
104  alphaPlus = alphaPlusDef->GetParticleName();
105  lowEnergyLimit[alphaPlus] = 1. * keV;
106  highEnergyLimit[alphaPlus] = 400. * MeV;
107 
108  helium = heliumDef->GetParticleName();
109  lowEnergyLimit[helium] = 1. * keV;
110  highEnergyLimit[helium] = 400. * MeV;
111 
112  //
113 
114  if (particle==hydrogenDef)
115  {
116  SetLowEnergyLimit(lowEnergyLimit[hydrogen]);
117  SetHighEnergyLimit(highEnergyLimit[hydrogen]);
118  }
119 
120  if (particle==alphaPlusDef)
121  {
122  SetLowEnergyLimit(lowEnergyLimit[alphaPlus]);
123  SetHighEnergyLimit(highEnergyLimit[alphaPlus]);
124  }
125 
126  if (particle==heliumDef)
127  {
128  SetLowEnergyLimit(lowEnergyLimit[helium]);
129  SetHighEnergyLimit(highEnergyLimit[helium]);
130  }
131 
132  // Final state
133 
134  //ALPHA+
135 
136  f0[0][0]=1.;
137  a0[0][0]=2.25;
138  a1[0][0]=-0.75;
139  b0[0][0]=-32.10;
140  c0[0][0]=0.600;
141  d0[0][0]=2.40;
142  x0[0][0]=4.60;
143 
144  x1[0][0]=-1.;
145  b1[0][0]=-1.;
146 
147  numberOfPartialCrossSections[0]=1;
148 
149  //HELIUM
150 
151  f0[0][1]=1.;
152  a0[0][1]=2.25;
153  a1[0][1]=-0.75;
154  b0[0][1]=-30.93;
155  c0[0][1]=0.590;
156  d0[0][1]=2.35;
157  x0[0][1]=4.29;
158 
159  f0[1][1]=1.;
160  a0[1][1]=2.25;
161  a1[1][1]=-0.75;
162  b0[1][1]=-32.61;
163  c0[1][1]=0.435;
164  d0[1][1]=2.70;
165  x0[1][1]=4.45;
166 
167  x1[0][1]=-1.;
168  b1[0][1]=-1.;
169 
170  x1[1][1]=-1.;
171  b1[1][1]=-1.;
172 
173  numberOfPartialCrossSections[1]=2;
174 
175  //
176 
177  if( verboseLevel>0 )
178  {
179  G4cout << "Dingfelder charge increase model is initialized " << G4endl
180  << "Energy range: "
181  << LowEnergyLimit() / keV << " keV - "
182  << HighEnergyLimit() / MeV << " MeV for "
183  << particle->GetParticleName()
184  << G4endl;
185  }
186 
187  // Initialize water density pointer
189 
190  if (isInitialised)
191  { return;}
193  isInitialised = true;
194 }
195 
196 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
197 
199  const G4ParticleDefinition* particleDefinition,
200  G4double k,
201  G4double,
202  G4double)
203 {
204  if (verboseLevel > 3)
205  {
206  G4cout
207  << "Calling CrossSectionPerVolume() of G4DNADingfelderChargeIncreaseModel"
208  << G4endl;
209  }
210 
211  // Calculate total cross section for model
212 
213  G4DNAGenericIonsManager *instance;
215 
216  if (
217  particleDefinition != instance->GetIon("hydrogen")
218  &&
219  particleDefinition != instance->GetIon("alpha+")
220  &&
221  particleDefinition != instance->GetIon("helium")
222  )
223 
224  return 0;
225 
226  G4double lowLim = 0;
227  G4double highLim = 0;
228  G4double totalCrossSection = 0.;
229 
230  G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
231 
232  if(waterDensity!= 0.0)
233  // if (material == nistwater || material->GetBaseMaterial() == nistwater)
234  {
235  const G4String& particleName = particleDefinition->GetParticleName();
236 
237  std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
238  pos1 = lowEnergyLimit.find(particleName);
239 
240  if (pos1 != lowEnergyLimit.end())
241  {
242  lowLim = pos1->second;
243  }
244 
245  std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
246  pos2 = highEnergyLimit.find(particleName);
247 
248  if (pos2 != highEnergyLimit.end())
249  {
250  highLim = pos2->second;
251  }
252 
253  if (k >= lowLim && k < highLim)
254  {
255  //HYDROGEN
256  if (particleDefinition == instance->GetIon("hydrogen"))
257  {
258  const G4double aa = 2.835;
259  const G4double bb = 0.310;
260  const G4double cc = 2.100;
261  const G4double dd = 0.760;
262  const G4double fac = 1.0e-18;
263  const G4double rr = 13.606 * eV;
264 
266  G4double x = t / rr;
267  G4double temp = 4.0 * pi * Bohr_radius/nm * Bohr_radius/nm * fac;
268  G4double sigmal = temp * cc * (std::pow(x,dd));
269  G4double sigmah = temp * (aa * std::log(1.0 + x) + bb) / x;
270  totalCrossSection = 1.0/(1.0/sigmal + 1.0/sigmah) *m*m;
271  }
272  else
273  {
274  totalCrossSection = Sum(k,particleDefinition);
275  }
276  }
277 
278  if (verboseLevel > 2)
279  {
280  G4cout << "__________________________________" << G4endl;
281  G4cout << "G4DNADingfelderChargeIncreaseModel - XS INFO START" << G4endl;
282  G4cout << "Kinetic energy(eV)=" << k/eV << " particle : " << particleName << G4endl;
283  G4cout << "Cross section per water molecule (cm^2)=" << totalCrossSection/cm/cm << G4endl;
284  G4cout << "Cross section per water molecule (cm^-1)=" << totalCrossSection*waterDensity/(1./cm) << G4endl;
285  // G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
286  G4cout << "G4DNADingfelderChargeIncreaseModel - XS INFO END" << G4endl;
287  }
288 
289  }
290 
291  return totalCrossSection*waterDensity;
292 // return totalCrossSection*material->GetAtomicNumDensityVector()[1];
293 
294 }
295 
296 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
297 
299  G4DynamicParticle*>* fvect,
300  const G4MaterialCutsCouple* /*couple*/,
301  const G4DynamicParticle* aDynamicParticle,
302  G4double,
303  G4double)
304 {
305  if (verboseLevel > 3)
306  {
307  G4cout
308  << "Calling SampleSecondaries() of G4DNADingfelderChargeIncreaseModel"
309  << G4endl;
310  }
311 
313 
314  G4ParticleDefinition* definition = aDynamicParticle->GetDefinition();
315 
316  G4double particleMass = definition->GetPDGMass();
317 
318  G4double inK = aDynamicParticle->GetKineticEnergy();
319 
320  G4int finalStateIndex = RandomSelect(inK,definition);
321 
322  G4int n = NumberOfFinalStates(definition,finalStateIndex);
323 
324  G4double outK = 0.;
325 
326  if (!statCode) outK = inK - IncomingParticleBindingEnergyConstant(definition,finalStateIndex);
327 
328  else outK = inK;
329 
330  if (statCode)
331  fParticleChangeForGamma->ProposeLocalEnergyDeposit(IncomingParticleBindingEnergyConstant(definition,finalStateIndex));
332 
334 
335  G4DNAGenericIonsManager* instance;
337 
338  G4double electronK;
339  if (definition == instance->GetIon("hydrogen")) electronK = inK*electron_mass_c2/proton_mass_c2;
340  else electronK = inK*electron_mass_c2/(particleMass);
341 
342  if (outK<0)
343  {
344  G4Exception("G4DNADingfelderChargeIncreaseModel::SampleSecondaries","em0004",
345  FatalException,"Final kinetic energy is negative.");
346  }
347 
348  G4DynamicParticle* dp = new G4DynamicParticle(OutgoingParticleDefinition(definition,finalStateIndex),
349  aDynamicParticle->GetMomentumDirection(),
350  outK);
351 
352  fvect->push_back(dp);
353 
354  n = n - 1;
355 
356  while (n>0)
357  {
358  n--;
359  fvect->push_back(new G4DynamicParticle
360  (G4Electron::Electron(), aDynamicParticle->GetMomentumDirection(), electronK) );
361  }
362 }
363 
364 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
365 
366 G4int G4DNADingfelderChargeIncreaseModel::NumberOfFinalStates(G4ParticleDefinition* particleDefinition,
367  G4int finalStateIndex)
368 
369 {
370  G4DNAGenericIonsManager* instance;
372  if (particleDefinition == instance->GetIon("hydrogen"))
373  return 2;
374  if (particleDefinition == instance->GetIon("alpha+"))
375  return 2;
376 
377  if (particleDefinition == instance->GetIon("helium"))
378  {
379  if (finalStateIndex == 0)
380  return 2;
381  return 3;
382  }
383  return 0;
384 }
385 
386 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
387 
388 G4ParticleDefinition* G4DNADingfelderChargeIncreaseModel::OutgoingParticleDefinition(G4ParticleDefinition* particleDefinition,
389  G4int finalStateIndex)
390 {
392  if (particleDefinition == instance->GetIon("hydrogen"))
393  return G4Proton::Proton();
394  if (particleDefinition == instance->GetIon("alpha+"))
395  return instance->GetIon("alpha++");
396 
397  if (particleDefinition == instance->GetIon("helium"))
398  {
399  if (finalStateIndex == 0)
400  return instance->GetIon("alpha+");
401  return instance->GetIon("alpha++");
402  }
403  return 0;
404 }
405 
406 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
407 
408 G4double G4DNADingfelderChargeIncreaseModel::IncomingParticleBindingEnergyConstant(G4ParticleDefinition* particleDefinition,
409  G4int finalStateIndex)
410 {
412  if (particleDefinition == instance->GetIon("hydrogen"))
413  return 13.6 * eV;
414 
415  if (particleDefinition == instance->GetIon("alpha+"))
416  {
417  // Binding energy for He+ -> He++ + e- 54.509 eV
418  // Binding energy for He -> He+ + e- 24.587 eV
419  return 54.509 * eV;
420  }
421 
422  if (particleDefinition == instance->GetIon("helium"))
423  {
424  // Binding energy for He+ -> He++ + e- 54.509 eV
425  // Binding energy for He -> He+ + e- 24.587 eV
426 
427  if (finalStateIndex == 0)
428  return 24.587 * eV;
429  return (54.509 + 24.587) * eV;
430  }
431 
432  return 0;
433 }
434 
435 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
436 
437 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
438 
439 G4double G4DNADingfelderChargeIncreaseModel::PartialCrossSection(G4double k,
440  G4int index,
441  const G4ParticleDefinition* particleDefinition)
442 {
443  G4int particleTypeIndex = 0;
444  G4DNAGenericIonsManager *instance;
446 
447  if (particleDefinition == instance->GetIon("alpha+"))
448  particleTypeIndex = 0;
449  if (particleDefinition == instance->GetIon("helium"))
450  particleTypeIndex = 1;
451 
452  //
453  // sigma(T) = f0 10 ^ y(log10(T/eV))
454  //
455  // / a0 x + b0 x < x0
456  // |
457  // y(x) = < a0 x + b0 - c0 (x - x0)^d0 x0 <= x < x1
458  // |
459  // \ a1 x + b1 x >= x1
460  //
461  //
462  // f0, a0, a1, b0, b1, c0, d0, x0, x1 are parameters that change for protons and helium (0, +, ++)
463  //
464  // f0 has been added to the code in order to manage partial (shell-dependent) cross sections (if no shell dependence is present. f0=1. Sum of f0 over the considered shells should give 1)
465  //
466  // From Rad. Phys. and Chem. 59 (2000) 255-275, M. Dingfelder et al.
467  // Inelastic-collision cross sections of liquid water for interactions of energetic proton
468  //
469 
470  if (x1[index][particleTypeIndex] < x0[index][particleTypeIndex])
471  {
472  //
473  // if x1 < x0 means that x1 and b1 will be calculated with the following formula (this piece of code is run on all alphas and not on protons)
474  //
475  // x1 = x0 + ((a0 - a1)/(c0 * d0)) ^ (1 / (d0 - 1))
476  //
477  // b1 = (a0 - a1) * x1 + b0 - c0 * (x1 - x0) ^ d0
478  //
479 
480  x1[index][particleTypeIndex] = x0[index][particleTypeIndex]
481  + std::pow((a0[index][particleTypeIndex] - a1[index][particleTypeIndex])
482  / (c0[index][particleTypeIndex]
483  * d0[index][particleTypeIndex]),
484  1. / (d0[index][particleTypeIndex] - 1.));
485  b1[index][particleTypeIndex] = (a0[index][particleTypeIndex]
486  - a1[index][particleTypeIndex]) * x1[index][particleTypeIndex]
487  + b0[index][particleTypeIndex]
488  - c0[index][particleTypeIndex]
489  * std::pow(x1[index][particleTypeIndex]
490  - x0[index][particleTypeIndex],
491  d0[index][particleTypeIndex]);
492  }
493 
494  G4double x(std::log10(k / eV));
495  G4double y;
496 
497  if (x < x0[index][particleTypeIndex])
498  y = a0[index][particleTypeIndex] * x + b0[index][particleTypeIndex];
499  else if (x < x1[index][particleTypeIndex])
500  y = a0[index][particleTypeIndex] * x + b0[index][particleTypeIndex]
501  - c0[index][particleTypeIndex]
502  * std::pow(x - x0[index][particleTypeIndex],
503  d0[index][particleTypeIndex]);
504  else
505  y = a1[index][particleTypeIndex] * x + b1[index][particleTypeIndex];
506 
507  return f0[index][particleTypeIndex] * std::pow(10., y) * m * m;
508 
509 }
510 
511 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
512 
513 G4int G4DNADingfelderChargeIncreaseModel::RandomSelect(G4double k,
514  const G4ParticleDefinition* particleDefinition)
515 {
516  G4int particleTypeIndex = 0;
517  G4DNAGenericIonsManager *instance;
519 
520  if (particleDefinition == instance->GetIon("hydrogen"))
521  return 0;
522  if (particleDefinition == instance->GetIon("alpha+"))
523  particleTypeIndex = 0;
524  if (particleDefinition == instance->GetIon("helium"))
525  particleTypeIndex = 1;
526 
527  const G4int n = numberOfPartialCrossSections[particleTypeIndex];
528  G4double* values(new G4double[n]);
529  G4double value = 0;
530  G4int i = n;
531 
532  while (i > 0)
533  {
534  i--;
535  values[i] = PartialCrossSection(k, i, particleDefinition);
536  value += values[i];
537  }
538 
539  value *= G4UniformRand();
540 
541  i = n;
542  while (i > 0)
543  {
544  i--;
545 
546  if (values[i] > value)
547  break;
548 
549  value -= values[i];
550  }
551 
552  delete[] values;
553 
554  return i;
555 }
556 
557 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
558 
559 G4double G4DNADingfelderChargeIncreaseModel::Sum(G4double k,
560  const G4ParticleDefinition* particleDefinition)
561 {
562  G4int particleTypeIndex = 0;
563  G4DNAGenericIonsManager *instance;
565 
566  if (particleDefinition == instance->GetIon("alpha+"))
567  particleTypeIndex = 0;
568  if (particleDefinition == instance->GetIon("helium"))
569  particleTypeIndex = 1;
570 
571  G4double totalCrossSection = 0.;
572 
573  for (G4int i = 0; i < numberOfPartialCrossSections[particleTypeIndex]; i++)
574  {
575  totalCrossSection += PartialCrossSection(k, i, particleDefinition);
576  }
577  return totalCrossSection;
578 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
static constexpr double proton_mass_c2
G4double GetKineticEnergy() const
static constexpr double Bohr_radius
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:602
size_t GetIndex() const
Definition: G4Material.hh:262
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
static constexpr double electron_mass_c2
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:724
static ulg bb
Definition: csz_inflate.cc:344
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
static constexpr double m
Definition: G4SIunits.hh:129
const XML_Char int const XML_Char * value
Definition: expat.h:331
const G4ThreeVector & GetMomentumDirection() const
static constexpr double cm
Definition: G4SIunits.hh:119
static G4Proton * Proton()
Definition: G4Proton.cc:93
static constexpr double eV
Definition: G4SIunits.hh:215
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
static G4DNAGenericIonsManager * Instance(void)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4DNAMolecularMaterial * Instance()
G4double GetPDGMass() const
static constexpr double nm
Definition: G4SIunits.hh:112
static G4Electron * Electron()
Definition: G4Electron.cc:94
static const G4double fac
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
static constexpr double pi
Definition: G4SIunits.hh:75
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:731
G4DNADingfelderChargeIncreaseModel(const G4ParticleDefinition *p=0, const G4String &nam="DNADingfelderChargeIncreaseModel")
static constexpr double keV
Definition: G4SIunits.hh:216
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132
G4ParticleDefinition * GetIon(const G4String &name)