p_excit_n_m Function

public function p_excit_n_m(eb, n, m) result(sigma)

Calculates the cross section for a proton-Hydrogen impact excitation transition from the n state to the m state at energy eb

Equation

References

Arguments

Type IntentOptional AttributesName
real(kind=Float64), intent(in) :: eb

Relative collision energy [keV/amu]

integer, intent(in) :: n

Initial atomic energy level/state

integer, intent(in) :: m

Final atomic energy level/state

Return Value real(kind=Float64)

Cross Section []


Calls

proc~~p_excit_n_m~~CallsGraph proc~p_excit_n_m p_excit_n_m proc~p_excit_n p_excit_n proc~p_excit_n_m->proc~p_excit_n proc~p_excit_1_janev p_excit_1_janev proc~p_excit_n->proc~p_excit_1_janev proc~p_excit_2_janev p_excit_2_janev proc~p_excit_n->proc~p_excit_2_janev proc~p_excit_3_janev p_excit_3_janev proc~p_excit_n->proc~p_excit_3_janev proc~p_excit_1_2_janev p_excit_1_2_janev proc~p_excit_1_janev->proc~p_excit_1_2_janev proc~p_excit_1_3_janev p_excit_1_3_janev proc~p_excit_1_janev->proc~p_excit_1_3_janev proc~p_excit_1_6_janev p_excit_1_6_janev proc~p_excit_1_janev->proc~p_excit_1_6_janev proc~p_excit_1_5_janev p_excit_1_5_janev proc~p_excit_1_janev->proc~p_excit_1_5_janev proc~p_excit_1_4_janev p_excit_1_4_janev proc~p_excit_1_janev->proc~p_excit_1_4_janev proc~p_excit_2_6_janev p_excit_2_6_janev proc~p_excit_2_janev->proc~p_excit_2_6_janev proc~p_excit_2_3_janev p_excit_2_3_janev proc~p_excit_2_janev->proc~p_excit_2_3_janev proc~p_excit_2_9_janev p_excit_2_9_janev proc~p_excit_2_janev->proc~p_excit_2_9_janev proc~p_excit_2_8_janev p_excit_2_8_janev proc~p_excit_2_janev->proc~p_excit_2_8_janev proc~p_excit_2_5_janev p_excit_2_5_janev proc~p_excit_2_janev->proc~p_excit_2_5_janev proc~p_excit_2_4_janev p_excit_2_4_janev proc~p_excit_2_janev->proc~p_excit_2_4_janev proc~p_excit_2_7_janev p_excit_2_7_janev proc~p_excit_2_janev->proc~p_excit_2_7_janev proc~p_excit_2_10_janev p_excit_2_10_janev proc~p_excit_2_janev->proc~p_excit_2_10_janev proc~p_excit_3_10_janev p_excit_3_10_janev proc~p_excit_3_janev->proc~p_excit_3_10_janev proc~p_excit_3_4_janev p_excit_3_4_janev proc~p_excit_3_janev->proc~p_excit_3_4_janev proc~p_excit_3_5_janev p_excit_3_5_janev proc~p_excit_3_janev->proc~p_excit_3_5_janev proc~p_excit_3_6_janev p_excit_3_6_janev proc~p_excit_3_janev->proc~p_excit_3_6_janev proc~p_excit_3_8_janev p_excit_3_8_janev proc~p_excit_3_janev->proc~p_excit_3_8_janev proc~p_excit_3_9_janev p_excit_3_9_janev proc~p_excit_3_janev->proc~p_excit_3_9_janev proc~p_excit_3_7_janev p_excit_3_7_janev proc~p_excit_3_janev->proc~p_excit_3_7_janev proc~p_excit_3_10_janev->proc~p_excit_3_6_janev proc~p_excit_2_6_janev->proc~p_excit_2_5_janev proc~p_excit_2_9_janev->proc~p_excit_2_5_janev proc~p_excit_2_8_janev->proc~p_excit_2_5_janev proc~p_excit_3_8_janev->proc~p_excit_3_6_janev proc~p_excit_3_9_janev->proc~p_excit_3_6_janev proc~p_excit_2_7_janev->proc~p_excit_2_5_janev proc~p_excit_3_7_janev->proc~p_excit_3_6_janev proc~p_excit_2_10_janev->proc~p_excit_2_5_janev

Contents

Source Code


Source Code

function p_excit_n_m(eb, n, m) result(sigma)
    !+Calculates the cross section for a proton-Hydrogen impact excitation transition from
    !+the `n` state to the `m` state at energy `eb`
    !+
    !+###Equation
    !+$$ H^+ + H(n) \rightarrow H^+ + H(m), m \gt n $$
    !+
    !+###References
    !+* Eq. 29.b and Table 4 in Ref. 2 for \(n = 1\) and \(m = 2\) [[atomic_tables(module)]]
    !+* Eq. 30 and Table 5 in Ref. 2 for \(n = 1\) and \(m = 3-6\) [[atomic_tables(module)]]
    !+* Eq. 31 and Table 5 in Ref. 2 for \(n = 1\) and \(m \gt 6\) [[atomic_tables(module)]]
    !+* Eq. 32 and Table 6 in Ref. 2 for \(n = 2\) and \(m \le 5\) [[atomic_tables(module)]]
    !+* Eq. 33 and Table 6 in Ref. 2 for \(n = 2\) and \(m = 6-10\) [[atomic_tables(module)]]
    !+* Eq. 34 and Table 6 in Ref. 2 for \(n = 2\) and \(m \gt 10\) [[atomic_tables(module)]]
    !+* Eq. 35 and Table 7 in Ref. 2 for \(n = 3\) and \(m \le 6\) [[atomic_tables(module)]]
    !+* Eq. 36 and Table 7 in Ref. 2 for \(n = 3\) and \(m = 7-10\) [[atomic_tables(module)]]
    !+* Eq. 37 and Table 7 in Ref. 2 for \(n = 3\) and \(m \gt 10\) [[atomic_tables(module)]]
    !+* Eq. 38-39 in Ref. 2 for \(n \gt 3\) and \(m \gt 4\) [[atomic_tables(module)]]
    !+
    real(Float64), intent(in) :: eb
        !+ Relative collision energy [keV/amu]
    integer, intent(in)       :: n
        !+ Initial atomic energy level/state
    integer, intent(in)       :: m
        !+ Final atomic energy level/state
    real(Float64)             :: sigma
        !+ Cross Section [\(cm^2\)]

    real(Float64), dimension(12) :: sigma_m

    sigma_m = p_excit_n(eb, n, 12)
    if(m.le.0) then
        sigma = sum(sigma_m)
    else
        sigma = sigma_m(m)
    endif

end function p_excit_n_m