initial_setup Subroutine

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

Initial setup of the propagator: call centrifugal matrix (kept throughotu the propagation) and other matrices at \(R_{\mathrm{min}}\)

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(inout) :: centrifugal_matrix_(number_of_channels_,number_of_channels_)

(R*2)centrifugal matrix - calculated once, will be used throughout the propagation

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

R-matrix at r_min


Calls

proc~~initial_setup~~CallsGraph proc~initial_setup initial_setup proc~calculate_pes_matrix calculate_pes_matrix proc~initial_setup->proc~calculate_pes_matrix proc~calculate_u_matrix calculate_u_matrix proc~initial_setup->proc~calculate_u_matrix proc~calculate_centrifugal_matrix calculate_centrifugal_matrix proc~initial_setup->proc~calculate_centrifugal_matrix proc~calculate_coupling_matrix calculate_coupling_matrix proc~initial_setup->proc~calculate_coupling_matrix proc~calculate_t_matrix calculate_t_matrix proc~initial_setup->proc~calculate_t_matrix proc~calculate_single_pes_matrix_element calculate_single_pes_matrix_element proc~calculate_pes_matrix->proc~calculate_single_pes_matrix_element interface~fill_symmetric_matrix fill_symmetric_matrix proc~calculate_pes_matrix->interface~fill_symmetric_matrix interface~invert_symmetric_matrix invert_symmetric_matrix 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~calculate_centrifugal_matrix->interface~fill_symmetric_matrix proc~delta_for_zero_omega delta_for_zero_omega proc~calculate_centrifugal_matrix->proc~delta_for_zero_omega proc~calculate_diagonal_centrifugal_element calculate_diagonal_centrifugal_element proc~calculate_centrifugal_matrix->proc~calculate_diagonal_centrifugal_element proc~calculate_offdiagonal_centrifugal_element calculate_offdiagonal_centrifugal_element proc~calculate_centrifugal_matrix->proc~calculate_offdiagonal_centrifugal_element 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~~initial_setup~~CalledByGraph proc~initial_setup initial_setup proc~numerov numerov proc~numerov->proc~initial_setup program~scattering SCATTERING program~scattering->proc~numerov

Contents

Source Code


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_) :: t_matrix_
real(kind=dp), private, dimension(number_of_channels_,number_of_channels_) :: u_matrix_

Source Code

      subroutine 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_)
         !! Initial setup of the propagator: call centrifugal matrix
         !! (kept throughotu the propagation) and other matrices
         !! at \\(R_{\mathrm{min}}\\)
         !---------------------------------------------------------------------!
         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(inout) :: centrifugal_matrix_(number_of_channels_,number_of_channels_)
            !! (R**2)*centrifugal matrix - calculated once,
            !! will be used throughout the propagation
         real(dp), intent(inout) :: r_matrix_(number_of_channels_,number_of_channels_)
            !! R-matrix at r_min
         !---------------------------------------------------------------------!   
         real(dp), dimension(number_of_channels_,number_of_channels_) ::       &
            pes_matrix_, coupling_matrix_, t_matrix_, u_matrix_
         !---------------------------------------------------------------------!   
         ! Calculate centrifugal matrix
         !---------------------------------------------------------------------!   
         call calculate_centrifugal_matrix(total_angular_momentum_,            &
            channel_indices_, channels_omega_values_, centrifugal_matrix_)
         !---------------------------------------------------------------------!   
         ! Calculate PES matrix at r_min
         !---------------------------------------------------------------------!   
         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 initial T-matrix and U-matrix
         !---------------------------------------------------------------------!
         call calculate_t_matrix(step_numerov_, coupling_matrix_, t_matrix_)
         call calculate_u_matrix(t_matrix_, u_matrix_)
         !---------------------------------------------------------------------!
         ! Initialize R-matrix: R-matrix at r_min + 1 = U-matrix at r_min
         !---------------------------------------------------------------------!
         r_matrix_ = u_matrix_
         !---------------------------------------------------------------------!
      end subroutine initial_setup