prepare_pes_matrix_elements Subroutine

private subroutine prepare_pes_matrix_elements(channel_indices, channels_omega_values, nonzero_terms_per_element, nonzero_legendre_indices, nonzero_algebraic_coefficients)

-- nonzero_terms_per_element - number of non-vanishing terms in the sum over \(\lambda\) in Eq. 1 in the "Coupling Matrix" section -- nonzero_legendre_indices - corresponding \(\lambda\) value for each non-vanishing coefficient is saved as an index to "legendre_indices" -- nonzero_algebraic_coefficients -- holds all non-vanishing \( g_{{\lambda},\gamma,\gamma'}^{Jp} \) coefficients

Arguments

Type IntentOptional Attributes Name
integer(kind=int32), intent(in) :: channel_indices(:)

holds the indices pointing to the basis arrays

integer(kind=int32), intent(in) :: channels_omega_values(:)

holds all values of \bar{\Omega}

integer(kind=int32), intent(inout) :: 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(inout) :: 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(inout) :: nonzero_algebraic_coefficients(:)

holds the values of the non-zero algebraic coefficients


Calls

proc~~prepare_pes_matrix_elements~~CallsGraph proc~prepare_pes_matrix_elements prepare_pes_matrix_elements proc~process_single_matrix_element process_single_matrix_element proc~prepare_pes_matrix_elements->proc~process_single_matrix_element proc~zero_projections_3j_condition zero_projections_3j_condition proc~process_single_matrix_element->proc~zero_projections_3j_condition proc~percival_seaton_coefficient percival_seaton_coefficient proc~process_single_matrix_element->proc~percival_seaton_coefficient proc~triangle_inequality_holds triangle_inequality_holds proc~zero_projections_3j_condition->proc~triangle_inequality_holds proc~is_sum_even is_sum_even proc~zero_projections_3j_condition->proc~is_sum_even fwig3jj fwig3jj proc~percival_seaton_coefficient->fwig3jj

Called by

proc~~prepare_pes_matrix_elements~~CalledByGraph proc~prepare_pes_matrix_elements prepare_pes_matrix_elements proc~initialize_pes_matrix initialize_pes_matrix proc~initialize_pes_matrix->proc~prepare_pes_matrix_elements program~scattering SCATTERING program~scattering->proc~initialize_pes_matrix

Contents


Variables

Type Visibility Attributes Name Initial
integer(kind=int32), private :: channel_index_1_
integer(kind=int32), private :: channel_index_2_
integer(kind=int32), private :: count_nonzero_algebraic_coefficients
integer(kind=int32), private :: count_nonzero_legendre_terms
integer(kind=int32), private :: count_nonzero_pes_matrix_elements
integer(kind=int32), private :: j_
integer(kind=int32), private :: j_prime_
integer(kind=int32), private :: lambda_
integer(kind=int32), private :: legendre_term_index_
integer(kind=int32), private :: omega_
integer(kind=int32), private :: omega_prime_
real(kind=dp), private :: pscoeff

Source Code

      subroutine prepare_pes_matrix_elements(channel_indices,                  &
         channels_omega_values, nonzero_terms_per_element,                     &
         nonzero_legendre_indices, nonzero_algebraic_coefficients)
         !! prepares:
         !! -- nonzero_terms_per_element - number of non-vanishing terms in
         !!    the sum over \\(\lambda\\) in Eq. 1 in the "Coupling Matrix" section
         !! -- nonzero_legendre_indices - corresponding \\(\lambda\\) value for
         !!    each non-vanishing coefficient is saved as an index to "legendre_indices"
         !! -- nonzero_algebraic_coefficients --  holds _all_ non-vanishing
         !!    \\( g\_{{\lambda},\gamma,\gamma'}^{Jp} \\) coefficients
         !---------------------------------------------------------------------!
         integer(int32), intent(in) :: channel_indices(:)
            !! holds the indices pointing to the basis arrays
         integer(int32), intent(in) :: channels_omega_values(:)
            !! holds all values of \bar{\Omega}
         integer(int32), intent(inout) :: 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(inout) :: 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(inout) :: nonzero_algebraic_coefficients(:)
            !! holds the values of the non-zero algebraic coefficients
         !---------------------------------------------------------------------!
         integer(int32) :: count_nonzero_pes_matrix_elements,                  &
            count_nonzero_algebraic_coefficients, count_nonzero_legendre_terms,&
            j_, j_prime_, omega_, omega_prime_, lambda_, channel_index_1_,     &
            channel_index_2_, legendre_term_index_
         real(dp) :: pscoeff
         !---------------------------------------------------------------------!
         nonzero_terms_per_element        = 0
         nonzero_legendre_indices         = 0
         nonzero_algebraic_coefficients    = 0
         count_nonzero_algebraic_coefficients     = 0
         count_nonzero_pes_matrix_elements  = 0
         !---------------------------------------------------------------------!
         do channel_index_1_ = 1, size(channel_indices)
            j_     = rot_levels(channel_indices(channel_index_1_))
            omega_ = channels_omega_values(channel_index_1_)
            do channel_index_2_ = 1, channel_index_1_
               j_prime_     = rot_levels(channel_indices(channel_index_2_))
               omega_prime_ = channels_omega_values(channel_index_2_)
               if (omega_  /= omega_prime_) cycle
               !---------------------------------------------------------------!
               ! passed \bar{\Omega} = \bar{\Omega}' condition
               !---------------------------------------------------------------!
               count_nonzero_pes_matrix_elements =                             &
                  count_nonzero_pes_matrix_elements + 1
               !---------------------------------------------------------------!
               ! process a single matrix element:
               ! determine non-zero terms in the sum over legendre polynomials
               ! for this element; these are saved to ...
               !---------------------------------------------------------------!
               call process_single_matrix_element(j_, j_prime_, omega_,        &
                  count_nonzero_algebraic_coefficients,                        &
                  count_nonzero_legendre_terms, nonzero_legendre_indices,      &
                  nonzero_algebraic_coefficients)
               
               nonzero_terms_per_element(count_nonzero_pes_matrix_elements)&
                  = count_nonzero_legendre_terms
            enddo
         enddo
         !---------------------------------------------------------------------!
      end subroutine prepare_pes_matrix_elements