general_propagation_step Subroutine

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

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

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


Calls

proc~~general_propagation_step~~CallsGraph proc~general_propagation_step 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 proc~float_to_character float_to_character proc~ispline->proc~float_to_character proc~write_warning write_warning proc~ispline->proc~write_warning proc~handle_lambda_index_error->proc~write_error proc~integer_to_character integer_to_character proc~handle_lambda_index_error->proc~integer_to_character proc~handle_coupling_index_error->proc~write_error proc~handle_coupling_index_error->proc~integer_to_character proc~write_message write_message proc~write_error->proc~write_message

Called by

proc~~general_propagation_step~~CalledByGraph proc~general_propagation_step general_propagation_step proc~numerov numerov proc~numerov->proc~general_propagation_step proc~handle_final_propagation_steps handle_final_propagation_steps proc~numerov->proc~handle_final_propagation_steps proc~handle_final_propagation_steps->proc~general_propagation_step 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, 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_

Source Code

      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