Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer(kind=int32), | intent(in) | :: | number_of_channels_ |
size of the basis |
||
real(kind=dp), | intent(in) | :: | step_numerov_ |
step of the propagator |
||
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_(number_of_channels_) |
holds the indices pointing to the basis arrays |
||
integer(kind=int32), | intent(in) | :: | channels_omega_values_(number_of_channels_) |
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(in) | :: | centrifugal_matrix_(number_of_channels_,number_of_channels_) |
(R*2)centrifugal matrix - |
||
real(kind=dp), | intent(inout) | :: | r_matrix_(number_of_channels_,number_of_channels_) |
on input: R-matrix at previous step on output: R-matrix at next step |
||
logical, | intent(in), | optional | :: | is_t_matrix_required_ | ||
real(kind=dp), | intent(out), | optional | :: | t_matrix_returned_(number_of_channels_,number_of_channels_) |
on input: R-matrix at previous step on output: R-matrix at next step |
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
real(kind=dp), | private, | dimension(number_of_channels_,number_of_channels_) | :: | coupling_matrix_ | |||
real(kind=dp), | private, | dimension(number_of_channels_,number_of_channels_) | :: | pes_matrix_ | |||
real(kind=dp), | private, | dimension(number_of_channels_,number_of_channels_) | :: | r_matrix_plus_ | |||
real(kind=dp), | private, | dimension(number_of_channels_,number_of_channels_) | :: | t_matrix_ | |||
real(kind=dp), | private, | dimension(number_of_channels_,number_of_channels_) | :: | u_matrix_ |
subroutine general_propagation_step(number_of_channels_, step_numerov_, &
total_angular_momentum_, intermolecular_distance_, channel_indices_, &
channels_omega_values_, nonzero_terms_per_element_, &
nonzero_legendre_indices_, nonzero_algebraic_coefficients_, &
centrifugal_matrix_, r_matrix_, is_t_matrix_required_, t_matrix_returned_)
!! ...
!---------------------------------------------------------------------!
integer(int32), intent(in) :: number_of_channels_
!! size of the basis
real(dp), intent(in) :: step_numerov_
!! step of the propagator
integer(int32), intent(in) :: total_angular_momentum_
!! total angular momentum
real(dp), intent(in) :: intermolecular_distance_
!! intermolecular distance
integer(int32), intent(in) :: channel_indices_(number_of_channels_)
!! holds the indices pointing to the basis arrays
integer(int32), intent(in) :: channels_omega_values_(number_of_channels_)
!! 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(in) :: centrifugal_matrix_(number_of_channels_,number_of_channels_)
!! (R**2)*centrifugal matrix -
real(dp), intent(inout) :: r_matrix_(number_of_channels_,number_of_channels_)
!! on input: R-matrix at previous step
!! on output: R-matrix at next step
logical, intent(in), optional :: is_t_matrix_required_
!! ...
real(dp), intent(out), optional :: t_matrix_returned_(number_of_channels_,number_of_channels_)
!! on input: R-matrix at previous step
!! on output: R-matrix at next step
!---------------------------------------------------------------------!
real(dp), dimension(number_of_channels_,number_of_channels_) :: &
pes_matrix_, coupling_matrix_, t_matrix_, u_matrix_, r_matrix_plus_
!---------------------------------------------------------------------!
! Calculate PES matrix at R
!---------------------------------------------------------------------!
call calculate_pes_matrix(total_angular_momentum_, &
intermolecular_distance_, channel_indices_, &
channels_omega_values_, nonzero_terms_per_element_, &
nonzero_legendre_indices_, nonzero_algebraic_coefficients_, &
pes_matrix_)
!---------------------------------------------------------------------!
! Merge centrifugal and PES matrix into Coupling matrix
!---------------------------------------------------------------------!
call calculate_coupling_matrix(intermolecular_distance_, pes_matrix_, &
centrifugal_matrix_, coupling_matrix_)
!---------------------------------------------------------------------!
! Calculate T-matrix and U-matrix
!---------------------------------------------------------------------!
call calculate_t_matrix(step_numerov_, coupling_matrix_, t_matrix_)
call calculate_u_matrix(t_matrix_, u_matrix_)
!---------------------------------------------------------------------!
! Invert R matrix from previous step
!---------------------------------------------------------------------!
call invert_symmetric_matrix(r_matrix_)
call fill_symmetric_matrix(r_matrix_, 'u')
!---------------------------------------------------------------------!
! R_{n+1} = U_{n} - R_{n}^{-1}
!---------------------------------------------------------------------!
r_matrix_plus_ = u_matrix_ - r_matrix_
!---------------------------------------------------------------------!
! Move R_{n+1} to R_{n}
!---------------------------------------------------------------------!
r_matrix_ = r_matrix_plus_
!---------------------------------------------------------------------!
if (present(is_t_matrix_required_) .and. is_t_matrix_required_) then
t_matrix_returned_ = t_matrix_
endif
!------------------------------------------------------------------------!
end subroutine general_propagation_step