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

Last change on this file since 1047 was 1007, checked in by garnier, 15 years ago

update to geant4.9.2

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