Calculates cross section for a proton-Hydrogen charge exchange interaction
from the state to the state at energy Erel
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=Float64), | intent(in) | :: | Erel | Relative collision energy [keV/amu] |
Cross Section []
function p_cx_2_2_adas(Erel) result(sigma)
!+Calculates cross section for a proton-Hydrogen charge exchange interaction
!+from the \(n=2\) state to the \(m=2\) state at energy `Erel`
!+
!+###Equation
!+ $$H^+ + H(2) \rightarrow H(2) + H^+$$
!+###References
!+* Ref. 4 [[atomic_tables(module)]]
real(Float64), intent(in) :: Erel
!+ Relative collision energy [keV/amu]
real(Float64) :: sigma
!+ Cross Section [\(cm^2\)]
real(Float64), dimension(11), parameter :: a2s = [-1.896015167d6, 4.431727330d6, &
-4.627815357d6, 2.843068107d6, &
-1.137952956d6, 3.100801094d5, &
-5.825744660d4, 7.452319142d3, &
-6.212350647d2, 3.047712749d1, &
-6.682658463d-1 ]
real(Float64), dimension(11), parameter :: a2p = [-1.614213508d5, 3.772469288d5, &
-3.924736424d5, 2.393127027d5, &
-9.470300966d4, 2.541276100d4, &
-4.682860453d3, 5.851219013d2, &
-4.744504549d1, 2.254460913d0, &
-4.767235839d-2 ]
real(Float64), parameter :: n = 2.d0
real(Float64) :: e, ee, fac, l, sigma2s, sigma2p
e = Erel * 1.d3 * n**2.0
if(Erel.le.1.5d2) then
ee = max(e, 1.d3)
fac = 1.d0
else
ee = 1.5e5 * n**2.d0
fac = 2.d15 * ((e*1.d-3)**(-5.5))
endif
l = log10(ee)
sigma2s = a2s(1) + a2s(2)*l + a2s(3)*l**2.0 + a2s(4)*l**3.0 + &
a2s(5)*l**4.0 + a2s(6)*l**5.0 + a2s(7)*l**6.0 + &
a2s(8)*l**7.0 + a2s(9)*l**8.0 + a2s(10)*l**9.0 + a2s(11)*l**10.0
sigma2s = 10.d0**(sigma2s)
sigma2p = a2p(1) + a2p(2)*l + a2p(3)*l**2.0 + a2p(4)*l**3.0 + &
a2p(5)*l**4.0 + a2p(6)*l**5.0 + a2p(7)*l**6.0 + &
a2p(8)*l**7.0 + a2p(9)*l**8.0 + a2p(10)*l**9.0 + a2p(11)*l**10.0
sigma2p = 10.d0**(sigma2p)
sigma = fac*(0.25*sigma2s + 0.75*sigma2p)
end function p_cx_2_2_adas