Renormalized Numerov propagator
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer(kind=int32), | intent(in) | :: | number_of_channels_ |
size of the basis |
||
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-zero terms in the sum (Eq. (6.21)) for each non-zero element of W/V |
||
integer(kind=int32), | intent(in) | :: | nonzero_legendre_indices_(:) |
holds the proper indices pointing to l1/l2/lltabs, which correspond to the non-vanishing elements of the sum (Eq. (6.21)) for each non-zero element of W/V |
||
real(kind=dp), | intent(in) | :: | nonzero_algebraic_coefficients_(:) |
holds the values of the non-zero algebraic coefficients |
||
integer(kind=int32), | intent(in) | :: | number_of_steps_ |
number of steps from r_min to r_max |
||
integer(kind=int32), | intent(in) | :: | total_angular_momentum_ |
total angular momentum |
||
real(kind=dp), | intent(inout) | :: | log_der_matrix_(number_of_channels_,number_of_channels_) |
resulting log-derivative matrix at r_max |
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
real(kind=dp), | private | :: | calculation_time_ | ||||
real(kind=dp), | private, | dimension(number_of_channels_,number_of_channels_) | :: | centrifugal_matrix_ | |||
integer(kind=int32), | private | :: | channel_index_1_ | ||||
integer(kind=int32), | private | :: | channel_index_2_ | ||||
real(kind=dp), | private | :: | finish | ||||
integer(kind=int32), | private | :: | i | ||||
real(kind=dp), | private | :: | intermolecular_distance_ | ||||
real(kind=dp), | private, | dimension(number_of_channels_,number_of_channels_) | :: | r_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_) | :: | r_matrix_r_max_ | |||
real(kind=dp), | private | :: | start | ||||
real(kind=dp), | private | :: | step_numerov_ | ||||
real(kind=dp), | private, | dimension(number_of_channels_,number_of_channels_) | :: | t_matrix_ | |||
real(kind=dp), | private, | dimension(number_of_channels_,number_of_channels_) | :: | t_matrix_minus_ | |||
real(kind=dp), | private, | dimension(number_of_channels_,number_of_channels_) | :: | t_matrix_plus_ |
subroutine numerov(number_of_channels_, channel_indices_, &
channels_omega_values_, nonzero_terms_per_element_, &
nonzero_legendre_indices_, nonzero_algebraic_coefficients_, &
number_of_steps_, total_angular_momentum_, log_der_matrix_)
!! Renormalized Numerov propagator
!---------------------------------------------------------------------!
integer(int32), intent(in) :: number_of_channels_
!! size of the basis
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-zero terms in the sum (Eq. (6.21)) for
!! each non-zero element of W/V
integer(int32), intent(in) :: nonzero_legendre_indices_(:)
!! holds the proper indices pointing to l1/l2/lltabs, which
!! correspond to the non-vanishing elements of the sum (Eq. (6.21))
!! for each non-zero element of W/V
real(dp), intent(in) :: nonzero_algebraic_coefficients_(:)
!! holds the values of the non-zero algebraic coefficients
integer(int32), intent(in) :: number_of_steps_
!! number of steps from r_min to r_max
integer(int32), intent(in) :: total_angular_momentum_
!! total angular momentum
real(dp), intent(inout) :: log_der_matrix_(number_of_channels_,number_of_channels_)
!! resulting log-derivative matrix at r_max
!---------------------------------------------------------------------!
integer(int32) :: i, channel_index_1_, channel_index_2_
real(dp) :: start, finish, intermolecular_distance_, step_numerov_, &
calculation_time_
real(dp), dimension(number_of_channels_,number_of_channels_) :: &
centrifugal_matrix_, &
t_matrix_minus_, t_matrix_, t_matrix_plus_, r_matrix_, &
r_matrix_r_max_, r_matrix_plus_
!---------------------------------------------------------------------!
CALL CPU_TIME(start)
step_numerov_ = (r_max - r_min)/dble(number_of_steps_ - 1)
intermolecular_distance_ = r_min
!---------------------------------------------------------------------!
! Initial setup: calculate centrifugal matrix and R_matrix at r_min + 1
!---------------------------------------------------------------------!
call initial_setup(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_)
!---------------------------------------------------------------------!
! Propagation loop
!---------------------------------------------------------------------!
do i=2, number_of_steps_ - 2
!------------------------------------------------------------------!
intermolecular_distance_ = r_min + (i-1)*step_numerov_
!------------------------------------------------------------------!
call 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_)
!------------------------------------------------------------------!
end do
!---------------------------------------------------------------------!
call handle_final_propagation_steps(number_of_channels_, &
step_numerov_, total_angular_momentum_, channel_indices_, &
channels_omega_values_, nonzero_terms_per_element_, &
nonzero_legendre_indices_, nonzero_algebraic_coefficients_, &
centrifugal_matrix_, r_matrix_, t_matrix_minus_, t_matrix_, &
t_matrix_plus_, r_matrix_r_max_, r_matrix_plus_)
!---------------------------------------------------------------------!
CALL CPU_TIME(finish)
!---------------------------------------------------------------------!
! Eq. (6.29)
!---------------------------------------------------------------------!
call calculate_log_der_matrix(step_numerov_,number_of_channels_, &
t_matrix_minus_,t_matrix_,t_matrix_plus_,r_matrix_r_max_, &
r_matrix_plus_,log_der_matrix_)
!---------------------------------------------------------------------!
call propagator_summary(r_min, r_max, number_of_steps_)
!---------------------------------------------------------------------!
if (print_level.ge.2) then
call time_count_summary(start, finish, calculation_time_, &
"Propagation completed in ")
endif
!---------------------------------------------------------------------!
end subroutine numerov