p_cx Function

public function p_cx(Erel, n_max, m_max) result(sigma)

Calculates a matrix of cross sections for proton-Hydrogen charge exchange interactions from the states at energy Erel

Equation

References

Arguments

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

Relative collision energy [keV/amu]

integer, intent(in) :: n_max

Number of initial atomic energy levels/states

integer, intent(in) :: m_max

Number of final atomic energy levels/states

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

Matrix of cross sections where the subscripts correspond to the transitions: p_cx[n,m] []


Calls

proc~~p_cx~~CallsGraph proc~p_cx p_cx proc~p_cx_n p_cx_n proc~p_cx->proc~p_cx_n proc~p_cx_janev p_cx_janev proc~p_cx_n->proc~p_cx_janev proc~m_spread m_spread proc~p_cx_n->proc~m_spread proc~p_cx_3_5_adas p_cx_3_5_adas proc~p_cx_n->proc~p_cx_3_5_adas proc~p_cx_3_4_adas p_cx_3_4_adas proc~p_cx_n->proc~p_cx_3_4_adas proc~p_cx_2 p_cx_2 proc~p_cx_n->proc~p_cx_2 proc~p_cx_3 p_cx_3 proc~p_cx_n->proc~p_cx_3 proc~p_cx_1 p_cx_1 proc~p_cx_n->proc~p_cx_1 proc~p_cx_1_4_adas p_cx_1_4_adas proc~p_cx_n->proc~p_cx_1_4_adas proc~p_cx_1_janev p_cx_1_janev proc~p_cx_janev->proc~p_cx_1_janev proc~p_cx_2_janev p_cx_2_janev proc~p_cx_janev->proc~p_cx_2_janev proc~p_cx_n_janev p_cx_n_janev proc~p_cx_janev->proc~p_cx_n_janev proc~p_cx_3_janev p_cx_3_janev proc~p_cx_janev->proc~p_cx_3_janev proc~p_cx_2->proc~p_cx_janev proc~p_cx_2->proc~m_spread proc~p_cx_1_2_adas p_cx_1_2_adas proc~p_cx_2->proc~p_cx_1_2_adas proc~p_cx_2_3_adas p_cx_2_3_adas proc~p_cx_2->proc~p_cx_2_3_adas proc~p_cx_2_2_adas p_cx_2_2_adas proc~p_cx_2->proc~p_cx_2_2_adas proc~p_cx_3->proc~p_cx_janev proc~p_cx_3->proc~m_spread proc~p_cx_3->proc~p_cx_3_5_adas proc~p_cx_3->proc~p_cx_3_4_adas proc~p_cx_3->proc~p_cx_1 proc~p_cx_3_2_adas p_cx_3_2_adas proc~p_cx_3->proc~p_cx_3_2_adas proc~p_cx_3_3_adas p_cx_3_3_adas proc~p_cx_3->proc~p_cx_3_3_adas proc~p_cx_3_6inf_adas p_cx_3_6inf_adas proc~p_cx_3->proc~p_cx_3_6inf_adas proc~p_cx_1_3_adas p_cx_1_3_adas proc~p_cx_3->proc~p_cx_1_3_adas proc~p_cx_1->proc~p_cx_janev proc~p_cx_1->proc~p_cx_1_4_adas proc~p_cx_1_1_adas p_cx_1_1_adas proc~p_cx_1->proc~p_cx_1_1_adas proc~p_cx_1->proc~p_cx_1_2_adas proc~p_cx_1->proc~p_cx_1_3_adas proc~p_cx_1_2_janev p_cx_1_2_janev proc~p_cx_1->proc~p_cx_1_2_janev proc~aljan1 aljan1 proc~p_cx_1_2_janev->proc~aljan1

Called by

proc~~p_cx~~CalledByGraph proc~p_cx p_cx proc~write_bb_h_h write_bb_H_H proc~write_bb_h_h->proc~p_cx program~generate_tables generate_tables program~generate_tables->proc~write_bb_h_h

Contents

Source Code


Source Code

function p_cx(Erel, n_max, m_max) result(sigma)
    !+Calculates a matrix of cross sections for proton-Hydrogen charge exchange interactions
    !+from the \(n=1..n_{max} \rightarrow m=1..m_{max}\) states at energy `Erel`
    !+
    !+###Equation
    !+ $$H^+ + H(n=1..n_{max}) \rightarrow H(m=1..m_{max}) + H^+$$
    !+###References
    !+* Ref. 2 [[atomic_tables(module)]]
    !+* Ref. 4 [[atomic_tables(module)]]
    !+* Ref. 8 [[atomic_tables(module)]]
    real(Float64), intent(in)             :: Erel
        !+ Relative collision energy [keV/amu]
    integer, intent(in)                   :: n_max
        !+ Number of initial atomic energy levels/states
    integer, intent(in)                   :: m_max
        !+ Number of final atomic energy levels/states
    real(Float64), dimension(n_max,m_max) :: sigma
        !+ Matrix of cross sections where the subscripts correspond
        !+ to the \(n \rightarrow m\) transitions: p_cx[n,m] [\(cm^2\)]

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

    integer :: n, m

    do n=1,12
        sigma_full(n,:) = p_cx_n(Erel, n, 12)
    enddo

    sigma = sigma_full(1:n_max,1:m_max)

end function p_cx