calculate_centrifugal_matrix Subroutine

public subroutine calculate_centrifugal_matrix(total_angular_momentum_, channel_indices_, channels_omega_values_, centrifugal_matrix_)

calculates the (R*2)centrifugal matrix from the second term of Eq. 3 in "What are coupled equations?" section; Matrix elements are given in Eq. 4 and 6 of "Coupling Matrix" secion

Arguments

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

total angular momentum

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}

real(kind=dp), intent(out) :: centrifugal_matrix_(:,:)

(output) - (R*2)centrifugal matrix


Calls

proc~~calculate_centrifugal_matrix~~CallsGraph proc~calculate_centrifugal_matrix calculate_centrifugal_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 interface~fill_symmetric_matrix fill_symmetric_matrix proc~calculate_centrifugal_matrix->interface~fill_symmetric_matrix proc~calculate_offdiagonal_centrifugal_element calculate_offdiagonal_centrifugal_element proc~calculate_centrifugal_matrix->proc~calculate_offdiagonal_centrifugal_element

Called by

proc~~calculate_centrifugal_matrix~~CalledByGraph proc~calculate_centrifugal_matrix calculate_centrifugal_matrix 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


Variables

Type Visibility Attributes Name Initial
real(kind=dp), private :: centtmp
integer(kind=int32), private :: channel_index_1_
integer(kind=int32), private :: channel_index_2_
real(kind=dp), private :: delta_1_
real(kind=dp), private :: delta_2_
integer(kind=int32), private :: j_
integer(kind=int32), private :: j_prime_
integer(kind=int32), private :: omega_
integer(kind=int32), private :: omega_prime_
integer(kind=int32), private :: v_
integer(kind=int32), private :: v_prime_

Source Code

      subroutine calculate_centrifugal_matrix(total_angular_momentum_,         &
         channel_indices_, channels_omega_values_, centrifugal_matrix_)    
         !! calculates the (R**2)*centrifugal matrix from the second term
         !! of Eq. 3 in "What are coupled equations?" section;
         !! Matrix elements are given in Eq. 4 and 6 of "Coupling Matrix" secion
         !---------------------------------------------------------------------!
         integer(int32), intent(in) :: total_angular_momentum_
            !! total angular momentum
         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}
         real(dp), intent(out) :: centrifugal_matrix_(:,:)
            !! (output) - (R**2)*centrifugal matrix
         !---------------------------------------------------------------------!
         integer(int32) :: omega_, omega_prime_, v_, j_, v_prime_, j_prime_,  &
            channel_index_1_, channel_index_2_
         real(dp) :: centtmp, delta_1_, delta_2_
         !---------------------------------------------------------------------!
         centrifugal_matrix_  = 0

         do channel_index_1_ = 1, size(channel_indices_)
            v_ = vib_levels(channel_indices_(channel_index_1_))
            j_ = rot_levels(channel_indices_(channel_index_1_))
            omega_ = channels_omega_values_(channel_index_1_)
            delta_1_ = delta_for_zero_omega(omega_)
            do channel_index_2_ = 1, channel_index_1_
               v_prime_ = vib_levels(channel_indices_(channel_index_2_))
               j_prime_ = rot_levels(channel_indices_(channel_index_2_))
               omega_prime_ = channels_omega_values_(channel_index_2_)
               delta_2_ = delta_for_zero_omega(omega_prime_)
               !---------------------------------------------------------------!
               if (v_ /= v_prime_ .or. j_ /= j_prime_ .or.                     &
                  abs(omega_ - omega_prime_) > 1) then
                  cycle
               endif
               !---------------------------------------------------------------!
               if (omega_ == omega_prime_) then
                  !------------------------------------------------------------!
                  ! Eq. 4 in "Coupling Matrix" section
                  !------------------------------------------------------------!
                  centrifugal_matrix_(channel_index_1_, channel_index_2_)      &
                     = calculate_diagonal_centrifugal_element(                 &
                     total_angular_momentum_, j_, omega_)
               else
                  !------------------------------------------------------------!
                  ! Eq. 5 in "Coupling Matrix" section
                  !------------------------------------------------------------!
                  centrifugal_matrix_(channel_index_1_, channel_index_2_)      &
                     = calculate_offdiagonal_centrifugal_element(              &
                     total_angular_momentum_, j_, omega_, omega_prime_,        &
                     delta_1_, delta_2_)
               endif
               !---------------------------------------------------------------!
            enddo
         enddo
         !---------------------------------------------------------------------!
         call fill_symmetric_matrix(centrifugal_matrix_, 'u')
         !---------------------------------------------------------------------!
      end subroutine calculate_centrifugal_matrix