source: trunk/source/processes/electromagnetic/muons/src/G4MuBremsstrahlungModel.cc@ 1006

Last change on this file since 1006 was 1005, checked in by garnier, 17 years ago

fichier oublies

File size: 15.4 KB
Line 
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.33 2009/02/20 14:48:16 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-02-ref-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.Ivanchenko)
49// 03-08-05 Angular correlations according to PRM (V.Ivanchenko)
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// 28-02-08 Use precomputed Z^1/3 and Log(A) (V.Ivanchenko)
54//
55
56//
57// Class Description:
58//
59//
60// -------------------------------------------------------------------
61//
62//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
63//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
64
65#include "G4MuBremsstrahlungModel.hh"
66#include "G4Gamma.hh"
67#include "G4MuonMinus.hh"
68#include "G4MuonPlus.hh"
69#include "Randomize.hh"
70#include "G4Material.hh"
71#include "G4Element.hh"
72#include "G4ElementVector.hh"
73#include "G4ProductionCutsTable.hh"
74#include "G4ParticleChangeForLoss.hh"
75
76//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
77//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
78
79using namespace std;
80
81G4MuBremsstrahlungModel::G4MuBremsstrahlungModel(const G4ParticleDefinition* p,
82 const G4String& nam)
83 : G4VEmModel(nam),
84 particle(0),
85 sqrte(sqrt(exp(1.))),
86 bh(202.4),
87 bh1(446.),
88 btf(183.),
89 btf1(1429.),
90 fParticleChange(0),
91 lowestKinEnergy(1.0*GeV),
92 minThreshold(1.0*keV)
93{
94 theGamma = G4Gamma::Gamma();
95 nist = G4NistManager::Instance();
96 if(p) SetParticle(p);
97}
98
99//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
100
101G4MuBremsstrahlungModel::~G4MuBremsstrahlungModel()
102{
103 size_t n = partialSumSigma.size();
104 if(n > 0) {
105 for(size_t i=0; i<n; i++) {
106 delete partialSumSigma[i];
107 }
108 }
109}
110
111//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
112
113G4double G4MuBremsstrahlungModel::MinEnergyCut(const G4ParticleDefinition*,
114 const G4MaterialCutsCouple*)
115{
116 return minThreshold;
117}
118
119//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
120
121void G4MuBremsstrahlungModel::Initialise(const G4ParticleDefinition* p,
122 const G4DataVector& cuts)
123{
124 if(p) SetParticle(p);
125
126 highKinEnergy = HighEnergyLimit();
127
128 // partial cross section is computed for fixed energy
129 G4double fixedEnergy = 0.5*highKinEnergy;
130
131 const G4ProductionCutsTable* theCoupleTable=
132 G4ProductionCutsTable::GetProductionCutsTable();
133 if(theCoupleTable) {
134 G4int numOfCouples = theCoupleTable->GetTableSize();
135
136 // clear old data
137 G4int nn = partialSumSigma.size();
138 G4int nc = cuts.size();
139 if(nn > 0) {
140 for (G4int ii=0; ii<nn; ii++){
141 G4DataVector* a = partialSumSigma[ii];
142 if ( a ) delete a;
143 }
144 partialSumSigma.clear();
145 }
146 // fill new data
147 if (numOfCouples>0) {
148 for (G4int i=0; i<numOfCouples; i++) {
149 G4double cute = DBL_MAX;
150
151 // protection for usage with extrapolator
152 if(i < nc) cute = cuts[i];
153
154 const G4MaterialCutsCouple* couple =
155 theCoupleTable->GetMaterialCutsCouple(i);
156 const G4Material* material = couple->GetMaterial();
157 G4DataVector* dv = ComputePartialSumSigma(material,fixedEnergy,cute);
158 partialSumSigma.push_back(dv);
159 }
160 }
161 }
162
163 // define pointer to G4ParticleChange
164 if(!fParticleChange) {
165 if(pParticleChange)
166 fParticleChange =
167 reinterpret_cast<G4ParticleChangeForLoss*>(pParticleChange);
168 else
169 fParticleChange = new G4ParticleChangeForLoss();
170 }
171}
172
173//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
174
175G4double G4MuBremsstrahlungModel::ComputeDEDXPerVolume(
176 const G4Material* material,
177 const G4ParticleDefinition*,
178 G4double kineticEnergy,
179 G4double cutEnergy)
180{
181 G4double dedx = 0.0;
182 if (kineticEnergy <= lowestKinEnergy) return dedx;
183
184 G4double tmax = kineticEnergy;
185 G4double cut = std::min(cutEnergy,tmax);
186 if(cut < minThreshold) cut = minThreshold;
187
188 const G4ElementVector* theElementVector = material->GetElementVector();
189 const G4double* theAtomicNumDensityVector =
190 material->GetAtomicNumDensityVector();
191
192 // loop for elements in the material
193 for (size_t i=0; i<material->GetNumberOfElements(); i++) {
194
195 G4double loss =
196 ComputMuBremLoss((*theElementVector)[i]->GetZ(), kineticEnergy, cut);
197
198 dedx += loss*theAtomicNumDensityVector[i];
199 }
200 // G4cout << "BR e= " << kineticEnergy << " dedx= " << dedx << G4endl;
201 if(dedx < 0.) dedx = 0.;
202 return dedx;
203}
204
205//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
206
207G4double G4MuBremsstrahlungModel::ComputMuBremLoss(G4double Z,
208 G4double tkin, G4double cut)
209{
210 G4double totalEnergy = mass + tkin;
211 G4double ak1 = 0.05;
212 G4int k2=5;
213 G4double xgi[]={0.03377,0.16940,0.38069,0.61931,0.83060,0.96623};
214 G4double wgi[]={0.08566,0.18038,0.23396,0.23396,0.18038,0.08566};
215 G4double loss = 0.;
216
217 G4double vcut = cut/totalEnergy;
218 G4double vmax = tkin/totalEnergy;
219
220 G4double aaa = 0.;
221 G4double bbb = vcut;
222 if(vcut>vmax) bbb=vmax ;
223 G4int kkk = (G4int)((bbb-aaa)/ak1)+k2 ;
224 G4double hhh=(bbb-aaa)/float(kkk) ;
225
226 G4double aa = aaa;
227 for(G4int l=0; l<kkk; l++)
228 {
229 for(G4int i=0; i<6; i++)
230 {
231 G4double ep = (aa + xgi[i]*hhh)*totalEnergy;
232 loss += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
233 }
234 aa += hhh;
235 }
236
237 loss *=hhh*totalEnergy ;
238
239 return loss;
240}
241
242//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
243
244G4double G4MuBremsstrahlungModel::ComputeMicroscopicCrossSection(
245 G4double tkin,
246 G4double Z,
247 G4double cut)
248{
249 G4double totalEnergy = tkin + mass;
250 G4double ak1 = 2.3;
251 G4int k2 = 4;
252 G4double xgi[]={0.03377,0.16940,0.38069,0.61931,0.83060,0.96623};
253 G4double wgi[]={0.08566,0.18038,0.23396,0.23396,0.18038,0.08566};
254 G4double cross = 0.;
255
256 if(cut >= tkin) return cross;
257
258 G4double vcut = cut/totalEnergy;
259 G4double vmax = tkin/totalEnergy;
260
261 G4double aaa = log(vcut);
262 G4double bbb = log(vmax);
263 G4int kkk = (G4int)((bbb-aaa)/ak1)+k2 ;
264 G4double hhh = (bbb-aaa)/G4double(kkk);
265
266 G4double aa = aaa;
267
268 for(G4int l=0; l<kkk; l++)
269 {
270 for(G4int i=0; i<6; i++)
271 {
272 G4double ep = exp(aa + xgi[i]*hhh)*totalEnergy;
273 cross += ep*wgi[i]*ComputeDMicroscopicCrossSection(tkin, Z, ep);
274 }
275 aa += hhh;
276 }
277
278 cross *=hhh;
279
280 //G4cout << "BR e= " << tkin<< " cross= " << cross/barn << G4endl;
281
282 return cross;
283}
284
285//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
286
287G4double G4MuBremsstrahlungModel::ComputeDMicroscopicCrossSection(
288 G4double tkin,
289 G4double Z,
290 G4double gammaEnergy)
291// differential cross section
292{
293 G4double dxsection = 0.;
294
295 if( gammaEnergy > tkin) return dxsection ;
296
297 G4double E = tkin + mass ;
298 G4double v = gammaEnergy/E ;
299 G4double delta = 0.5*mass*mass*v/(E-gammaEnergy) ;
300 G4double rab0=delta*sqrte ;
301
302 G4int iz = G4int(Z);
303 if(iz < 1) iz = 1;
304
305 G4double z13 = 1.0/nist->GetZ13(iz);
306 G4double dn = 1.54*nist->GetA27(iz);
307
308 G4double b,b1,dnstar ;
309
310 if(1 == iz)
311 {
312 b = bh;
313 b1 = bh1;
314 dnstar = dn;
315 }
316 else
317 {
318 b = btf;
319 b1 = btf1;
320 dnstar = dn/std::pow(dn, 1./Z);
321 }
322
323 // nucleus contribution logarithm
324 G4double rab1=b*z13;
325 G4double fn=log(rab1/(dnstar*(electron_mass_c2+rab0*rab1))*
326 (mass+delta*(dnstar*sqrte-2.))) ;
327 if(fn <0.) fn = 0. ;
328 // electron contribution logarithm
329 G4double epmax1=E/(1.+0.5*mass*rmass/E) ;
330 G4double fe=0.;
331 if(gammaEnergy<epmax1)
332 {
333 G4double rab2=b1*z13*z13 ;
334 fe=log(rab2*mass/((1.+delta*rmass/(electron_mass_c2*sqrte))*
335 (electron_mass_c2+rab0*rab2))) ;
336 if(fe<0.) fe=0. ;
337 }
338
339 dxsection = coeff*(1.-v*(1. - 0.75*v))*Z*(fn*Z + fe)/gammaEnergy;
340
341 return dxsection;
342}
343
344//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
345
346G4double G4MuBremsstrahlungModel::ComputeCrossSectionPerAtom(
347 const G4ParticleDefinition*,
348 G4double kineticEnergy,
349 G4double Z, G4double,
350 G4double cutEnergy,
351 G4double maxEnergy)
352{
353 G4double cross = 0.0;
354 if (kineticEnergy <= lowestKinEnergy) return cross;
355 G4double tmax = std::min(maxEnergy, kineticEnergy);
356 G4double cut = std::min(cutEnergy, kineticEnergy);
357 if(cut < minThreshold) cut = minThreshold;
358 if (cut >= tmax) return cross;
359
360 cross = ComputeMicroscopicCrossSection (kineticEnergy, Z, cut);
361 if(tmax < kineticEnergy) {
362 cross -= ComputeMicroscopicCrossSection(kineticEnergy, Z, tmax);
363 }
364 return cross;
365}
366
367//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
368
369G4DataVector* G4MuBremsstrahlungModel::ComputePartialSumSigma(
370 const G4Material* material,
371 G4double kineticEnergy,
372 G4double cut)
373
374// Build the table of cross section per element.
375// The table is built for material
376// This table is used to select randomly an element in the material.
377{
378 G4int nElements = material->GetNumberOfElements();
379 const G4ElementVector* theElementVector = material->GetElementVector();
380 const G4double* theAtomNumDensityVector =
381 material->GetAtomicNumDensityVector();
382
383 G4DataVector* dv = new G4DataVector();
384
385 G4double cross = 0.0;
386
387 for (G4int i=0; i<nElements; i++ ) {
388 cross += theAtomNumDensityVector[i]
389 * ComputeMicroscopicCrossSection(kineticEnergy,
390 (*theElementVector)[i]->GetZ(), cut);
391 dv->push_back(cross);
392 }
393 return dv;
394}
395
396//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
397
398void G4MuBremsstrahlungModel::SampleSecondaries(
399 std::vector<G4DynamicParticle*>* vdp,
400 const G4MaterialCutsCouple* couple,
401 const G4DynamicParticle* dp,
402 G4double minEnergy,
403 G4double maxEnergy)
404{
405 G4double kineticEnergy = dp->GetKineticEnergy();
406 // check against insufficient energy
407 G4double tmax = std::min(kineticEnergy, maxEnergy);
408 G4double tmin = std::min(kineticEnergy, minEnergy);
409 if(tmin < minThreshold) tmin = minThreshold;
410 if(tmin >= tmax) return;
411
412 // ===== sampling of energy transfer ======
413
414 G4ParticleMomentum partDirection = dp->GetMomentumDirection();
415
416 // select randomly one element constituing the material
417 const G4Element* anElement = SelectRandomAtom(couple);
418 G4double Z = anElement->GetZ();
419
420 G4double totalEnergy = kineticEnergy + mass;
421 G4double totalMomentum = sqrt(kineticEnergy*(kineticEnergy + 2.0*mass));
422
423 G4double func1 = tmin*
424 ComputeDMicroscopicCrossSection(kineticEnergy,Z,tmin);
425
426 G4double lnepksi, epksi;
427 G4double func2;
428
429 do {
430 lnepksi = log(tmin) + G4UniformRand()*log(kineticEnergy/tmin);
431 epksi = exp(lnepksi);
432 func2 = epksi*ComputeDMicroscopicCrossSection(kineticEnergy,Z,epksi);
433
434 } while(func2 < func1*G4UniformRand());
435
436 G4double gEnergy = epksi;
437
438 // ===== sample angle =====
439
440 G4double gam = totalEnergy/mass;
441 G4double rmax = gam*std::min(1.0, totalEnergy/gEnergy - 1.0);
442 G4double rmax2= rmax*rmax;
443 G4double x = G4UniformRand()*rmax2/(1.0 + rmax2);
444
445 G4double theta = sqrt(x/(1.0 - x))/gam;
446 G4double sint = sin(theta);
447 G4double phi = twopi * G4UniformRand() ;
448 G4double dirx = sint*cos(phi), diry = sint*sin(phi), dirz = cos(theta) ;
449
450 G4ThreeVector gDirection(dirx, diry, dirz);
451 gDirection.rotateUz(partDirection);
452
453 partDirection *= totalMomentum;
454 partDirection -= gEnergy*gDirection;
455 partDirection = partDirection.unit();
456
457 // primary change
458 kineticEnergy -= gEnergy;
459 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
460 fParticleChange->SetProposedMomentumDirection(partDirection);
461
462 // save secondary
463 G4DynamicParticle* aGamma =
464 new G4DynamicParticle(theGamma,gDirection,gEnergy);
465 vdp->push_back(aGamma);
466}
467
468//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
469
470const G4Element* G4MuBremsstrahlungModel::SelectRandomAtom(
471 const G4MaterialCutsCouple* couple) const
472{
473 // select randomly 1 element within the material
474
475 const G4Material* material = couple->GetMaterial();
476 G4int nElements = material->GetNumberOfElements();
477 const G4ElementVector* theElementVector = material->GetElementVector();
478 if(1 == nElements) return (*theElementVector)[0];
479 else if(1 > nElements) return 0;
480
481 G4DataVector* dv = partialSumSigma[couple->GetIndex()];
482 G4double rval = G4UniformRand()*((*dv)[nElements-1]);
483 for (G4int i=0; i<nElements; i++) {
484 if (rval <= (*dv)[i]) return (*theElementVector)[i];
485 }
486 return (*theElementVector)[nElements-1];
487}
488
489//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracBrowser for help on using the repository browser.