MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
mod_bremsstrahlung.t
Go to the documentation of this file.
1 ! Several bremsstrahlung cross sections are provided in this module
3  use mod_physics
4  use mod_constants
5  implicit none
6 
7  double precision, parameter :: r0=2.8179e-13,alpha=7.2992d-3
8  double precision, parameter :: mec2=511.d0
9 
10 contains
11 
12  subroutine bremcross_kramers(Z,Ee,Eph,sigmaB)
13  ! Kramers bremsstrahlung cross section
14  ! Z:ion change; Ee:electron energy (keV); Eph:photon energy (keV)
15  ! sigmaB: cross section (cm^-2 keV^-1)
16 
17  integer :: Z
18  double precision :: Eph,Ee,sigmaB
19 
20  double precision :: E1,E2,p1,p2,k
21  double precision :: sigma0
22 
23  sigmab=0.d0
24 
25  if (ee>eph) then
26  sigma0=(16.d0/3)*alpha*r0**2
27  e1=ee/mec2 ! initial electron kinetic energy in units of mec2
28  k=eph/mec2 ! photon energy in units of mec2
29  e2=e1-k ! kinetic energy of scattered electron
30  p1=sqrt(e1*(2+e1))
31  p2=sqrt(e2*(2+e2))
32  sigmab=sigma0*z**2/(k*p1**2)
33  sigmab=sigmab/mec2
34  endif
35 
36  end subroutine bremcross_kramers
37 
38  subroutine bremcross_betheheitler(Z,Ee,Eph,sigmaB)
39  ! Bethe-Heitler bremsstrahlung cross section
40  ! Formula 3BN in Koch, Motz, Reviews of Modern Physics 31, 920 (1959)
41  ! Z:ion change; Ee:electron energy (keV); Eph:photon energy (keV)
42  ! sigmaB: cross section (cm^-2 keV^-1)
43 
44  integer :: Z
45  double precision :: Eph,Ee,sigmaB
46 
47  double precision :: E1,E2,p1,p2,k
48  double precision :: sigma0
49 
50  sigmab=0.d0
51 
52  if (ee>eph) then
53  sigma0=(16.d0/3)*alpha*r0**2
54  e1=ee/mec2 ! initial electron kinetic energy in units of mec2
55  k=eph/mec2 ! photon energy in units of mec2
56  e2=e1-k ! kinetic energy of scattered electron
57  p1=sqrt(e1*(2+e1))
58  p2=sqrt(e2*(2+e2))
59  sigmab=sigma0*z**2/(k*p1**2)*log((p1+p2)/(p1-p2))
60  sigmab=sigmab/mec2
61  endif
62 
63  end subroutine bremcross_betheheitler
64 
65  subroutine bremcross_haug(Z,Ee,Eph,sigmaB)
66  ! Haug bremsstrahlung cross section
67  ! Equation (4) and (5) in Haug, Astron. Astrophys. 326, 417 (1997)
68  ! Z:ion change; Ee:electron energy (keV); Eph:photon energy (keV)
69  ! sigmaB: cross section (cm^-2 keV^-1)
70 
71  integer :: Z
72  double precision :: Eph,Ee,sigmaB
73 
74  double precision :: E1,E2,Et1,Et2,p1,p2,k
75  double precision :: sigma0
76  double precision :: term1,term2,term3,a1,a2,fE
77 
78  sigmab=0.d0
79 
80  if (ee>eph) then
81  sigma0=2.d0*alpha*r0**2
82  e1=ee/mec2 ! initial electron kinetic energy in units of mec2
83  k=eph/mec2 ! photon energy in units of mec2
84  e2=e1-k ! kinetic energy of scattered electron
85  p1=sqrt(e1*(2+e1))
86  p2=sqrt(e2*(2+e2))
87  et1=e1+1.d0 ! initial electron total energy
88  et2=e2+1.d0 ! total energy of scattered electron
89 
90  term1=(4.d0/3)*et1*et2+k**2-(7.d0/15)*(k**2/(et1*et2))
91  term1=term1-(11.d0/70)*(k**2)*(p1**2+p2**2)/(et1*et2)**4
92  term2=2*log((et1*et2+p1*p2-1)/k)
93  term3=1+(1/(et1*et2))+(7.d0/20)*(p1*2+p2*2)/(et1*et2)**3
94  term3=term3+((9.d0/28)*k**2+(263.d0/210)*(p1**2)*(p2**2))/(et1*et2)**3
95  term3=term3*(p1*p2)/(et1*et2)
96 
97  a1=alpha*z*et1/p1
98  a2=alpha*z*et2/p2
99  fe=(a2/a1)*(1.d0-exp(-2*dpi*a1))/(1.d0-exp(-2*dpi*a2))
100 
101  sigmab=sigma0*z**2/(k*p1**2)*term1*(term2-term3)*fe
102  sigmab=sigmab/mec2
103  endif
104 
105  end subroutine bremcross_haug
106 
107 end module mod_bremsstrahlung
double precision, parameter alpha
subroutine bremcross_betheheitler(Z, Ee, Eph, sigmaB)
subroutine bremcross_kramers(Z, Ee, Eph, sigmaB)
double precision, parameter r0
double precision, parameter mec2
subroutine bremcross_haug(Z, Ee, Eph, sigmaB)
Module for physical and numeric constants.
Definition: mod_constants.t:2
double precision, parameter dpi
Pi.
Definition: mod_constants.t:22
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition: mod_physics.t:4