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: G4MuBremsstrahlungModel.cc,v 1.24 2007/11/08 11:48:28 vnivanch Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-01-patch-02 $ |
---|
28 | // |
---|
29 | // ------------------------------------------------------------------- |
---|
30 | // |
---|
31 | // GEANT4 Class file |
---|
32 | // |
---|
33 | // |
---|
34 | // File name: G4MuBremsstrahlungModel |
---|
35 | // |
---|
36 | // Author: Vladimir Ivanchenko on base of Laszlo Urban code |
---|
37 | // |
---|
38 | // Creation date: 24.06.2002 |
---|
39 | // |
---|
40 | // Modifications: |
---|
41 | // |
---|
42 | // 04-12-02 Change G4DynamicParticle constructor in PostStepDoIt (V.Ivanchenko) |
---|
43 | // 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko) |
---|
44 | // 24-01-03 Fix for compounds (V.Ivanchenko) |
---|
45 | // 27-01-03 Make models region aware (V.Ivanchenko) |
---|
46 | // 13-02-03 Add name (V.Ivanchenko) |
---|
47 | // 10-02-04 Add lowestKinEnergy (V.Ivanchenko) |
---|
48 | // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko) |
---|
49 | // 03-08-05 Angular correlations according to PRM (V.Ivantchenko) |
---|
50 | // 13-02-06 add ComputeCrossSectionPerAtom (mma) |
---|
51 | // 21-03-06 Fix problem of initialisation in case when cuts are not defined (VI) |
---|
52 | // 07-11-07 Improve sampling of final state (A.Bogdanov) |
---|
53 | // |
---|
54 | |
---|
55 | // |
---|
56 | // Class Description: |
---|
57 | // |
---|
58 | // |
---|
59 | // ------------------------------------------------------------------- |
---|
60 | // |
---|
61 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
62 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
63 | |
---|
64 | #include "G4MuBremsstrahlungModel.hh" |
---|
65 | #include "G4Gamma.hh" |
---|
66 | #include "G4MuonMinus.hh" |
---|
67 | #include "G4MuonPlus.hh" |
---|
68 | #include "Randomize.hh" |
---|
69 | #include "G4Material.hh" |
---|
70 | #include "G4Element.hh" |
---|
71 | #include "G4ElementVector.hh" |
---|
72 | #include "G4ProductionCutsTable.hh" |
---|
73 | #include "G4ParticleChangeForLoss.hh" |
---|
74 | |
---|
75 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
76 | |
---|
77 | // static members |
---|
78 | // |
---|
79 | G4double G4MuBremsstrahlungModel::zdat[]={1., 4., 13., 29., 92.}; |
---|
80 | G4double G4MuBremsstrahlungModel::adat[]={1.01, 9.01, 26.98, 63.55, 238.03}; |
---|
81 | G4double G4MuBremsstrahlungModel::tdat[]={1.e3, 1.e4, 1.e5, 1.e6, 1.e7, 1.e8, |
---|
82 | 1.e9, 1.e10}; |
---|
83 | |
---|
84 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
85 | |
---|
86 | using namespace std; |
---|
87 | |
---|
88 | G4MuBremsstrahlungModel::G4MuBremsstrahlungModel(const G4ParticleDefinition* p, |
---|
89 | const G4String& nam) |
---|
90 | : G4VEmModel(nam), |
---|
91 | particle(0), |
---|
92 | lowestKinEnergy(1.0*GeV), |
---|
93 | minThreshold(1.0*keV), |
---|
94 | nzdat(5), |
---|
95 | ntdat(8), |
---|
96 | NBIN(1000), |
---|
97 | cutFixed(0.98*keV), |
---|
98 | ignoreCut(false), |
---|
99 | samplingTablesAreFilled(false) |
---|
100 | { |
---|
101 | theGamma = G4Gamma::Gamma(); |
---|
102 | if(p) SetParticle(p); |
---|
103 | } |
---|
104 | |
---|
105 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
106 | |
---|
107 | G4MuBremsstrahlungModel::~G4MuBremsstrahlungModel() |
---|
108 | { |
---|
109 | size_t n = partialSumSigma.size(); |
---|
110 | if(n > 0) { |
---|
111 | for(size_t i=0; i<n; i++) { |
---|
112 | delete partialSumSigma[i]; |
---|
113 | } |
---|
114 | } |
---|
115 | } |
---|
116 | |
---|
117 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
118 | |
---|
119 | G4double G4MuBremsstrahlungModel::MinEnergyCut(const G4ParticleDefinition*, |
---|
120 | const G4MaterialCutsCouple*) |
---|
121 | { |
---|
122 | return minThreshold; |
---|
123 | } |
---|
124 | |
---|
125 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
126 | |
---|
127 | void G4MuBremsstrahlungModel::SetParticle(const G4ParticleDefinition* p) |
---|
128 | { |
---|
129 | if(!particle) { |
---|
130 | particle = p; |
---|
131 | mass = particle->GetPDGMass(); |
---|
132 | } |
---|
133 | } |
---|
134 | |
---|
135 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
136 | |
---|
137 | void G4MuBremsstrahlungModel::Initialise(const G4ParticleDefinition* p, |
---|
138 | const G4DataVector& cuts) |
---|
139 | { |
---|
140 | if(p) SetParticle(p); |
---|
141 | |
---|
142 | highKinEnergy = HighEnergyLimit(); |
---|
143 | |
---|
144 | G4double fixedEnergy = 0.5*highKinEnergy; |
---|
145 | |
---|
146 | const G4ProductionCutsTable* theCoupleTable= |
---|
147 | G4ProductionCutsTable::GetProductionCutsTable(); |
---|
148 | if(theCoupleTable) { |
---|
149 | G4int numOfCouples = theCoupleTable->GetTableSize(); |
---|
150 | |
---|
151 | G4int nn = partialSumSigma.size(); |
---|
152 | G4int nc = cuts.size(); |
---|
153 | if(nn > 0) { |
---|
154 | for (G4int ii=0; ii<nn; ii++){ |
---|
155 | G4DataVector* a=partialSumSigma[ii]; |
---|
156 | if ( a ) delete a; |
---|
157 | } |
---|
158 | partialSumSigma.clear(); |
---|
159 | } |
---|
160 | if (numOfCouples>0) { |
---|
161 | for (G4int i=0; i<numOfCouples; i++) { |
---|
162 | G4double cute = DBL_MAX; |
---|
163 | if(i < nc) cute = cuts[i]; |
---|
164 | if(cute < cutFixed || ignoreCut) cute = cutFixed; |
---|
165 | const G4MaterialCutsCouple* couple = |
---|
166 | theCoupleTable->GetMaterialCutsCouple(i); |
---|
167 | const G4Material* material = couple->GetMaterial(); |
---|
168 | G4DataVector* dv = ComputePartialSumSigma(material,fixedEnergy,cute); |
---|
169 | partialSumSigma.push_back(dv); |
---|
170 | } |
---|
171 | } |
---|
172 | } |
---|
173 | if(!samplingTablesAreFilled) MakeSamplingTables(); |
---|
174 | if(pParticleChange) |
---|
175 | fParticleChange = |
---|
176 | reinterpret_cast<G4ParticleChangeForLoss*>(pParticleChange); |
---|
177 | else |
---|
178 | fParticleChange = new G4ParticleChangeForLoss(); |
---|
179 | } |
---|
180 | |
---|
181 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
182 | |
---|
183 | G4double G4MuBremsstrahlungModel::ComputeDEDXPerVolume( |
---|
184 | const G4Material* material, |
---|
185 | const G4ParticleDefinition*, |
---|
186 | G4double kineticEnergy, |
---|
187 | G4double cutEnergy) |
---|
188 | { |
---|
189 | G4double dedx = 0.0; |
---|
190 | if (kineticEnergy <= lowestKinEnergy || ignoreCut) return dedx; |
---|
191 | |
---|
192 | G4double tmax = kineticEnergy; |
---|
193 | G4double cut = min(cutEnergy,tmax); |
---|
194 | if(cut < cutFixed) cut = cutFixed; |
---|
195 | |
---|
196 | const G4ElementVector* theElementVector = material->GetElementVector(); |
---|
197 | const G4double* theAtomicNumDensityVector = |
---|
198 | material->GetAtomicNumDensityVector(); |
---|
199 | |
---|
200 | // loop for elements in the material |
---|
201 | for (size_t i=0; i<material->GetNumberOfElements(); i++) { |
---|
202 | |
---|
203 | G4double Z = (*theElementVector)[i]->GetZ(); |
---|
204 | G4double A = (*theElementVector)[i]->GetA()/(g/mole) ; |
---|
205 | |
---|
206 | G4double loss = ComputMuBremLoss(Z, A, kineticEnergy, cut); |
---|
207 | |
---|
208 | dedx += loss*theAtomicNumDensityVector[i]; |
---|
209 | } |
---|
210 | if(dedx < 0.) dedx = 0.; |
---|
211 | return dedx; |
---|
212 | } |
---|
213 | |
---|
214 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
215 | |
---|
216 | G4double G4MuBremsstrahlungModel::ComputMuBremLoss(G4double Z, G4double A, |
---|
217 | G4double tkin, G4double cut) |
---|
218 | { |
---|
219 | G4double totalEnergy = mass + tkin; |
---|
220 | G4double ak1 = 0.05; |
---|
221 | G4int k2=5; |
---|
222 | G4double xgi[]={0.03377,0.16940,0.38069,0.61931,0.83060,0.96623}; |
---|
223 | G4double wgi[]={0.08566,0.18038,0.23396,0.23396,0.18038,0.08566}; |
---|
224 | G4double loss = 0.; |
---|
225 | |
---|
226 | G4double vcut = cut/totalEnergy; |
---|
227 | G4double vmax = tkin/totalEnergy; |
---|
228 | |
---|
229 | G4double aaa = 0.; |
---|
230 | G4double bbb = vcut; |
---|
231 | if(vcut>vmax) bbb=vmax ; |
---|
232 | G4int kkk = (G4int)((bbb-aaa)/ak1)+k2 ; |
---|
233 | G4double hhh=(bbb-aaa)/float(kkk) ; |
---|
234 | |
---|
235 | G4double aa = aaa; |
---|
236 | for(G4int l=0; l<kkk; l++) |
---|
237 | { |
---|
238 | for(G4int i=0; i<6; i++) |
---|
239 | { |
---|
240 | G4double ep = (aa + xgi[i]*hhh)*totalEnergy; |
---|
241 | loss += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, A, ep); |
---|
242 | } |
---|
243 | aa += hhh; |
---|
244 | } |
---|
245 | |
---|
246 | loss *=hhh*totalEnergy ; |
---|
247 | |
---|
248 | return loss; |
---|
249 | } |
---|
250 | |
---|
251 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
252 | |
---|
253 | G4double G4MuBremsstrahlungModel::ComputeMicroscopicCrossSection( |
---|
254 | G4double tkin, |
---|
255 | G4double Z, |
---|
256 | G4double A, |
---|
257 | G4double cut) |
---|
258 | { |
---|
259 | G4double totalEnergy = tkin + mass; |
---|
260 | G4double ak1 = 2.3; |
---|
261 | G4int k2 = 4; |
---|
262 | G4double xgi[]={0.03377,0.16940,0.38069,0.61931,0.83060,0.96623}; |
---|
263 | G4double wgi[]={0.08566,0.18038,0.23396,0.23396,0.18038,0.08566}; |
---|
264 | G4double cross = 0.; |
---|
265 | |
---|
266 | if(cut >= tkin) return cross; |
---|
267 | |
---|
268 | G4double vcut = cut/totalEnergy; |
---|
269 | G4double vmax = tkin/totalEnergy; |
---|
270 | |
---|
271 | G4double aaa = log(vcut); |
---|
272 | G4double bbb = log(vmax); |
---|
273 | G4int kkk = (G4int)((bbb-aaa)/ak1)+k2 ; |
---|
274 | G4double hhh = (bbb-aaa)/float(kkk); |
---|
275 | |
---|
276 | G4double aa = aaa; |
---|
277 | |
---|
278 | for(G4int l=0; l<kkk; l++) |
---|
279 | { |
---|
280 | for(G4int i=0; i<6; i++) |
---|
281 | { |
---|
282 | G4double ep = exp(aa + xgi[i]*hhh)*totalEnergy; |
---|
283 | cross += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, A, ep); |
---|
284 | } |
---|
285 | aa += hhh; |
---|
286 | } |
---|
287 | |
---|
288 | cross *=hhh; |
---|
289 | |
---|
290 | return cross; |
---|
291 | } |
---|
292 | |
---|
293 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
294 | |
---|
295 | G4double G4MuBremsstrahlungModel::ComputeDMicroscopicCrossSection( |
---|
296 | G4double tkin, |
---|
297 | G4double Z, |
---|
298 | G4double A, |
---|
299 | G4double gammaEnergy) |
---|
300 | // differential cross section |
---|
301 | { |
---|
302 | static const G4double sqrte=sqrt(exp(1.)) ; |
---|
303 | static const G4double bh=202.4,bh1=446.,btf=183.,btf1=1429. ; |
---|
304 | static const G4double rmass=mass/electron_mass_c2 ; |
---|
305 | static const G4double cc=classic_electr_radius/rmass ; |
---|
306 | static const G4double coeff= 16.*fine_structure_const*cc*cc/3. ; |
---|
307 | |
---|
308 | G4double dxsection = 0.; |
---|
309 | |
---|
310 | if( gammaEnergy > tkin) return dxsection ; |
---|
311 | |
---|
312 | G4double E = tkin + mass ; |
---|
313 | G4double v = gammaEnergy/E ; |
---|
314 | G4double delta = 0.5*mass*mass*v/(E-gammaEnergy) ; |
---|
315 | G4double rab0=delta*sqrte ; |
---|
316 | |
---|
317 | G4double z13 = exp(-log(Z)/3.) ; |
---|
318 | G4double dn = 1.54*exp(0.27*log(A)) ; |
---|
319 | |
---|
320 | G4double b,b1,dnstar ; |
---|
321 | |
---|
322 | if(Z<1.5) |
---|
323 | { |
---|
324 | b=bh; |
---|
325 | b1=bh1; |
---|
326 | dnstar=dn ; |
---|
327 | } |
---|
328 | else |
---|
329 | { |
---|
330 | b=btf; |
---|
331 | b1=btf1; |
---|
332 | dnstar = exp((1.-1./Z)*log(dn)) ; |
---|
333 | } |
---|
334 | |
---|
335 | // nucleus contribution logarithm |
---|
336 | G4double rab1=b*z13; |
---|
337 | G4double fn=log(rab1/(dnstar*(electron_mass_c2+rab0*rab1))* |
---|
338 | (mass+delta*(dnstar*sqrte-2.))) ; |
---|
339 | if(fn <0.) fn = 0. ; |
---|
340 | // electron contribution logarithm |
---|
341 | G4double epmax1=E/(1.+0.5*mass*rmass/E) ; |
---|
342 | G4double fe=0.; |
---|
343 | if(gammaEnergy<epmax1) |
---|
344 | { |
---|
345 | G4double rab2=b1*z13*z13 ; |
---|
346 | fe=log(rab2*mass/((1.+delta*rmass/(electron_mass_c2*sqrte))* |
---|
347 | (electron_mass_c2+rab0*rab2))) ; |
---|
348 | if(fe<0.) fe=0. ; |
---|
349 | } |
---|
350 | |
---|
351 | dxsection = coeff*(1.-v*(1. - 0.75*v))*Z*(fn*Z + fe)/gammaEnergy; |
---|
352 | |
---|
353 | return dxsection; |
---|
354 | } |
---|
355 | |
---|
356 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
357 | |
---|
358 | G4double G4MuBremsstrahlungModel::ComputeCrossSectionPerAtom( |
---|
359 | const G4ParticleDefinition*, |
---|
360 | G4double kineticEnergy, |
---|
361 | G4double Z, G4double A, |
---|
362 | G4double cutEnergy, |
---|
363 | G4double) |
---|
364 | { |
---|
365 | G4double cut = min(cutEnergy, kineticEnergy); |
---|
366 | if(cut < cutFixed || ignoreCut) cut = cutFixed; |
---|
367 | G4double cross = |
---|
368 | ComputeMicroscopicCrossSection (kineticEnergy, Z, A/(g/mole), cut); |
---|
369 | return cross; |
---|
370 | } |
---|
371 | |
---|
372 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
373 | |
---|
374 | G4double G4MuBremsstrahlungModel::CrossSectionPerVolume( |
---|
375 | const G4Material* material, |
---|
376 | const G4ParticleDefinition*, |
---|
377 | G4double kineticEnergy, |
---|
378 | G4double cutEnergy, |
---|
379 | G4double maxEnergy) |
---|
380 | { |
---|
381 | G4double cross = 0.0; |
---|
382 | if (cutEnergy >= maxEnergy || kineticEnergy <= lowestKinEnergy) return cross; |
---|
383 | |
---|
384 | G4double tmax = min(maxEnergy, kineticEnergy); |
---|
385 | G4double cut = min(cutEnergy, tmax); |
---|
386 | if(cut < cutFixed || ignoreCut) cut = cutFixed; |
---|
387 | |
---|
388 | const G4ElementVector* theElementVector = material->GetElementVector(); |
---|
389 | const G4double* theAtomNumDensityVector = |
---|
390 | material->GetAtomicNumDensityVector(); |
---|
391 | |
---|
392 | for (size_t i=0; i<material->GetNumberOfElements(); i++) { |
---|
393 | |
---|
394 | G4double Z = (*theElementVector)[i]->GetZ(); |
---|
395 | G4double A = (*theElementVector)[i]->GetA()/(g/mole); |
---|
396 | |
---|
397 | G4double cr = ComputeMicroscopicCrossSection(kineticEnergy, Z, A, cut); |
---|
398 | |
---|
399 | if(tmax < kineticEnergy) { |
---|
400 | cr -= ComputeMicroscopicCrossSection(kineticEnergy, Z, A, tmax); |
---|
401 | } |
---|
402 | cross += theAtomNumDensityVector[i] * cr; |
---|
403 | } |
---|
404 | |
---|
405 | return cross; |
---|
406 | } |
---|
407 | |
---|
408 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
409 | |
---|
410 | G4DataVector* G4MuBremsstrahlungModel::ComputePartialSumSigma( |
---|
411 | const G4Material* material, |
---|
412 | G4double kineticEnergy, |
---|
413 | G4double cut) |
---|
414 | |
---|
415 | // Build the table of cross section per element. The table is built for MATERIAL |
---|
416 | // This table is used by DoIt to select randomly an element in the material. |
---|
417 | { |
---|
418 | G4int nElements = material->GetNumberOfElements(); |
---|
419 | const G4ElementVector* theElementVector = material->GetElementVector(); |
---|
420 | const G4double* theAtomNumDensityVector = |
---|
421 | material->GetAtomicNumDensityVector(); |
---|
422 | |
---|
423 | G4DataVector* dv = new G4DataVector(); |
---|
424 | |
---|
425 | G4double cross = 0.0; |
---|
426 | |
---|
427 | for (G4int i=0; i<nElements; i++ ) { |
---|
428 | |
---|
429 | G4double Z = (*theElementVector)[i]->GetZ(); |
---|
430 | G4double A = (*theElementVector)[i]->GetA()/(g/mole) ; |
---|
431 | cross += theAtomNumDensityVector[i] |
---|
432 | * ComputeMicroscopicCrossSection(kineticEnergy, Z, A, cut); |
---|
433 | dv->push_back(cross); |
---|
434 | } |
---|
435 | return dv; |
---|
436 | } |
---|
437 | |
---|
438 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
439 | |
---|
440 | void G4MuBremsstrahlungModel::MakeSamplingTables() |
---|
441 | { |
---|
442 | |
---|
443 | G4double AtomicNumber,AtomicWeight,KineticEnergy, |
---|
444 | TotalEnergy,Maxep; |
---|
445 | |
---|
446 | for (G4int iz=0; iz<nzdat; iz++) |
---|
447 | { |
---|
448 | AtomicNumber = zdat[iz]; |
---|
449 | AtomicWeight = adat[iz]*g/mole ; |
---|
450 | |
---|
451 | for (G4int it=0; it<ntdat; it++) |
---|
452 | { |
---|
453 | KineticEnergy = tdat[it]; |
---|
454 | TotalEnergy = KineticEnergy + mass; |
---|
455 | Maxep = KineticEnergy ; |
---|
456 | |
---|
457 | G4double CrossSection = 0.0 ; |
---|
458 | |
---|
459 | // calculate the differential cross section |
---|
460 | // numerical integration in |
---|
461 | // log ............... |
---|
462 | G4double c = log(Maxep/cutFixed) ; |
---|
463 | G4double ymin = -5. ; |
---|
464 | G4double ymax = 0. ; |
---|
465 | G4double dy = (ymax-ymin)/NBIN ; |
---|
466 | |
---|
467 | G4double y = ymin - 0.5*dy ; |
---|
468 | G4double yy = ymin - dy ; |
---|
469 | G4double x = exp(y); |
---|
470 | G4double fac = exp(dy); |
---|
471 | G4double dx = exp(yy)*(fac - 1.0); |
---|
472 | |
---|
473 | for (G4int i=0 ; i<NBIN; i++) |
---|
474 | { |
---|
475 | y += dy ; |
---|
476 | x *= fac; |
---|
477 | dx*= fac; |
---|
478 | G4double ep = cutFixed*exp(c*x) ; |
---|
479 | |
---|
480 | CrossSection += ep*dx*ComputeDMicroscopicCrossSection( |
---|
481 | KineticEnergy,AtomicNumber, |
---|
482 | AtomicWeight,ep) ; |
---|
483 | ya[i]=y ; |
---|
484 | proba[iz][it][i] = CrossSection ; |
---|
485 | |
---|
486 | } |
---|
487 | |
---|
488 | proba[iz][it][NBIN] = CrossSection ; |
---|
489 | ya[NBIN] = 0. ; // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
490 | |
---|
491 | if(CrossSection > 0.) |
---|
492 | { |
---|
493 | for(G4int ib=0; ib<=NBIN; ib++) |
---|
494 | { |
---|
495 | proba[iz][it][ib] /= CrossSection ; |
---|
496 | } |
---|
497 | } |
---|
498 | } |
---|
499 | } |
---|
500 | samplingTablesAreFilled = true; |
---|
501 | } |
---|
502 | |
---|
503 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
504 | |
---|
505 | void G4MuBremsstrahlungModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp, |
---|
506 | const G4MaterialCutsCouple* couple, |
---|
507 | const G4DynamicParticle* dp, |
---|
508 | G4double minEnergy, |
---|
509 | G4double maxEnergy) |
---|
510 | { |
---|
511 | G4double kineticEnergy = dp->GetKineticEnergy(); |
---|
512 | // check against insufficient energy |
---|
513 | G4double tmax = min(kineticEnergy, maxEnergy); |
---|
514 | G4double tmin = min(kineticEnergy, minEnergy); |
---|
515 | if(tmin < cutFixed || ignoreCut) tmin = cutFixed; |
---|
516 | if(tmin >= tmax) return; |
---|
517 | |
---|
518 | // ===== the begining of a new code ====== |
---|
519 | // ===== sampling of energy transfer ====== |
---|
520 | |
---|
521 | G4ParticleMomentum partDirection = dp->GetMomentumDirection(); |
---|
522 | |
---|
523 | // select randomly one element constituing the material |
---|
524 | const G4Element* anElement = SelectRandomAtom(couple); |
---|
525 | |
---|
526 | G4double totalEnergy = kineticEnergy + mass; |
---|
527 | G4double totalMomentum = sqrt(kineticEnergy*(kineticEnergy + 2.0*mass)); |
---|
528 | |
---|
529 | G4double AtomicNumber = anElement->GetZ(); |
---|
530 | G4double AtomicWeight = anElement->GetA()/(g/mole); |
---|
531 | |
---|
532 | G4double func1 = tmin*ComputeDMicroscopicCrossSection( |
---|
533 | kineticEnergy,AtomicNumber, |
---|
534 | AtomicWeight,tmin); |
---|
535 | |
---|
536 | G4double lnepksi, epksi; |
---|
537 | G4double func2; |
---|
538 | G4double ksi2; |
---|
539 | |
---|
540 | do { |
---|
541 | lnepksi = log(tmin) + G4UniformRand()*log(kineticEnergy/tmin); |
---|
542 | epksi = exp(lnepksi); |
---|
543 | func2 = epksi*ComputeDMicroscopicCrossSection( |
---|
544 | kineticEnergy,AtomicNumber, |
---|
545 | AtomicWeight,epksi); |
---|
546 | ksi2 = G4UniformRand(); |
---|
547 | |
---|
548 | } while(func2/func1 < ksi2); |
---|
549 | |
---|
550 | // ===== the end of a new code ===== |
---|
551 | |
---|
552 | // create G4DynamicParticle object for the Gamma |
---|
553 | G4double gEnergy = epksi; |
---|
554 | |
---|
555 | // sample angle |
---|
556 | G4double gam = totalEnergy/mass; |
---|
557 | G4double rmax = gam*min(1.0, totalEnergy/gEnergy - 1.0); |
---|
558 | rmax *= rmax; |
---|
559 | G4double x = G4UniformRand()*rmax/(1.0 + rmax); |
---|
560 | |
---|
561 | G4double theta = sqrt(x/(1.0 - x))/gam; |
---|
562 | G4double sint = sin(theta); |
---|
563 | G4double phi = twopi * G4UniformRand() ; |
---|
564 | G4double dirx = sint*cos(phi), diry = sint*sin(phi), dirz = cos(theta) ; |
---|
565 | |
---|
566 | G4ThreeVector gDirection(dirx, diry, dirz); |
---|
567 | gDirection.rotateUz(partDirection); |
---|
568 | |
---|
569 | partDirection *= totalMomentum; |
---|
570 | partDirection -= gEnergy*gDirection; |
---|
571 | partDirection = partDirection.unit(); |
---|
572 | |
---|
573 | // primary change |
---|
574 | kineticEnergy -= gEnergy; |
---|
575 | fParticleChange->SetProposedKineticEnergy(kineticEnergy); |
---|
576 | fParticleChange->SetProposedMomentumDirection(partDirection); |
---|
577 | |
---|
578 | // save secondary |
---|
579 | G4DynamicParticle* aGamma = new G4DynamicParticle(theGamma,gDirection,gEnergy); |
---|
580 | vdp->push_back(aGamma); |
---|
581 | } |
---|
582 | |
---|
583 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|
584 | |
---|
585 | const G4Element* G4MuBremsstrahlungModel::SelectRandomAtom( |
---|
586 | const G4MaterialCutsCouple* couple) const |
---|
587 | { |
---|
588 | // select randomly 1 element within the material |
---|
589 | |
---|
590 | const G4Material* material = couple->GetMaterial(); |
---|
591 | G4int nElements = material->GetNumberOfElements(); |
---|
592 | const G4ElementVector* theElementVector = material->GetElementVector(); |
---|
593 | if(1 == nElements) return (*theElementVector)[0]; |
---|
594 | else if(1 > nElements) return 0; |
---|
595 | |
---|
596 | G4DataVector* dv = partialSumSigma[couple->GetIndex()]; |
---|
597 | G4double rval = G4UniformRand()*((*dv)[nElements-1]); |
---|
598 | for (G4int i=0; i<nElements; i++) { |
---|
599 | if (rval <= (*dv)[i]) return (*theElementVector)[i]; |
---|
600 | } |
---|
601 | return (*theElementVector)[nElements-1]; |
---|
602 | } |
---|
603 | |
---|
604 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... |
---|