e_excit_n Function

public function e_excit_n(eb, n, m_max) result(sigma)

Calculates an array of cross sections for a electron-Hydrogen impact excitation transition from the n state to the state at energy eb

Equation

References

Arguments

TypeIntentOptionalAttributesName
real(kind=Float64), intent(in) :: eb

Collision energy [keV]

integer, intent(in) :: n

Initial atomic energy level/state

integer, intent(in) :: m_max

Number of m states to calculate

Return Value real(kind=Float64),dimension(m_max)

Array of cross sections where the m'th index refers to an excitation from the n m state []


Calls

proc~~e_excit_n~~CallsGraph proc~e_excit_n e_excit_n proc~e_excit_1_janev e_excit_1_janev proc~e_excit_n->proc~e_excit_1_janev proc~e_excit_2_3_janev e_excit_2_3_janev proc~e_excit_n->proc~e_excit_2_3_janev proc~e_excit_f e_excit_f proc~e_excit_n->proc~e_excit_f proc~e_excit_1_janev->proc~e_excit_f proc~e_excit_1_2_janev e_excit_1_2_janev proc~e_excit_1_janev->proc~e_excit_1_2_janev proc~e_excit_1_3_janev e_excit_1_3_janev proc~e_excit_1_janev->proc~e_excit_1_3_janev proc~e_excit_1_4_janev e_excit_1_4_janev proc~e_excit_1_janev->proc~e_excit_1_4_janev proc~e_excit_1_5_janev e_excit_1_5_janev proc~e_excit_1_janev->proc~e_excit_1_5_janev

Called by

proc~~e_excit_n~~CalledByGraph proc~e_excit_n e_excit_n proc~e_excit_n_m e_excit_n_m proc~e_excit_n_m->proc~e_excit_n proc~e_excit e_excit proc~e_excit->proc~e_excit_n proc~write_bb_h_e write_bb_H_e proc~write_bb_h_e->proc~e_excit program~generate_tables generate_tables program~generate_tables->proc~write_bb_h_e

Contents

Source Code


Source Code

function e_excit_n(eb, n, m_max) result(sigma)
    !+Calculates an array of cross sections for a electron-Hydrogen impact excitation transition from
    !+the `n` state to the \(m=1..m_{max}\) state at energy `eb`
    !+
    !+###Equation
    !+$$ e + H(n) \rightarrow e + H(m=1..m_{max}), m \gt n $$
    !+
    !+###References
    !+* Eq. 5 and Table 2 in Ref. 2 [[atomic_tables(module)]]
    !+* Eqs. 6-7 in Ref. 2 [[atomic_tables(module)]]
    !+* Section 2.1.1 B in Ref. 2 [[atomic_tables(module)]]
    !+
    real(Float64), intent(in)       :: eb
        !+ Collision energy [keV]
    integer, intent(in)             :: n
        !+ Initial atomic energy level/state
    integer, intent(in)             :: m_max
        !+ Number of m states to calculate
    real(Float64), dimension(m_max) :: sigma
        !+ Array of cross sections where the m'th index refers to
        !+an excitation from the `n`\(\rightarrow\)`m` state [\(cm^2\)]

    integer :: m
    real(Float64) :: nf, mf
    real(Float64) :: x, y, A, B, bn, r, deltaE

    nf = real(n)

    if(n.eq.1) then
        sigma = e_excit_1_janev(eb, m_max)
    else
        m_loop: do m=1,m_max
            mf = real(m)
            if(n.ge.m) then
                sigma(m) = 0.d0
                cycle m_loop
            endif
            if((n.eq.2).and.(m.eq.3)) then
                sigma(m) = e_excit_2_3_janev(eb)
            else
                deltaE=13.6*(1.0/nf**2 - 1.0/mf**2)
                x = (eb*1.d3)/deltaE
                y = 1.0 - (nf/mf)**2
                r = 1.94/nf**1.57

                A = 2.0 * nf**2 * e_excit_f(n,m)/y
                bn = 1.0/nf*(4.0 - 18.63/nf + 36.24/nf**2 - 28.09/nf**3)
                B = 4.0 * nf**4/(mf**3*y**2)*(1.0 + 4.0/(3.0*y) + bn/y**2.0)

                sigma(m) = 1.76e-16*nf**2/(y*x)*(1.0 - exp(-r*y*x))* &
                           (A*(log(x) + 1.0/(2.0*x)) + (B - A*log(2.0*n**2.0/y))* &
                           (1.0 - 1.0/x))

                if(x.le.1.0) sigma(m) = 0.d0
            endif
        enddo m_loop
    endif

end function e_excit_n