calculate_single_pes_matrix_element Function

private function calculate_single_pes_matrix_element(intermolecular_distance_, channel_index_1_, channel_index_2_, channel_indices_, count_nonzero_legendre_terms, count_nonzero_algebraic_coefficients_, nonzero_legendre_indices_, nonzero_algebraic_coefficients_) result(matrix_element_)

Implementation of Eq. 1 in "Coupling Matrix" section; diagonal contribution from wavevectors (see the last term in Eq. 3 of "What are coupled equations" section) is added

Arguments

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

intermolecular distance

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

indices identifying matrix element

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

indices identifying matrix element

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

holds the indices pointing to the basis arrays

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

number of non-zero terms in Legendre expansion for a single matrix element of the interaction potential

integer(kind=int32), intent(inout) :: count_nonzero_algebraic_coefficients_

running index counting non-zero coupling coefficients, \( g_{{\lambda},\gamma,\gamma'}^{Jp} \) in the whole matrix; incremented in this subroutine

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

holds indices pointing to legendre_indices, which correspond to the non-vanishing elements of the sum over \(\lambda\) for each non-zero element of the pes matrix;

real(kind=dp), intent(in) :: nonzero_algebraic_coefficients_(:)

holds the values of the non-zero algebraic coefficients

Return Value real(kind=dp)

matirx element of the interaction potential contribution to the coupling matrix


Calls

proc~~calculate_single_pes_matrix_element~~CallsGraph proc~calculate_single_pes_matrix_element calculate_single_pes_matrix_element proc~get_radial_coupling_term_value get_radial_coupling_term_value proc~calculate_single_pes_matrix_element->proc~get_radial_coupling_term_value proc~wavevector_squared_from_energy wavevector_squared_from_energy proc~calculate_single_pes_matrix_element->proc~wavevector_squared_from_energy 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~total_energy total_energy proc~wavevector_squared_from_energy->proc~total_energy proc~write_error write_error proc~wavevector_squared_from_energy->proc~write_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~handle_lambda_index_error->proc~write_error proc~integer_to_character integer_to_character proc~handle_lambda_index_error->proc~integer_to_character proc~handle_coupling_index_error->proc~write_error proc~handle_coupling_index_error->proc~integer_to_character proc~write_message write_message proc~write_error->proc~write_message proc~write_warning->proc~write_message

Called by

proc~~calculate_single_pes_matrix_element~~CalledByGraph proc~calculate_single_pes_matrix_element calculate_single_pes_matrix_element 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
real(kind=dp), private :: algebraic_coefficient_
real(kind=dp), private :: internal_energy_
integer(kind=int32), private :: j_
integer(kind=int32), private :: j_prime_
integer(kind=int32), private :: lambda_
integer(kind=int32), private :: lambda_index_
real(kind=dp), private :: radial_term_
real(kind=dp), private :: sum_over_lambda_
integer(kind=int32), private :: v_
integer(kind=int32), private :: v_prime_

Source Code

      function calculate_single_pes_matrix_element(intermolecular_distance_,   &
         channel_index_1_, channel_index_2_, channel_indices_,                 &
         count_nonzero_legendre_terms, count_nonzero_algebraic_coefficients_,  &
         nonzero_legendre_indices_, nonzero_algebraic_coefficients_)           &
         result(matrix_element_)
         !! Implementation of Eq. 1 in "Coupling Matrix" section;
         !! diagonal contribution from wavevectors (see the last term in
         !! Eq. 3 of "What are coupled equations" section) is added
         !---------------------------------------------------------------------!
         real(dp), intent(in) :: intermolecular_distance_
            !! intermolecular distance
         integer(int32), intent(in) :: channel_index_1_, channel_index_2_
            !! indices identifying matrix element
         integer(int32), intent(in) :: channel_indices_(:)
            !! holds the indices pointing to the basis arrays
         integer(int32), intent(in) :: count_nonzero_legendre_terms
            !! number of non-zero terms in Legendre expansion for a single
            !! matrix element of the interaction potential
         integer(int32), intent(in) :: nonzero_legendre_indices_(:)
            !! holds indices pointing to legendre_indices, which correspond to
            !! the non-vanishing elements of the sum over \\(\lambda\\)
            !! for each non-zero element of the pes matrix;
         real(dp), intent(in) :: nonzero_algebraic_coefficients_(:)
            !! holds the values of the non-zero algebraic coefficients
         integer(int32), intent(inout) :: count_nonzero_algebraic_coefficients_
            !! running index counting non-zero coupling coefficients,
            !! \\( g\_{{\lambda},\gamma,\gamma'}^{Jp} \\) in the whole matrix;
            !! incremented in this subroutine
         real(dp) :: matrix_element_
            !! matirx element of the interaction potential contribution
            !! to the coupling matrix
         !---------------------------------------------------------------------!
         integer(int32) :: v_, j_, v_prime_, j_prime_, lambda_, lambda_index_
         real(dp) :: internal_energy_, sum_over_lambda_,                       &
            algebraic_coefficient_, radial_term_
         !---------------------------------------------------------------------!
         v_ = vib_levels(channel_indices_(channel_index_1_))
         j_ = rot_levels(channel_indices_(channel_index_1_))
         v_prime_ = vib_levels(channel_indices_(channel_index_2_))
         j_prime_ = rot_levels(channel_indices_(channel_index_2_))
         internal_energy_ = internal_energies(channel_indices_(channel_index_1_))
         !---------------------------------------------------------------------!
         sum_over_lambda_ = 0.0_dp
         do lambda_index_ = 1, count_nonzero_legendre_terms
            !------------------------------------------------------------------!
            count_nonzero_algebraic_coefficients_                              &
               = count_nonzero_algebraic_coefficients_ + 1
            !------------------------------------------------------------------!
            lambda_ = legendre_indices(nonzero_legendre_indices_(              &
               count_nonzero_algebraic_coefficients_))
            algebraic_coefficient_ = nonzero_algebraic_coefficients_(          &
               count_nonzero_algebraic_coefficients_)
            !------------------------------------------------------------------!
            call get_radial_coupling_term_value(intermolecular_distance_,      &
               lambda_, v_, j_, v_prime_, j_prime_, radial_term_)
            !------------------------------------------------------------------!
            sum_over_lambda_ = sum_over_lambda_                                &
               + algebraic_coefficient_ * radial_term_
            !------------------------------------------------------------------!
         enddo
         matrix_element_ =  2.0_dp * reduced_mass * sum_over_lambda_
         !---------------------------------------------------------------------!
         if (channel_index_1_ == channel_index_2_) then
            matrix_element_ = matrix_element_                                  &
               - wavevector_squared_from_energy(internal_energy_)
         endif
         !---------------------------------------------------------------------!
      end function calculate_single_pes_matrix_element