source: trunk/source/processes/hadronic/stopping/src/G4StopElementSelector.cc@ 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.4 KB
RevLine 
[819]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: G4StopElementSelector.cc,v 1.16 2007/10/02 18:27:43 vnivanch Exp $
[1007]27// GEANT4 tag $Name: geant4-09-02 $
[819]28//
29// File: G4StopElementSelector
30//
31// Author: V.Ivanchenko (Vladimir.Ivanchenko@cern.ch)
32//
33// Creation date: 2 April 2000
34//
35// Modifications:
36// 18/08/2000 V.Ivanchenko Update description
37// 17/05/2006 V.Ivanchenko Cleanup
38// 02/10/2007 V.Ivanchenko Fixed typo in computation of Lambda-factor
39// proposed by Victor Pec
40//
41//---------------------------------------------------------------------
42
43#include "G4StopElementSelector.hh"
44#include "Randomize.hh"
45#include "G4Material.hh"
46
47//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
48
49// constructor
50G4StopElementSelector::G4StopElementSelector()
51{ }
52
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
54
55// destructor
56G4StopElementSelector::~G4StopElementSelector()
57{ }
58
59
60//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
61
62G4Element* G4StopElementSelector::GetElement(const G4Material* aMaterial)
63{
64 // Fermi-Teller Z-low of mu- capture and exceptions
65 // for halogens and oxigen.
66 // N.C.Mukhopadhyay Phys. Rep. 30 (1977) 1.
67 G4int i;
68 G4double Z;
69 const G4int numberOfElements = aMaterial->GetNumberOfElements();
70 const G4ElementVector* theElementVector = aMaterial->GetElementVector();
71
72 if(1 == numberOfElements) return (*theElementVector)[0];
73
74 const G4double* theAtomicNumberDensity = aMaterial->GetAtomicNumDensityVector();
75
76 G4double sum = 0.0;
77 for ( i=0; i < numberOfElements; i++ ) {
78
79 Z = (*theElementVector)[i]->GetZ();
80
81 // Halogens
82 if( (9.0 == Z) || (17.0 == Z) || (35.0 == Z) || (53.0 == Z) || (85.0 == Z) ) {
83 sum += 0.66 * Z * theAtomicNumberDensity[i] ;
84
85 // Oxigen
86 } else if( 8.0 == Z ) {
87 sum += 0.56 * Z * theAtomicNumberDensity[i] ;
88
89 // Others
90 } else {
91 sum += Z * theAtomicNumberDensity[i] ;
92 }
93 }
94
95 G4double random = G4UniformRand() * sum;
96 sum = 0.0 ;
97 i = -1;
98
99 // Selection of element
100 do {
101 i++;
102 Z = (*theElementVector)[i]->GetZ();
103
104 // Galogens
105 if( (9.0 == Z) || (17.0 == Z) || (35.0 == Z) || (53.0 == Z) || (85.0 == Z) ) {
106 sum += 0.66 * Z * theAtomicNumberDensity[i] ;
107
108 // Oxigen
109 } else if( 8.0 == Z ) {
110 sum += 0.56 * Z * theAtomicNumberDensity[i] ;
111
112 // Others
113 } else {
114 sum += Z * theAtomicNumberDensity[i] ;
115 }
116 } while ( (sum < random) && (i < numberOfElements - 1) );
117
118 return (*theElementVector)[i];
119}
120
121//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
122
123G4double G4StopElementSelector::GetMuonCaptureRate(G4double Z, G4double A)
124{
125 // Initialized data
126
127 // static std::vector<G4double> zeff(100);
128 static G4double zeff[100] = {
129 1.,1.98,2.95,3.89,4.8,5.72,6.61,7.49,8.32,9.12,9.95,10.69,11.48,12.22,
130 12.91,13.64,14.24,14.89,15.53,16.15,16.75,17.38,18.04,18.49,
131 19.06,19.59,20.1,20.66,21.12,21.61,22.02,22.43,22.84,23.24,
132 23.65,24.06,24.47,24.85,25.23,25.61,25.99,26.37,26.69,27.,
133 27.32,27.63,27.95,28.2,28.42,28.64,28.79,29.03,29.27,29.51,
134 29.75,29.99,30.2,30.36,30.53,30.69,30.85,31.01,31.18,31.34,
135 31.48,31.62,31.76,31.9,32.05,32.19,32.33,32.47,32.61,32.76,
136 32.94,33.11,33.29,33.46,33.64,33.81,34.21,34.18,34.,34.1,
137 34.21,34.31,34.42,34.52,34.63,34.73,34.84,34.94,35.04,35.15,
138 35.25,35.36,35.46,35.57,35.67,35.78 };
139
140 // Mu- capture data from B.B.Balashov, G.Ya.Korenman, P.A.Eramgan
141 // Atomizdat, 1978. (Experimental capture velocities)
142
143 const size_t ListZE = 65;
144 static G4int ListZExp[ListZE] = {
145 3, 4, 5, 6, 7, 8, 9, 10, 11, 12,
146 13, 14, 15, 16, 17, 18, 19, 20, 22, 23,
147 24, 25, 26, 27, 28, 31, 32, 33, 34, 37,
148 38, 39, 40, 41, 42, 45, 46, 47, 48, 49,
149 50, 51, 52, 53, 55, 56, 57, 58, 59, 60,
150 62, 64, 65, 67, 72, 73, 74, 80, 81, 82,
151 83, 90, 92, 93};
152
153 static G4double ListCaptureVel[ListZE] = {
154 0.0057, 0.010, 0.0258, 0.0371, 0.0644,
155 0.0974, 0.144, 0.250, 0.386, 0.479,
156 0.700, 0.849, 1.119, 1.338, 1.40,
157 1.30, 1.98, 2.45, 2.60, 3.19,
158 3.29, 3.91, 4.41, 4.96, 5.74,
159 5.68, 5.53, 6.06, 5.69, 6.89,
160 7.25, 7.89, 8.59, 10.40, 9.22,
161 10.01, 10.00, 10.88, 10.62, 11.37,
162 10.68, 10.49, 9.06, 11.20, 10.98,
163 10.18, 10.71, 11.44, 13.45, 12.32,
164 12.22, 12.09, 12.73, 12.95, 13.03,
165 12.86, 13.13, 13.39, 12.74, 13.78,
166 13.02, 13.26, 13.10, 14.00, 14.70};
167
168
169 // Local variables
170 G4double zeff2, xmu, a2ze, r1, r2;
171 G4double lambda;
172
173 // == Effective charges from Ford and Wills Nucl Phys 35(1962)295.
174 // == Untabulated charges are interpolated.
175 // == Mu capture lifetime (Goulard and Primakoff PRC10(1974)2034.
176
177 G4int i = G4int(Z) - 1 ;
178 if(i > 99) i = 99;
179
180 const G4double b0a = -.03;
181 const G4double b0b = -.25;
182 const G4double b0c = 3.24;
183 const G4double t1 = 875.e-10;
184 r1 = zeff[i];
185 zeff2 = r1 * r1;
186 // ^-4 -> ^-5 suggested by user
187 xmu = zeff2 * 2.663e-5;
188 a2ze = 0.5 * A / Z;
189 r2 = 1.0 - xmu;
190 lambda = t1 * zeff2 * zeff2 * (r2 * r2) * (1.0 - (1.0 - xmu) * .75704) *
191 (a2ze * b0a + 1.0 - (a2ze - 1.0) * b0b -
192 (2.0 * (A - Z) + std::abs(a2ze - 1.) ) * b0c / (A * 4.) );
193
194 // == Mu capture data are taken if exist
195 for (unsigned int j = 0; j < ListZE; j++) {
196 if( ListZExp[j] == i + 1) {
197 lambda = ListCaptureVel[j] / microsecond;
198 break;
199 }
200 }
201
202 return lambda;
203}
204
205//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
206
207G4double G4StopElementSelector::GetMuonDecayRate(G4double Z, G4double /* A */)
208{
209 // Decay time on K-shell
210 // N.C.Mukhopadhyay Phys. Rep. 30 (1977) 1.
211
212 G4double lambda = 1.0 - 2.5 * Z * Z / (137.0*137.0);
213 if( 0.5 > lambda ) lambda = 0.5;
214 return lambda * 0.445 / microsecond;
215}
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
Note: See TracBrowser for help on using the repository browser.