get_radial_coupling_term_value Subroutine

public subroutine get_radial_coupling_term_value(intermolecular_distance, lambda_, v_, j_, v_prime_, j_prime_, radial_term_value_)

Returns the interpolated value of a specific radial coupling term at a given distance.

Arguments

Type IntentOptional Attributes Name
real(kind=dp), intent(in) :: intermolecular_distance

Intermolecular distance, \(R\)

integer(kind=int32), intent(in) :: lambda_

Legendre expansion index

integer(kind=int32), intent(in) :: v_

pre-collisional vibrational quantum number

integer(kind=int32), intent(in) :: j_

pre-collisional rotational quantum number

integer(kind=int32), intent(in) :: v_prime_

post-collisional vibrational quantum number

integer(kind=int32), intent(in) :: j_prime_

post-collisional rotational quantum number

real(kind=dp), intent(out) :: radial_term_value_

Value of the radial coupling coefficient


Calls

proc~~get_radial_coupling_term_value~~CallsGraph proc~get_radial_coupling_term_value get_radial_coupling_term_value proc~find_lambda_index find_lambda_index proc~get_radial_coupling_term_value->proc~find_lambda_index proc~find_coupling_index find_coupling_index proc~get_radial_coupling_term_value->proc~find_coupling_index proc~ispline ispline proc~get_radial_coupling_term_value->proc~ispline proc~handle_lambda_index_error handle_lambda_index_error proc~get_radial_coupling_term_value->proc~handle_lambda_index_error proc~handle_coupling_index_error handle_coupling_index_error proc~get_radial_coupling_term_value->proc~handle_coupling_index_error proc~float_to_character float_to_character proc~ispline->proc~float_to_character proc~write_warning write_warning proc~ispline->proc~write_warning proc~integer_to_character integer_to_character proc~handle_lambda_index_error->proc~integer_to_character proc~write_error write_error proc~handle_lambda_index_error->proc~write_error proc~handle_coupling_index_error->proc~integer_to_character proc~handle_coupling_index_error->proc~write_error proc~write_message write_message proc~write_error->proc~write_message proc~write_warning->proc~write_message

Called by

proc~~get_radial_coupling_term_value~~CalledByGraph proc~get_radial_coupling_term_value get_radial_coupling_term_value proc~calculate_single_pes_matrix_element calculate_single_pes_matrix_element proc~calculate_single_pes_matrix_element->proc~get_radial_coupling_term_value proc~calculate_pes_matrix calculate_pes_matrix proc~calculate_pes_matrix->proc~calculate_single_pes_matrix_element proc~initial_setup initial_setup proc~initial_setup->proc~calculate_pes_matrix proc~general_propagation_step general_propagation_step proc~general_propagation_step->proc~calculate_pes_matrix proc~numerov numerov proc~numerov->proc~initial_setup proc~numerov->proc~general_propagation_step proc~handle_final_propagation_steps handle_final_propagation_steps proc~numerov->proc~handle_final_propagation_steps proc~handle_final_propagation_steps->proc~general_propagation_step program~scattering SCATTERING program~scattering->proc~numerov

Contents


Variables

Type Visibility Attributes Name Initial
integer(kind=int32), private :: coupling_index
integer(kind=int32), private :: lambda_index

Source Code

      subroutine get_radial_coupling_term_value(intermolecular_distance,       &
         lambda_, v_, j_, v_prime_, j_prime_, radial_term_value_)
         !! Returns the interpolated value of a specific radial coupling term
         !! at a given distance.
         !---------------------------------------------------------------------!
         real(dp), intent(in) :: intermolecular_distance
            !! Intermolecular distance, \\(R\\)
         integer(int32), intent(in) :: lambda_
            !! Legendre expansion index
         integer(int32), intent(in) :: v_
            !! pre-collisional vibrational quantum number
         integer(int32), intent(in) :: j_
            !! pre-collisional rotational quantum number
         integer(int32), intent(in) :: v_prime_
            !! post-collisional vibrational quantum number
         integer(int32), intent(in) :: j_prime_
            !! post-collisional rotational quantum number
         real(dp), intent(out) :: radial_term_value_
            !! Value of the radial coupling coefficient
         !---------------------------------------------------------------------!
         integer(int32) :: lambda_index, coupling_index
         !---------------------------------------------------------------------!
         lambda_index = find_lambda_index(lambda_)
         !---------------------------------------------------------------------!
         if (lambda_index == 0) then
            call handle_lambda_index_error(lambda_)
            return
         endif
         !---------------------------------------------------------------------!
         coupling_index = find_coupling_index(v_, j_, v_prime_, j_prime_)
         !---------------------------------------------------------------------!
         if (coupling_index == 0) then
            call handle_coupling_index_error(v_, j_, v_prime_, j_prime_)
            return
         endif
         !---------------------------------------------------------------------!
         radial_term_value_ = ISPLINE(intermolecular_distance,                 &
            number_of_r_points, r_grid, coupling_terms(:, lambda_index,        &
            coupling_index), coupling_terms_b_coeffs(:, lambda_index,          &
            coupling_index), coupling_terms_c_coeffs(:, lambda_index,          &
            coupling_index), coupling_terms_d_coeffs(:, lambda_index,          &
            coupling_index))
         !---------------------------------------------------------------------!
      end subroutine get_radial_coupling_term_value