calculate_offdiagonal_centrifugal_element Function

private function calculate_offdiagonal_centrifugal_element(total_angular_momentum_, j_, omega_, omega_prime_, delta_1_, delta_2_) result(offdiagonal_element_)

calculates off-diagonal element of the centrifgual matrix, see Eq. 5 in "Coupling Matrix" section

Arguments

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

total angular momentum

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

rotational angular momentum

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

\(\bar{\Omega}\)

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

\(\bar{\Omega}'\)

real(kind=dp), intent(in) :: delta_1_

Kronecker delta functions determined earlier

real(kind=dp), intent(in) :: delta_2_

Kronecker delta functions determined earlier

Return Value real(kind=dp)

(output) value of the off-diagonal element of the centrifgual matrix


Called by

proc~~calculate_offdiagonal_centrifugal_element~~CalledByGraph proc~calculate_offdiagonal_centrifugal_element calculate_offdiagonal_centrifugal_element proc~calculate_centrifugal_matrix calculate_centrifugal_matrix proc~calculate_centrifugal_matrix->proc~calculate_offdiagonal_centrifugal_element proc~initial_setup initial_setup proc~initial_setup->proc~calculate_centrifugal_matrix proc~numerov numerov proc~numerov->proc~initial_setup program~scattering SCATTERING program~scattering->proc~numerov

Contents


Source Code

      function calculate_offdiagonal_centrifugal_element(                      &
         total_angular_momentum_, j_, omega_, omega_prime_, delta_1_, delta_2_)&
         result(offdiagonal_element_)
         !! calculates off-diagonal element of the centrifgual matrix, see
         !! Eq. 5 in "Coupling Matrix" section
         !---------------------------------------------------------------------!
         integer(int32), intent(in) :: total_angular_momentum_
            !! total angular momentum
         integer(int32), intent(in) :: j_
            !! rotational angular momentum
         integer(int32), intent(in) :: omega_
            !! \\(\bar{\Omega}\\)
         integer(int32), intent(in) :: omega_prime_
            !! \\(\bar{\Omega}'\\)
         real(dp), intent(in) :: delta_1_, delta_2_
            !! Kronecker delta functions determined earlier
         real(dp) :: offdiagonal_element_
            !! (output) value of the off-diagonal element of
            !! the centrifgual matrix
         !---------------------------------------------------------------------!
         offdiagonal_element_ = - sqrt(real(                                   &
            (total_angular_momentum_ * (total_angular_momentum_ + 1)           &
            - omega_ * omega_prime_) * (j_ * (j_ + 1) - omega_ * omega_prime_),&
            dp)) * sqrt( (1.0_dp + delta_1_) * (1.0_dp + delta_2_) )
         !---------------------------------------------------------------------!
      end function calculate_offdiagonal_centrifugal_element