source: trunk/source/processes/electromagnetic/highenergy/src/G4GammaConversionToMuons.cc@ 830

Last change on this file since 830 was 819, checked in by garnier, 17 years ago

import all except CVS

File size: 14.9 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//
27// $Id: G4GammaConversionToMuons.cc,v 1.4 2006/06/29 19:32:40 gunter Exp $
28// GEANT4 tag $Name: geant4-09-01-patch-02 $
29//
30// ------------ G4GammaConversionToMuons physics process ------
31// by H.Burkhardt, S. Kelner and R. Kokoulin, April 2002
32//
33//
34// 07-08-02: missprint in OR condition in DoIt : f1<0 || f1>f1_max ..etc ...
35// 25-10-04: migrade to new interfaces of ParticleChange (vi)
36// ---------------------------------------------------------------------------
37
38#include "G4GammaConversionToMuons.hh"
39#include "G4EnergyLossTables.hh"
40#include "G4UnitsTable.hh"
41#include "G4MuonPlus.hh"
42#include "G4MuonMinus.hh"
43
44//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
45
46using namespace std;
47
48G4GammaConversionToMuons::G4GammaConversionToMuons(const G4String& processName,
49 G4ProcessType type):G4VDiscreteProcess (processName, type),
50 LowestEnergyLimit (4*G4MuonPlus::MuonPlus()->GetPDGMass()), // 4*Mmuon
51 HighestEnergyLimit(1e21*eV), // ok to 1e21eV=1e12GeV, then LPM suppression
52 CrossSecFactor(1.)
53{ }
54
55//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
56
57// destructor
58
59G4GammaConversionToMuons::~G4GammaConversionToMuons() // (empty) destructor
60{ }
61
62//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
63
64G4bool G4GammaConversionToMuons::IsApplicable(
65 const G4ParticleDefinition& particle)
66{
67 return ( &particle == G4Gamma::Gamma() );
68}
69
70//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
71
72void G4GammaConversionToMuons::BuildPhysicsTable(const G4ParticleDefinition&)
73// Build cross section and mean free path tables
74{ //here no tables, just calling PrintInfoDefinition
75 PrintInfoDefinition();
76}
77
78//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
79
80G4double G4GammaConversionToMuons::GetMeanFreePath(const G4Track& aTrack,
81 G4double, G4ForceCondition*)
82
83// returns the photon mean free path in GEANT4 internal units
84// (MeanFreePath is a private member of the class)
85
86{
87 const G4DynamicParticle* aDynamicGamma = aTrack.GetDynamicParticle();
88 G4double GammaEnergy = aDynamicGamma->GetKineticEnergy();
89 G4Material* aMaterial = aTrack.GetMaterial();
90
91 if (GammaEnergy < LowestEnergyLimit)
92 MeanFreePath = DBL_MAX;
93 else
94 MeanFreePath = ComputeMeanFreePath(GammaEnergy,aMaterial);
95
96 return MeanFreePath;
97}
98
99//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
100
101G4double G4GammaConversionToMuons::ComputeMeanFreePath(G4double GammaEnergy,
102 G4Material* aMaterial)
103
104// computes and returns the photon mean free path in GEANT4 internal units
105{
106 const G4ElementVector* theElementVector = aMaterial->GetElementVector();
107 const G4double* NbOfAtomsPerVolume = aMaterial->GetVecNbOfAtomsPerVolume();
108
109 G4double SIGMA = 0 ;
110
111 for ( size_t i=0 ; i < aMaterial->GetNumberOfElements() ; i++ )
112 {
113 G4double AtomicZ = (*theElementVector)[i]->GetZ();
114 G4double AtomicA = (*theElementVector)[i]->GetA()/(g/mole);
115 SIGMA += NbOfAtomsPerVolume[i] *
116 ComputeCrossSectionPerAtom(GammaEnergy,AtomicZ,AtomicA);
117 }
118 return SIGMA > DBL_MIN ? 1./SIGMA : DBL_MAX;
119}
120
121//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
122
123G4double G4GammaConversionToMuons::GetCrossSectionPerAtom(
124 const G4DynamicParticle* aDynamicGamma,
125 G4Element* anElement)
126
127// gives the total cross section per atom in GEANT4 internal units
128{
129 G4double GammaEnergy = aDynamicGamma->GetKineticEnergy();
130 G4double AtomicZ = anElement->GetZ();
131 G4double AtomicA = anElement->GetA()/(g/mole);
132 G4double crossSection =
133 ComputeCrossSectionPerAtom(GammaEnergy,AtomicZ,AtomicA);
134 return crossSection;
135}
136
137//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
138
139G4double G4GammaConversionToMuons::ComputeCrossSectionPerAtom(
140 G4double Egam, G4double Z, G4double A)
141
142// Calculates the microscopic cross section in GEANT4 internal units.
143// Total cross section parametrisation from H.Burkhardt
144// It gives a good description at any energy (from 0 to 10**21 eV)
145{ static const G4double Mmuon=G4MuonPlus::MuonPlus()->GetPDGMass();
146 static const G4double Mele=electron_mass_c2;
147 static const G4double Rc=elm_coupling/Mmuon; // classical particle radius
148 static const G4double sqrte=sqrt(exp(1.));
149 static const G4double PowSat=-0.88;
150
151 static G4double CrossSection = 0.0 ;
152
153 if ( A < 1. ) return 0;
154 if ( Egam < 4*Mmuon ) return 0 ; // below threshold return 0
155
156 static G4double EgamLast=0,Zlast=0,PowThres,Ecor,B,Dn,Zthird,Winfty,WMedAppr,
157 Wsatur,sigfac;
158
159 if(Zlast==Z && Egam==EgamLast) return CrossSection; // already calculated
160 EgamLast=Egam;
161
162 if(Zlast!=Z) // new element
163 { Zlast=Z;
164 if(Z==1) // special case of Hydrogen
165 { B=202.4;
166 Dn=1.49;
167 }
168 else
169 { B=183.;
170 Dn=1.54*pow(A,0.27);
171 }
172 Zthird=pow(Z,-1./3.); // Z**(-1/3)
173 Winfty=B*Zthird*Mmuon/(Dn*Mele);
174 WMedAppr=1./(4.*Dn*sqrte*Mmuon);
175 Wsatur=Winfty/WMedAppr;
176 sigfac=4.*fine_structure_const*Z*Z*Rc*Rc;
177 PowThres=1.479+0.00799*Dn;
178 Ecor=-18.+4347./(B*Zthird);
179 }
180 G4double CorFuc=1.+.04*log(1.+Ecor/Egam);
181 G4double Eg=pow(1.-4.*Mmuon/Egam,PowThres)*pow( pow(Wsatur,PowSat)+
182 pow(Egam,PowSat),1./PowSat); // threshold and saturation
183 CrossSection=7./9.*sigfac*log(1.+WMedAppr*CorFuc*Eg);
184 CrossSection*=CrossSecFactor; // increase the CrossSection by (by default 1)
185 return CrossSection;
186}
187
188//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
189
190void G4GammaConversionToMuons::SetCrossSecFactor(G4double fac)
191// Set the factor to artificially increase the cross section
192{ CrossSecFactor=fac;
193 G4cout << "The cross section for GammaConversionToMuons is artificially "
194 << "increased by the CrossSecFactor=" << CrossSecFactor << G4endl;
195}
196
197//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
198
199G4VParticleChange* G4GammaConversionToMuons::PostStepDoIt(
200 const G4Track& aTrack,
201 const G4Step& aStep)
202//
203// generation of gamma->mu+mu-
204//
205{
206 aParticleChange.Initialize(aTrack);
207 G4Material* aMaterial = aTrack.GetMaterial();
208
209 static const G4double Mmuon=G4MuonPlus::MuonPlus()->GetPDGMass();
210 static const G4double Mele=electron_mass_c2;
211 static const G4double sqrte=sqrt(exp(1.));
212
213 // current Gamma energy and direction, return if energy too low
214 const G4DynamicParticle *aDynamicGamma = aTrack.GetDynamicParticle();
215 G4double Egam = aDynamicGamma->GetKineticEnergy();
216 if (Egam < 4*Mmuon) return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
217 G4ParticleMomentum GammaDirection = aDynamicGamma->GetMomentumDirection();
218
219 // select randomly one element constituting the material
220 const G4Element& anElement = *SelectRandomAtom(aDynamicGamma, aMaterial);
221 G4double Z = anElement.GetZ();
222 G4double A = anElement.GetA()/(g/mole);
223
224 static G4double Zlast=0,B,Dn,Zthird,Winfty,A027,C1Num2,C2Term2;
225 if(Zlast!=Z) // the element has changed
226 { Zlast=Z;
227 if(Z==1) // special case of Hydrogen
228 { B=202.4;
229 Dn=1.49;
230 }
231 else
232 { B=183.;
233 Dn=1.54*pow(A,0.27);
234 }
235 Zthird=pow(Z,-1./3.); // Z**(-1/3)
236 Winfty=B*Zthird*Mmuon/(Dn*Mele);
237 A027=pow(A,0.27);
238 G4double C1Num=0.35*A027;
239 C1Num2=C1Num*C1Num;
240 C2Term2=Mele/(183.*Zthird*Mmuon);
241 }
242
243 G4double GammaMuonInv=Mmuon/Egam;
244 G4double sqrtx=sqrt(.25-GammaMuonInv);
245 G4double xmax=.5+sqrtx;
246 G4double xmin=.5-sqrtx;
247
248 // generate xPlus according to the differential cross section by rejection
249 G4double Ds2=(Dn*sqrte-2.);
250 G4double sBZ=sqrte*B*Zthird/Mele;
251 G4double LogWmaxInv=1./log(Winfty*(1.+2.*Ds2*GammaMuonInv)
252 /(1.+2.*sBZ*Mmuon*GammaMuonInv));
253 G4double xPlus,xMinus,xPM,result,W;
254 do
255 { xPlus=xmin+G4UniformRand()*(xmax-xmin);
256 xMinus=1.-xPlus;
257 xPM=xPlus*xMinus;
258 G4double del=Mmuon*Mmuon/(2.*Egam*xPM);
259 W=Winfty*(1.+Ds2*del/Mmuon)/(1.+sBZ*del);
260 if(W<1.) W=1.; // to avoid negative cross section at xmin
261 G4double xxp=1.-4./3.*xPM; // the main xPlus dependence
262 result=xxp*log(W)*LogWmaxInv;
263 if(result>1.)
264 { G4cout << "error in dSigxPlusGen, result=" << result << " is >1" << '\n';
265 exit(10);
266 }
267 }
268 while (G4UniformRand() > result);
269
270 // now generate the angular variables via the auxilary variables t,psi,rho
271 G4double t;
272 G4double psi;
273 G4double rho;
274
275 G4double thetaPlus,thetaMinus,phiHalf; // final angular variables
276
277 do // t, psi, rho generation start (while angle < pi)
278 {
279 //generate t by the rejection method
280 G4double C1=C1Num2* GammaMuonInv/xPM;
281 G4double f1_max=(1.-xPM) / (1.+C1);
282 G4double f1; // the probability density
283 do
284 { t=G4UniformRand();
285 f1=(1.-2.*xPM+4.*xPM*t*(1.-t)) / (1.+C1/(t*t));
286 if(f1<0 || f1> f1_max) // should never happend
287 { G4cout << "outside allowed range f1=" << f1 << G4endl;
288 exit(1);
289 }
290 }
291 while ( G4UniformRand()*f1_max > f1);
292 // generate psi by the rejection method
293 G4double f2_max=1.-2.*xPM*(1.-4.*t*(1.-t));
294
295 // long version
296 G4double f2;
297 do
298 { psi=2.*pi*G4UniformRand();
299 f2=1.-2.*xPM+4.*xPM*t*(1.-t)*(1.+cos(2.*psi));
300 if(f2<0 || f2> f2_max) // should never happend
301 { G4cout << "outside allowed range f2=" << f2 << G4endl;
302 exit(1);
303 }
304 }
305 while ( G4UniformRand()*f2_max > f2);
306
307 // generate rho by direct transformation
308 G4double C2Term1=GammaMuonInv/(2.*xPM*t);
309 G4double C2=4./sqrt(xPM)*pow(C2Term1*C2Term1+C2Term2*C2Term2,2.);
310 G4double rhomax=1.9/A027*(1./t-1.);
311 G4double beta=log( (C2+pow(rhomax,4.))/C2 );
312 rho=pow(C2 *( exp(beta*G4UniformRand())-1. ) ,0.25);
313
314 //now get from t and psi the kinematical variables
315 G4double u=sqrt(1./t-1.);
316 G4double xiHalf=0.5*rho*cos(psi);
317 phiHalf=0.5*rho/u*sin(psi);
318
319 thetaPlus =GammaMuonInv*(u+xiHalf)/xPlus;
320 thetaMinus=GammaMuonInv*(u-xiHalf)/xMinus;
321
322 } while ( std::abs(thetaPlus)>pi || std::abs(thetaMinus) >pi);
323
324 // now construct the vectors
325 // azimuthal symmetry, take phi0 at random between 0 and 2 pi
326 G4double phi0=2.*pi*G4UniformRand();
327 G4double EPlus=xPlus*Egam;
328 G4double EMinus=xMinus*Egam;
329
330 // mu+ mu- directions for gamma in z-direction
331 G4ThreeVector MuPlusDirection ( sin(thetaPlus) *cos(phi0+phiHalf),
332 sin(thetaPlus) *sin(phi0+phiHalf), cos(thetaPlus) );
333 G4ThreeVector MuMinusDirection (-sin(thetaMinus)*cos(phi0-phiHalf),
334 -sin(thetaMinus) *sin(phi0-phiHalf), cos(thetaMinus) );
335 // rotate to actual gamma direction
336 MuPlusDirection.rotateUz(GammaDirection);
337 MuMinusDirection.rotateUz(GammaDirection);
338 aParticleChange.SetNumberOfSecondaries(2);
339 // create G4DynamicParticle object for the particle1
340 G4DynamicParticle* aParticle1= new G4DynamicParticle(
341 G4MuonPlus::MuonPlus(),MuPlusDirection,EPlus-Mmuon);
342 aParticleChange.AddSecondary(aParticle1);
343 // create G4DynamicParticle object for the particle2
344 G4DynamicParticle* aParticle2= new G4DynamicParticle(
345 G4MuonMinus::MuonMinus(),MuMinusDirection,EMinus-Mmuon);
346 aParticleChange.AddSecondary(aParticle2);
347 //
348 // Kill the incident photon
349 //
350 aParticleChange.ProposeMomentumDirection( 0., 0., 0. ) ;
351 aParticleChange.ProposeEnergy( 0. ) ;
352 aParticleChange.ProposeTrackStatus( fStopAndKill ) ;
353 // Reset NbOfInteractionLengthLeft and return aParticleChange
354 return G4VDiscreteProcess::PostStepDoIt( aTrack, aStep );
355}
356
357//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
358
359G4Element* G4GammaConversionToMuons::SelectRandomAtom(
360 const G4DynamicParticle* aDynamicGamma,
361 G4Material* aMaterial)
362{
363 // select randomly 1 element within the material, invoked by PostStepDoIt
364
365 const G4int NumberOfElements = aMaterial->GetNumberOfElements();
366 const G4ElementVector* theElementVector = aMaterial->GetElementVector();
367 if (NumberOfElements == 1) return (*theElementVector)[0];
368
369 const G4double* NbOfAtomsPerVolume = aMaterial->GetVecNbOfAtomsPerVolume();
370
371 G4double PartialSumSigma = 0. ;
372 G4double rval = G4UniformRand()/MeanFreePath;
373
374
375 for ( G4int i=0 ; i < NumberOfElements ; i++ )
376 { PartialSumSigma += NbOfAtomsPerVolume[i] *
377 GetCrossSectionPerAtom(aDynamicGamma, (*theElementVector)[i]);
378 if (rval <= PartialSumSigma) return ((*theElementVector)[i]);
379 }
380 G4cout << " WARNING !!! - The Material '"<< aMaterial->GetName()
381 << "' has no elements, NULL pointer returned." << G4endl;
382 return NULL;
383}
384
385//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
386
387void G4GammaConversionToMuons::PrintInfoDefinition()
388{
389 G4String comments ="gamma->mu+mu- Bethe Heitler process.\n";
390 G4cout << G4endl << GetProcessName() << ": " << comments
391 << " good cross section parametrization from "
392 << G4BestUnit(LowestEnergyLimit,"Energy")
393 << " to " << HighestEnergyLimit/GeV << " GeV for all Z." << G4endl;
394}
395
396//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
Note: See TracBrowser for help on using the repository browser.