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
Type | Intent | Optional | 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 |
matirx element of the interaction potential contribution to the coupling matrix
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_ |
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