handle_final_propagation_steps Subroutine

private subroutine 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_)

Handles propagation at the last two grid points: R_{N-1} and R_{N}: provides T-matrix at N-1, N and N+1 points and the Ratio matrix at N and N+1 points

Arguments

Type IntentOptional 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

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} \cdot\) centrifugal matrix -

real(kind=dp), intent(inout) :: r_matrix_(number_of_channels_,number_of_channels_)

Ratio matrix at N-1 step

real(kind=dp), intent(out) :: t_matrix_minus_(number_of_channels_,number_of_channels_)

T-matrix at N-1 step

real(kind=dp), intent(out) :: t_matrix_(number_of_channels_,number_of_channels_)

T-matrix at N step

real(kind=dp), intent(out) :: t_matrix_plus_(number_of_channels_,number_of_channels_)

T-matrix at N+1 step

real(kind=dp), intent(out) :: r_matrix_r_max_(number_of_channels_,number_of_channels_)

Ratio matrix at N step

real(kind=dp), intent(out) :: r_matrix_plus_(number_of_channels_,number_of_channels_)

Ratio matrix at N+1 step


Calls

proc~~handle_final_propagation_steps~~CallsGraph proc~handle_final_propagation_steps handle_final_propagation_steps proc~general_propagation_step general_propagation_step proc~handle_final_propagation_steps->proc~general_propagation_step proc~calculate_pes_matrix calculate_pes_matrix proc~general_propagation_step->proc~calculate_pes_matrix interface~invert_symmetric_matrix invert_symmetric_matrix proc~general_propagation_step->interface~invert_symmetric_matrix interface~fill_symmetric_matrix fill_symmetric_matrix proc~general_propagation_step->interface~fill_symmetric_matrix proc~calculate_u_matrix calculate_u_matrix proc~general_propagation_step->proc~calculate_u_matrix proc~calculate_coupling_matrix calculate_coupling_matrix proc~general_propagation_step->proc~calculate_coupling_matrix proc~calculate_t_matrix calculate_t_matrix proc~general_propagation_step->proc~calculate_t_matrix proc~calculate_pes_matrix->interface~fill_symmetric_matrix proc~calculate_single_pes_matrix_element calculate_single_pes_matrix_element proc~calculate_pes_matrix->proc~calculate_single_pes_matrix_element proc~calculate_u_matrix->interface~invert_symmetric_matrix proc~calculate_u_matrix->interface~fill_symmetric_matrix interface~add_scalar_to_diagonal add_scalar_to_diagonal proc~calculate_u_matrix->interface~add_scalar_to_diagonal 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

Called by

proc~~handle_final_propagation_steps~~CalledByGraph proc~handle_final_propagation_steps handle_final_propagation_steps proc~numerov numerov proc~numerov->proc~handle_final_propagation_steps program~scattering SCATTERING program~scattering->proc~numerov

Contents


Variables

Type Visibility Attributes Name Initial
real(kind=dp), private, dimension(number_of_channels_,number_of_channels_) :: coupling_matrix_
real(kind=dp), private :: intermolecular_distance_
logical, private :: is_t_matrix_required_
real(kind=dp), private, dimension(number_of_channels_,number_of_channels_) :: pes_matrix_
real(kind=dp), private, dimension(number_of_channels_,number_of_channels_) :: u_matrix_

Source Code

      subroutine 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_)
         !! Handles propagation at the last two grid points:
         !! R_{N-1} and R_{N}: provides T-matrix at N-1, N and N+1 points
         !! and the Ratio matrix at N and N+1 points
         !---------------------------------------------------------------------!
         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
         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} \cdot\\) centrifugal matrix -
         real(dp), intent(inout) :: r_matrix_(number_of_channels_,number_of_channels_)
            !! Ratio matrix at N-1 step
         real(dp), intent(out) :: t_matrix_minus_(number_of_channels_,number_of_channels_)
            !! T-matrix at N-1 step
         real(dp), intent(out) :: t_matrix_(number_of_channels_,number_of_channels_)
            !! T-matrix at N step
         real(dp), intent(out) :: t_matrix_plus_(number_of_channels_,number_of_channels_)
            !! T-matrix at N+1 step
         real(dp), intent(out) :: r_matrix_r_max_(number_of_channels_,number_of_channels_)
            !! Ratio matrix at N step
         real(dp), intent(out) :: r_matrix_plus_(number_of_channels_,number_of_channels_)
            !! Ratio matrix at N+1 step
         !---------------------------------------------------------------------!
         logical :: is_t_matrix_required_
         real(dp) :: intermolecular_distance_
         real(dp), dimension(number_of_channels_,number_of_channels_) ::       &
            pes_matrix_, coupling_matrix_, u_matrix_
         !---------------------------------------------------------------------!
         is_t_matrix_required_ = .true.
         !---------------------------------------------------------------------!
         ! N - 1 step
         !---------------------------------------------------------------------!
         intermolecular_distance_ = r_max - 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_,   &
            is_t_matrix_required_, t_matrix_minus_)
         r_matrix_r_max_ = r_matrix_
         !---------------------------------------------------------------------!
         ! N  step
         !---------------------------------------------------------------------!
         intermolecular_distance_ = intermolecular_distance_ + 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_,   &
            is_t_matrix_required_, t_matrix_)
         r_matrix_plus_ = r_matrix_
         !---------------------------------------------------------------------!
         ! N + 1 step
         !---------------------------------------------------------------------!
         intermolecular_distance_ = intermolecular_distance_ + 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_,   &
            is_t_matrix_required_, t_matrix_plus_)
         !---------------------------------------------------------------------!
      end subroutine handle_final_propagation_steps