source: trunk/source/processes/electromagnetic/lowenergy/include/G4IonParametrisedLossModel.icc@ 1036

Last change on this file since 1036 was 1007, checked in by garnier, 17 years ago

update to geant4.9.2

File size: 7.8 KB
RevLine 
[968]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//
28// ===========================================================================
29// GEANT4 class
30//
31// Class: G4IonParametrisedLossModel
32//
33// Base class: G4VEmModel (utils)
34//
35// Author: Anton Lechner (Anton.Lechner@cern.ch)
36//
37// First implementation: 10. 11. 2008
38//
[1007]39// Modifications:
[968]40//
41//
42// Class description:
43// Model for computing the energy loss of ions by employing a
44// parameterisation of dE/dx tables (default ICRU 73 tables). For
45// ion-material combinations and/or projectile energies not covered
46// by this model, the G4BraggIonModel and G4BetheBloch models are
47// employed.
48//
49// Comments:
50//
51// ===========================================================================
52
53
54inline G4double G4IonParametrisedLossModel::DeltaRayMeanEnergyTransferRate(
55 const G4Material* material,
56 const G4ParticleDefinition* particle,
57 G4double kineticEnergy,
58 G4double cutEnergy) {
59
60 // ############## Mean energy transferred to delta-rays ###################
61 // Computes the mean energy transfered to delta-rays per unit length,
62 // considering only delta-rays with energies above the energy threshold
63 // (energy cut)
64 //
65 // The mean energy transfer rate is derived by using the differential
66 // cross section given in the references below.
67 //
68 // See Geant4 physics reference manual (version 9.1), section 9.1.3
69 //
70 // Ref.: W.M. Yao et al, Jour. of Phys. G 33 (2006) 1.
71 // B. Rossi, High energy particles, New York, NY: Prentice-Hall (1952).
72 //
73 // (Implementation adapted from G4BraggIonModel)
74
75
76 // *** Variables:
77 // kineticEnergy = kinetic energy of projectile
78 // totEnergy = total energy of projectile, i.e. kinetic energy
79 // plus rest energy (Mc^2)
80 // betaSquared = beta of projectile squared, calculated as
81 // beta^2 = 1 - 1 / (E/Mc^2)^2
82 // = T * ( E + Mc^2 ) / E^2
83 // where T = kineticEnergy, E = totEnergy
84 // cutEnergy = energy threshold for secondary particle production
85 // i.e. energy cut, below which energy transfered to
86 // electrons is treated as continuous loss of projectile
87 // maxKinEnergy = maximum energy transferable to secondary electrons
88 // meanRate = mean kinetic energy of delta ray (per unit length)
89 // (above cutEnergy)
90
91 G4double meanRate = 0.0;
92
93 G4double maxKinEnergy = MaxSecondaryEnergy(particle, kineticEnergy);
94
95 if (cutEnergy < maxKinEnergy) {
96
97 G4double totalEnergy = kineticEnergy + cacheMass;
98 G4double betaSquared = kineticEnergy *
99 (totalEnergy + cacheMass) / (totalEnergy * totalEnergy);
100
101 G4double cutMaxEnergyRatio = cutEnergy / maxKinEnergy;
102
103 meanRate =
104 (- std::log(cutMaxEnergyRatio) - (1.0 - cutMaxEnergyRatio) * betaSquared) *
105 twopi_mc2_rcl2 *
106 (material->GetTotNbOfElectPerVolume()) / betaSquared;
107
108 meanRate *= GetChargeSquareRatio(particle, material, kineticEnergy);
109 }
110
111 return meanRate;
112}
113
114
115inline
116G4double G4IonParametrisedLossModel::MaxSecondaryEnergy(
117 const G4ParticleDefinition* particle,
118 G4double kineticEnergy) {
119
120 // ############## Maximum energy of secondaries ##########################
121 // Function computes maximum energy of secondary electrons which are
122 // released by an ion
123 //
124 // See Geant4 physics reference manual (version 9.1), section 9.1.1
125 //
126 // Ref.: W.M. Yao et al, Jour. of Phys. G 33 (2006) 1.
127 // C.Caso et al. (Part. Data Group), Europ. Phys. Jour. C 3 1 (1998).
128 // B. Rossi, High energy particles, New York, NY: Prentice-Hall (1952).
129 //
130 // (Implementation adapted from G4BraggIonModel)
131
132 if(particle != cacheParticle) UpdateCache(particle);
133
134 G4double tau = kineticEnergy/cacheMass;
135 G4double tmax = 2.0 * electron_mass_c2 * tau * (tau + 2.) /
136 (1. + 2.0 * (tau + 1.) * cacheElecMassRatio +
137 cacheElecMassRatio * cacheElecMassRatio);
138
139 return tmax;
140}
141
142
143inline
144void G4IonParametrisedLossModel::UpdateCache(
145 const G4ParticleDefinition* particle) {
146
147 cacheParticle = particle;
148 cacheMass = particle -> GetPDGMass();
149 cacheElecMassRatio = electron_mass_c2 / cacheMass;
150 G4double q = particle -> GetPDGCharge() / eplus;
151 cacheChargeSquare = q * q;
152}
153
154
155inline
156G4double G4IonParametrisedLossModel::GetChargeSquareRatio(
157 const G4ParticleDefinition* particle,
158 const G4Material* material,
159 G4double kineticEnergy) { // Kinetic energy
160
161 G4double chargeSquareRatio = corrections ->
162 EffectiveChargeSquareRatio(particle,
163 material,
164 kineticEnergy);
165 corrFactor = chargeSquareRatio *
166 corrections -> EffectiveChargeCorrection(particle,
167 material,
168 kineticEnergy);
169 return corrFactor;
170}
171
172
173inline
174G4double G4IonParametrisedLossModel::GetParticleCharge(
175 const G4ParticleDefinition* particle,
176 const G4Material* material,
177 G4double kineticEnergy) { // Kinetic energy
178
179 return corrections -> GetParticleCharge(particle, material, kineticEnergy);
180}
181
182
183inline
184LossTableList::iterator G4IonParametrisedLossModel::IsApplicable(
185 const G4ParticleDefinition* particle, // Projectile (ion)
186 const G4Material* material) { // Target material
187
[1007]188 LossTableList::iterator iter = lossTableList.begin();
[968]189 LossTableList::iterator iterTables = lossTableList.begin();
190 LossTableList::iterator iterTables_end = lossTableList.end();
191
192 for(;iterTables != iterTables_end; iterTables++) {
193 G4bool isApplicable = (*iterTables) ->
194 IsApplicable(particle, material);
195 if(isApplicable) {
196 iter = iterTables;
197 break;
198 }
199 }
200
201 return iter;
202}
Note: See TracBrowser for help on using the repository browser.