Aq_excit_1_janev Function

public function Aq_excit_1_janev(eb, q, m_max) result(sigma)

Calculates an array of the excitation cross sections for a neutral Hydrogen atom transitioning from the state to the m=1..m_max states due to a collision an ion with charge q at energy eb

Equation

References

Arguments

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

Collision energy [keV]

integer, intent(in) :: q

Ion charge

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 state to the m'th state []


Calls

proc~~aq_excit_1_janev~~CallsGraph proc~aq_excit_1_janev Aq_excit_1_janev proc~aq_excit_1_2_janev Aq_excit_1_2_janev proc~aq_excit_1_janev->proc~aq_excit_1_2_janev proc~aq_excit_1_3_janev Aq_excit_1_3_janev proc~aq_excit_1_janev->proc~aq_excit_1_3_janev proc~aq_excit_1_4_janev Aq_excit_1_4_janev proc~aq_excit_1_janev->proc~aq_excit_1_4_janev proc~aq_excit_1_5_janev Aq_excit_1_5_janev proc~aq_excit_1_janev->proc~aq_excit_1_5_janev

Called by

proc~~aq_excit_1_janev~~CalledByGraph proc~aq_excit_1_janev Aq_excit_1_janev proc~aq_excit_n Aq_excit_n proc~aq_excit_n->proc~aq_excit_1_janev proc~aq_excit_n_m Aq_excit_n_m proc~aq_excit_n_m->proc~aq_excit_n proc~aq_excit Aq_excit proc~aq_excit->proc~aq_excit_n proc~write_bb_h_aq write_bb_H_Aq proc~write_bb_h_aq->proc~aq_excit program~generate_tables generate_tables program~generate_tables->proc~write_bb_h_aq

Contents

Source Code


Source Code

function Aq_excit_1_janev(eb, q, m_max) result(sigma)
    !+Calculates an array of the excitation cross sections for a neutral Hydrogen atom transitioning from
    !+the \(n=1\) state to the m=1..`m_max` states due to a collision an ion with charge `q` at energy `eb`
    !+
    !+###Equation
    !+$$ A^{q+} + H(1) \rightarrow A^{q+} + H(m=1..m_{max}), q \gt 4 $$
    !+
    !+###References
    !+* Page 132 in Ref. 5 [[atomic_tables(module)]]
    !+* Page 134 in Ref. 5 [[atomic_tables(module)]]
    !+* Page 136 in Ref. 5 [[atomic_tables(module)]]
    !+
    real(Float64), intent(in)       :: eb
        !+ Collision energy [keV]
    integer, intent(in)             :: q
        !+ Ion charge
    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=1\) state to the m'th state [\(cm^2\)]

    integer :: m

    sigma = 0.d0
    do m=1,m_max
        select case (m)
            case (1)
                sigma(1) = 0.d0
            case (2)
                sigma(2) = Aq_excit_1_2_janev(eb, q)
            case (3)
                sigma(3) = Aq_excit_1_3_janev(eb, q)
            case (4)
                sigma(4) = Aq_excit_1_4_janev(eb, q)
            case (5)
                sigma(5) = Aq_excit_1_5_janev(eb, q)
            case DEFAULT
                sigma(m) = sigma(5)*(5.0/m)**3.0
        end select
    enddo

end function Aq_excit_1_janev