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
Type | Intent | Optional | 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 |
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_ |
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