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 | // $Id: G4CrossSectionChargeIncrease.cc,v 1.3 2007/12/10 16:31:21 gunter Exp $ |
---|
28 | // GEANT4 tag $Name: $ |
---|
29 | // |
---|
30 | // Contact Author: Maria Grazia Pia (Maria.Grazia.Pia@cern.ch) |
---|
31 | // |
---|
32 | // Reference: TNS Geant4-DNA paper |
---|
33 | // Reference for implementation model: NIM. 155, pp. 145-156, 1978 |
---|
34 | |
---|
35 | // History: |
---|
36 | // ----------- |
---|
37 | // Date Name Modification |
---|
38 | // 28 Apr 2007 M.G. Pia Created in compliance with design described in TNS paper |
---|
39 | // |
---|
40 | // ------------------------------------------------------------------- |
---|
41 | |
---|
42 | // Class description: |
---|
43 | // Geant4-DNA Cross total cross section for electron elastic scattering in water |
---|
44 | // Reference: TNS Geant4-DNA paper |
---|
45 | // S. Chauvie et al., Geant4 physics processes for microdosimetry simulation: |
---|
46 | // design foundation and implementation of the first set of models, |
---|
47 | // IEEE Trans. Nucl. Sci., vol. 54, no. 6, Dec. 2007. |
---|
48 | // Further documentation available from http://www.ge.infn.it/geant4/dna |
---|
49 | |
---|
50 | // ------------------------------------------------------------------- |
---|
51 | |
---|
52 | |
---|
53 | #include "G4CrossSectionChargeIncrease.hh" |
---|
54 | #include "G4Track.hh" |
---|
55 | #include "G4DynamicParticle.hh" |
---|
56 | #include "G4ParticleDefinition.hh" |
---|
57 | #include "G4DNAGenericIonsManager.hh" |
---|
58 | |
---|
59 | G4CrossSectionChargeIncrease::G4CrossSectionChargeIncrease() |
---|
60 | { |
---|
61 | // Default energy limits (defined for protection against anomalous behaviour only) |
---|
62 | name = "ChargeIncrease"; |
---|
63 | lowEnergyLimitDefault = 1 * keV; |
---|
64 | highEnergyLimitDefault = 10 * MeV; |
---|
65 | |
---|
66 | G4DNAGenericIonsManager *instance; |
---|
67 | instance = G4DNAGenericIonsManager::Instance(); |
---|
68 | G4ParticleDefinition* hydrogenDef = instance->GetIon("hydrogen"); |
---|
69 | G4ParticleDefinition* alphaPlusDef = instance->GetIon("alpha+"); |
---|
70 | G4ParticleDefinition* heliumDef = instance->GetIon("helium"); |
---|
71 | |
---|
72 | G4String hydrogen; |
---|
73 | G4String alphaPlus; |
---|
74 | G4String helium; |
---|
75 | |
---|
76 | if (hydrogenDef != 0) |
---|
77 | { |
---|
78 | hydrogen = hydrogenDef->GetParticleName(); |
---|
79 | lowEnergyLimit[hydrogen] = 1. * keV; |
---|
80 | highEnergyLimit[hydrogen] = 10. * MeV; |
---|
81 | } |
---|
82 | else |
---|
83 | { |
---|
84 | G4Exception("G4CrossSectionChargeIncrease Constructor: hydrogen is not defined"); |
---|
85 | } |
---|
86 | |
---|
87 | if (alphaPlusDef != 0) |
---|
88 | { |
---|
89 | alphaPlus = alphaPlusDef->GetParticleName(); |
---|
90 | lowEnergyLimit[alphaPlus] = 1. * keV; |
---|
91 | highEnergyLimit[alphaPlus] = 10. * MeV; |
---|
92 | } |
---|
93 | else |
---|
94 | { |
---|
95 | G4Exception("G4CrossSectionChargeIncrease Constructor: alphaPlus is not defined"); |
---|
96 | } |
---|
97 | |
---|
98 | if (heliumDef != 0) |
---|
99 | { |
---|
100 | helium = heliumDef->GetParticleName(); |
---|
101 | lowEnergyLimit[helium] = 1. * keV; |
---|
102 | highEnergyLimit[helium] = 10. * MeV; |
---|
103 | } |
---|
104 | else |
---|
105 | { |
---|
106 | G4Exception("G4CrossSectionChargeIncrease Constructor: helium is not defined"); |
---|
107 | } |
---|
108 | |
---|
109 | } |
---|
110 | |
---|
111 | |
---|
112 | G4CrossSectionChargeIncrease::~G4CrossSectionChargeIncrease() |
---|
113 | {} |
---|
114 | |
---|
115 | |
---|
116 | G4double G4CrossSectionChargeIncrease::CrossSection(const G4Track& track) |
---|
117 | { |
---|
118 | G4double lowLim = lowEnergyLimitDefault; |
---|
119 | G4double highLim = highEnergyLimitDefault; |
---|
120 | |
---|
121 | const G4DynamicParticle* particle = track.GetDynamicParticle(); |
---|
122 | G4double k = particle->GetKineticEnergy(); |
---|
123 | |
---|
124 | const G4String& particleName = particle->GetDefinition()->GetParticleName(); |
---|
125 | |
---|
126 | G4DNAGenericIonsManager *instance; |
---|
127 | instance = G4DNAGenericIonsManager::Instance(); |
---|
128 | |
---|
129 | const G4ParticleDefinition* particleDefinition = track.GetDefinition(); |
---|
130 | |
---|
131 | if ( |
---|
132 | particleDefinition != instance->GetIon("hydrogen") |
---|
133 | && |
---|
134 | particleDefinition != instance->GetIon("alpha+") |
---|
135 | && |
---|
136 | particleDefinition != instance->GetIon("helium") |
---|
137 | ) |
---|
138 | |
---|
139 | G4Exception("G4CrossSectionChargeIncrease: attempting to calculate cross section for wrong particle"); |
---|
140 | |
---|
141 | |
---|
142 | // Retrieve energy limits for the current particle type |
---|
143 | |
---|
144 | std::map< G4String,G4double,std::less<G4String> >::iterator pos1; |
---|
145 | pos1 = lowEnergyLimit.find(particleName); |
---|
146 | |
---|
147 | // Lower limit |
---|
148 | if (pos1 != lowEnergyLimit.end()) |
---|
149 | { |
---|
150 | lowLim = pos1->second; |
---|
151 | } |
---|
152 | |
---|
153 | // Upper limit |
---|
154 | std::map< G4String,G4double,std::less<G4String> >::iterator pos2; |
---|
155 | pos2 = highEnergyLimit.find(particleName); |
---|
156 | |
---|
157 | if (pos2 != highEnergyLimit.end()) |
---|
158 | { |
---|
159 | highLim = pos2->second; |
---|
160 | } |
---|
161 | |
---|
162 | G4double totalCrossSection = 0.; |
---|
163 | |
---|
164 | if (k >= lowLim && k <= highLim) |
---|
165 | { |
---|
166 | //HYDROGEN |
---|
167 | if (particleDefinition == instance->GetIon("hydrogen")) |
---|
168 | { |
---|
169 | const G4double aa = 2.835; |
---|
170 | const G4double bb = 0.310; |
---|
171 | const G4double cc = 2.100; |
---|
172 | const G4double dd = 0.760; |
---|
173 | const G4double fac = 1.0e-18; |
---|
174 | const G4double rr = 13.606 * eV; |
---|
175 | |
---|
176 | G4double t = k / (proton_mass_c2/electron_mass_c2); |
---|
177 | G4double x = t / rr; |
---|
178 | G4double temp = 4.0 * pi * Bohr_radius/nm * Bohr_radius/nm * fac; |
---|
179 | G4double sigmal = temp * cc * (std::pow(x,dd)); |
---|
180 | G4double sigmah = temp * (aa * std::log(1.0 + x) + bb) / x; |
---|
181 | totalCrossSection = 1.0/(1.0/sigmal + 1.0/sigmah) *m*m; |
---|
182 | } |
---|
183 | else |
---|
184 | { |
---|
185 | totalCrossSection = partialCrossSection.Sum(k,particleDefinition); |
---|
186 | } |
---|
187 | } |
---|
188 | |
---|
189 | return totalCrossSection; |
---|
190 | } |
---|
191 | |
---|