calculate_pes_matrix Subroutine

public subroutine calculate_pes_matrix(total_angular_momentum_, intermolecular_distance_, channel_indices_, channels_omega_values_, nonzero_terms_per_element_, nonzero_legendre_indices_, nonzero_algebraic_coefficients_, vmatrix)

calculates the contribution to the coupling matrix from the the interaction potential (PES); see 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
integer(kind=int32), intent(in) :: total_angular_momentum_

total angular momentum

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

intermolecular distance

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

holds the indices pointing to the basis arrays

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

holds all values of \bar{\Omega}

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

keeps the number of non-vanishing elements of the sum over \(\lambda\) for each non-zero element of the pes matrix

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

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

(output) - the interaction potential contribution to the pes matrix


Calls

proc~~calculate_pes_matrix~~CallsGraph proc~calculate_pes_matrix calculate_pes_matrix proc~calculate_single_pes_matrix_element calculate_single_pes_matrix_element proc~calculate_pes_matrix->proc~calculate_single_pes_matrix_element interface~fill_symmetric_matrix fill_symmetric_matrix proc~calculate_pes_matrix->interface~fill_symmetric_matrix 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_pes_matrix~~CalledByGraph proc~calculate_pes_matrix calculate_pes_matrix 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

Source Code


Variables

Type Visibility Attributes Name Initial
integer(kind=int32), private :: channel_index_1_
integer(kind=int32), private :: channel_index_2_
integer(kind=int32), private :: count_nonzero_algebraic_coefficients_
integer(kind=int32), private :: count_nonzero_coupling_matrix_elements
integer(kind=int32), private :: count_nonzero_legendre_terms
integer(kind=int32), private :: omega_
integer(kind=int32), private :: omega_prime_

Source Code

      subroutine calculate_pes_matrix(total_angular_momentum_,                 &
         intermolecular_distance_, channel_indices_, channels_omega_values_,   &
         nonzero_terms_per_element_, nonzero_legendre_indices_,                &
         nonzero_algebraic_coefficients_, vmatrix)
         !! calculates the contribution to the coupling matrix
         !! from the the interaction potential (PES);
         !! see 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
         !---------------------------------------------------------------------!
         integer(int32), intent(in) :: total_angular_momentum_
            !! total angular momentum
         real(dp), intent(in) :: intermolecular_distance_
            !! intermolecular distance
         integer(int32), intent(in) :: channel_indices_(:)
            !! holds the indices pointing to the basis arrays
         integer(int32), intent(in) :: channels_omega_values_(:)
            !! holds all values of \bar{\Omega}
         integer(int32), intent(in) :: nonzero_terms_per_element_(:)
            !! keeps the number of non-vanishing elements of the sum
            !! over \\(\lambda\\) for each non-zero element of the pes matrix
         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
         real(dp), intent(out) :: vmatrix(:,:)
            !! (output) - the interaction potential contribution
            !! to the pes matrix
         !---------------------------------------------------------------------!
         integer(int32) :: count_nonzero_algebraic_coefficients_,              &
            count_nonzero_coupling_matrix_elements,                            &
            count_nonzero_legendre_terms, channel_index_1_, channel_index_2_,  &
            omega_, omega_prime_
         !---------------------------------------------------------------------!
         vmatrix   = 0
         count_nonzero_algebraic_coefficients_    = 0
         count_nonzero_coupling_matrix_elements = 0
         !---------------------------------------------------------------------!
         do channel_index_1_ = 1, size(channel_indices_)
            omega_ = channels_omega_values_(channel_index_1_)
            do channel_index_2_ = 1, channel_index_1_
               omega_prime_ = channels_omega_values_(channel_index_2_)
               !---------------------------------------------------------------!
               if (omega_ /= omega_prime_) cycle
               !---------------------------------------------------------------!
               count_nonzero_coupling_matrix_elements                          &
                  = count_nonzero_coupling_matrix_elements + 1
               !---------------------------------------------------------------!
               ! process a single matrix element:
               ! get number of  non-zero terms in Legendre expansion for this
               ! matrix element
               !---------------------------------------------------------------!
               count_nonzero_legendre_terms                                    &
                  = nonzero_terms_per_element_(count_nonzero_coupling_matrix_elements)
               !---------------------------------------------------------------!
               ! implementation of Eq. 1 in "Coupling Matrix" section
               !---------------------------------------------------------------!
               vmatrix(channel_index_1_, channel_index_2_)                     &
                  = 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_)
               !---------------------------------------------------------------!
             enddo
         enddo
         !---------------------------------------------------------------------!
         call fill_symmetric_matrix(vmatrix,'u')
         !---------------------------------------------------------------------!
      end subroutine calculate_pes_matrix