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