Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4StopElementSelector.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$
27//
28// File: G4StopElementSelector
29//
30// Author: V.Ivanchenko ([email protected])
31//
32// Creation date: 2 April 2000
33//
34// Modifications:
35// 18/08/2000 V.Ivanchenko Update description
36// 17/05/2006 V.Ivanchenko Cleanup
37// 02/10/2007 V.Ivanchenko Fixed typo in computation of Lambda-factor
38// proposed by Victor Pec
39//
40//---------------------------------------------------------------------
41
44#include "G4SystemOfUnits.hh"
45#include "Randomize.hh"
46#include "G4Material.hh"
47
48//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
49
50// constructor
52{ }
53
54//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
55
56// destructor
58{ }
59
60
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
62
64{
65 // Fermi-Teller Z-low of mu- capture and exceptions
66 // for halogens and oxigen.
67 // N.C.Mukhopadhyay Phys. Rep. 30 (1977) 1.
68 G4int i;
69 G4double Z;
70 const G4int numberOfElements = aMaterial->GetNumberOfElements();
71 const G4ElementVector* theElementVector = aMaterial->GetElementVector();
72
73 if(1 == numberOfElements) return (*theElementVector)[0];
74
75 const G4double* theAtomicNumberDensity = aMaterial->GetAtomicNumDensityVector();
76
77 G4double sum = 0.0;
78 for ( i=0; i < numberOfElements; i++ ) {
79
80 Z = (*theElementVector)[i]->GetZ();
81
82 // Halogens
83 if( (9.0 == Z) || (17.0 == Z) || (35.0 == Z) || (53.0 == Z) || (85.0 == Z) ) {
84 sum += 0.66 * Z * theAtomicNumberDensity[i] ;
85
86 // Oxigen
87 } else if( 8.0 == Z ) {
88 sum += 0.56 * Z * theAtomicNumberDensity[i] ;
89
90 // Others
91 } else {
92 sum += Z * theAtomicNumberDensity[i] ;
93 }
94 }
95
96 G4double random = G4UniformRand() * sum;
97 sum = 0.0 ;
98 i = -1;
99
100 // Selection of element
101 do {
102 i++;
103 Z = (*theElementVector)[i]->GetZ();
104
105 // Galogens
106 if( (9.0 == Z) || (17.0 == Z) || (35.0 == Z) || (53.0 == Z) || (85.0 == Z) ) {
107 sum += 0.66 * Z * theAtomicNumberDensity[i] ;
108
109 // Oxigen
110 } else if( 8.0 == Z ) {
111 sum += 0.56 * Z * theAtomicNumberDensity[i] ;
112
113 // Others
114 } else {
115 sum += Z * theAtomicNumberDensity[i] ;
116 }
117 } while ( (sum < random) && (i < numberOfElements - 1) );
118
119 return (*theElementVector)[i];
120}
121
122//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
123
125{
126 // Initialized data
127
128 // static std::vector<G4double> zeff(100);
129 static G4double zeff[100] = {
130 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,
131 12.91,13.64,14.24,14.89,15.53,16.15,16.75,17.38,18.04,18.49,
132 19.06,19.59,20.1,20.66,21.12,21.61,22.02,22.43,22.84,23.24,
133 23.65,24.06,24.47,24.85,25.23,25.61,25.99,26.37,26.69,27.,
134 27.32,27.63,27.95,28.2,28.42,28.64,28.79,29.03,29.27,29.51,
135 29.75,29.99,30.2,30.36,30.53,30.69,30.85,31.01,31.18,31.34,
136 31.48,31.62,31.76,31.9,32.05,32.19,32.33,32.47,32.61,32.76,
137 32.94,33.11,33.29,33.46,33.64,33.81,34.21,34.18,34.,34.1,
138 34.21,34.31,34.42,34.52,34.63,34.73,34.84,34.94,35.04,35.15,
139 35.25,35.36,35.46,35.57,35.67,35.78 };
140
141 // Mu- capture data from B.B.Balashov, G.Ya.Korenman, P.A.Eramgan
142 // Atomizdat, 1978. (Experimental capture velocities)
143 // Data for Hydrogen from Phys. Rev. Lett. 99(2007)032002
144
145 const size_t ListZE = 66;
146 static G4int ListZExp[ListZE] = {1,
147 3, 4, 5, 6, 7, 8, 9, 10, 11, 12,
148 13, 14, 15, 16, 17, 18, 19, 20, 22, 23,
149 24, 25, 26, 27, 28, 31, 32, 33, 34, 37,
150 38, 39, 40, 41, 42, 45, 46, 47, 48, 49,
151 50, 51, 52, 53, 55, 56, 57, 58, 59, 60,
152 62, 64, 65, 67, 72, 73, 74, 80, 81, 82,
153 83, 90, 92, 93};
154
155 static G4double ListCaptureVel[ListZE] = {0.000725,
156 0.0057, 0.010, 0.0258, 0.0371, 0.0644,
157 0.0974, 0.144, 0.250, 0.386, 0.479,
158 0.700, 0.849, 1.119, 1.338, 1.40,
159 1.30, 1.98, 2.45, 2.60, 3.19,
160 3.29, 3.91, 4.41, 4.96, 5.74,
161 5.68, 5.53, 6.06, 5.69, 6.89,
162 7.25, 7.89, 8.59, 10.40, 9.22,
163 10.01, 10.00, 10.88, 10.62, 11.37,
164 10.68, 10.49, 9.06, 11.20, 10.98,
165 10.18, 10.71, 11.44, 13.45, 12.32,
166 12.22, 12.09, 12.73, 12.95, 13.03,
167 12.86, 13.13, 13.39, 12.74, 13.78,
168 13.02, 13.26, 13.10, 14.00, 14.70};
169
170
171 // Local variables
172 G4double zeff2, xmu, a2ze, r1, r2;
173 G4double lambda;
174
175 // == Effective charges from Ford and Wills Nucl Phys 35(1962)295.
176 // == Untabulated charges are interpolated.
177 // == Mu capture lifetime (Goulard and Primakoff PRC10(1974)2034.
178
179 G4int i = G4int(Z) - 1 ;
180 if(i > 99) i = 99;
181
182 const G4double b0a = -.03;
183 const G4double b0b = -.25;
184 const G4double b0c = 3.24;
185 const G4double t1 = 875.e-10;
186 r1 = zeff[i];
187 zeff2 = r1 * r1;
188 // ^-4 -> ^-5 suggested by user
189 xmu = zeff2 * 2.663e-5;
190 a2ze = 0.5 * A / Z;
191 r2 = 1.0 - xmu;
192 lambda = t1 * zeff2 * zeff2 * (r2 * r2) * (1.0 - (1.0 - xmu) * .75704) *
193 (a2ze * b0a + 1.0 - (a2ze - 1.0) * b0b -
194 (2.0 * (A - Z) + std::abs(a2ze - 1.) ) * b0c / (A * 4.) );
195
196 // == Mu capture data are taken if exist
197 for (unsigned int j = 0; j < ListZE; j++) {
198 if( ListZExp[j] == i + 1) {
199 lambda = ListCaptureVel[j] / microsecond;
200 break;
201 }
202 }
203
204 return lambda;
205}
206
207//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
208
210{
211 // Decay time on K-shell
212 // N.C.Mukhopadhyay Phys. Rep. 30 (1977) 1.
213
214 G4double lambda = 1.0;
215 if(Z > 1) {
216 G4double x = Z*fine_structure_const;
217 lambda -= 2.5 * x * x;
218 if( 0.5 > lambda ) { lambda = 0.5; }
219 }
220 return lambda * 0.445164 / microsecond;
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
284
285
286
287
288
289
std::vector< G4Element * > G4ElementVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4UniformRand()
Definition: Randomize.hh:53
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:189
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:215
G4Element * GetElement(const G4Material *aMaterial)
G4double GetMuonDecayRate(G4double Z, G4double A)
G4double GetMuonCaptureRate(G4double Z, G4double A)