calculate_sf_matrix_from_bf_matrix Subroutine

public subroutine calculate_sf_matrix_from_bf_matrix(number_of_channels, total_angular_momentum_, channel_indices, channels_omega_values, channel_l_values, bf_matrix, sf_matrix)

takes as an input matrix in the body-fixed frame and transforms it to the spec-fixed frame; iterates over all matrix elements and calls calculate_single_SF_element

Arguments

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

size of the basis

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

total angular momentum

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) :: channel_l_values(number_of_channels)

holds all values of l

real(kind=dp), intent(in) :: bf_matrix(number_of_channels,number_of_channels)

matrix in the BF frame

real(kind=dp), intent(inout) :: sf_matrix(number_of_channels,number_of_channels)

(output) matrix in the SF frame


Calls

proc~~calculate_sf_matrix_from_bf_matrix~~CallsGraph proc~calculate_sf_matrix_from_bf_matrix calculate_sf_matrix_from_bf_matrix proc~calculate_single_sf_element calculate_single_SF_element proc~calculate_sf_matrix_from_bf_matrix->proc~calculate_single_sf_element proc~p_coeff p_coeff proc~calculate_single_sf_element->proc~p_coeff fwig3jj fwig3jj proc~p_coeff->fwig3jj

Called by

proc~~calculate_sf_matrix_from_bf_matrix~~CalledByGraph proc~calculate_sf_matrix_from_bf_matrix calculate_sf_matrix_from_bf_matrix program~scattering SCATTERING program~scattering->proc~calculate_sf_matrix_from_bf_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 :: j1_
integer(kind=int32), private :: j1p_
integer(kind=int32), private :: l_
integer(kind=int32), private :: lp_
integer(kind=int32), private :: omega_
integer(kind=int32), private :: omegap_
real(kind=dp), private :: single_sf_element
integer(kind=int32), private :: v1_
integer(kind=int32), private :: v1p_

Source Code

      subroutine calculate_sf_matrix_from_bf_matrix(number_of_channels,        &
         total_angular_momentum_, channel_indices, channels_omega_values,      &
         channel_l_values, bf_matrix, sf_matrix)
         !! takes as an input matrix in the body-fixed frame and transforms it 
         !! to the spec-fixed frame; iterates over all matrix elements
         !! and calls calculate_single_SF_element
         !---------------------------------------------------------------------!
         integer(int32), intent(in) :: number_of_channels
            !! size of the basis
         integer(int32), intent(in) :: total_angular_momentum_
            !! total angular momentum
         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) :: channel_l_values(number_of_channels)
            !! holds all values of l
         real(dp), intent(in) :: bf_matrix(number_of_channels, number_of_channels)
            !! matrix in the BF frame
         real(dp), intent(inout) :: sf_matrix(number_of_channels, number_of_channels)
            !! (output) matrix in the SF frame
         !---------------------------------------------------------------------!
         integer(int32) :: l_, lp_, omega_, omegap_, v1_, j1_, v1p_, j1p_,     &
            channel_index_1_, channel_index_2_
         real(dp) :: single_sf_element
         !---------------------------------------------------------------------!
         do channel_index_1_ = 1, number_of_channels
            v1_ = vib_levels(channel_indices(channel_index_1_))
            j1_ = rot_levels(channel_indices(channel_index_1_))
            l_  = channel_l_values(channel_index_1_)
            do channel_index_2_ = 1, number_of_channels
               v1p_ = vib_levels(channel_indices(channel_index_2_))
               j1p_ = rot_levels(channel_indices(channel_index_2_))
               lp_  = channel_l_values(channel_index_2_)
               call calculate_single_SF_element(number_of_channels,            &
                  total_angular_momentum_, v1_, j1_, v1p_, j1p_, l_, lp_,      &
                  channel_indices, channels_omega_values, bf_matrix,    &
                  single_sf_element)
               sf_matrix(channel_index_1_, channel_index_2_) = single_sf_element
            enddo
         enddo
         !---------------------------------------------------------------------!         
      end subroutine calculate_sf_matrix_from_bf_matrix