calculate_s_matrix Subroutine

public subroutine calculate_s_matrix(number_of_open_channels, k_matrix, s_matrix_real, s_matrix_imag)

calculates S-matrix from open-open portion of the K-matrix using Eq. (12) in "Solution of coupled equations"

Arguments

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

number of open channels

real(kind=dp), intent(in) :: k_matrix(number_of_open_channels,number_of_open_channels)

K-matrix

real(kind=dp), intent(inout) :: s_matrix_real(number_of_open_channels,number_of_open_channels)

(output) real part of the S-matrix

real(kind=dp), intent(inout) :: s_matrix_imag(number_of_open_channels,number_of_open_channels)

(output) imaginary part of the S-matrix


Calls

proc~~calculate_s_matrix~~CallsGraph proc~calculate_s_matrix calculate_s_matrix dgemm dgemm proc~calculate_s_matrix->dgemm interface~invert_symmetric_matrix invert_symmetric_matrix proc~calculate_s_matrix->interface~invert_symmetric_matrix interface~fill_symmetric_matrix fill_symmetric_matrix proc~calculate_s_matrix->interface~fill_symmetric_matrix

Called by

proc~~calculate_s_matrix~~CalledByGraph proc~calculate_s_matrix calculate_s_matrix program~scattering SCATTERING program~scattering->proc~calculate_s_matrix

Contents

Source Code


Variables

Type Visibility Attributes Name Initial
integer(kind=int32), private :: open_channel_index_1_
integer(kind=int32), private :: open_channel_index_2_
real(kind=dp), private :: s_tmp_matrix(number_of_open_channels,number_of_open_channels)

Source Code

      subroutine calculate_s_matrix(number_of_open_channels, k_matrix,         &
         s_matrix_real, s_matrix_imag)
         !! calculates S-matrix from open-open portion of the K-matrix using
         !! Eq. (12) in "Solution of coupled equations" 
         !---------------------------------------------------------------------!
         integer(int32), intent(in) :: number_of_open_channels
            !! number of open channels
         real(dp), intent(in) :: k_matrix(number_of_open_channels,number_of_open_channels)
            !! K-matrix
         real(dp), intent(inout) :: s_matrix_real(number_of_open_channels,number_of_open_channels)
            !! (output) real part of the S-matrix
         real(dp), intent(inout) :: s_matrix_imag(number_of_open_channels,number_of_open_channels)
            !! (output) imaginary part of the S-matrix
         !---------------------------------------------------------------------!
         integer(int32) :: open_channel_index_1_, open_channel_index_2_
         real(dp) :: s_tmp_matrix(number_of_open_channels,number_of_open_channels)
         !---------------------------------------------------------------------!
         s_matrix_real = 0
         s_matrix_imag = 0
         !---------------------------------------------------------------------!
         !---------------------------------------------------------------------!
         call DGEMM('N','N',number_of_open_channels,number_of_open_channels,   &
            number_of_open_channels,0.5_dp,k_matrix,number_of_open_channels,   &
            k_matrix,number_of_open_channels,0.0_dp,s_tmp_matrix,number_of_open_channels)
         !---------------------------------------------------------------------!
         do open_channel_index_1_ = 1, number_of_open_channels
            s_tmp_matrix(open_channel_index_1_, open_channel_index_1_) =       &
               s_tmp_matrix(open_channel_index_1_, open_channel_index_1_)      &
               + 0.5_dp
         enddo
         !---------------------------------------------------------------------!
         call invert_symmetric_matrix(s_tmp_matrix)
         call fill_symmetric_matrix(s_tmp_matrix, 'u')
         !---------------------------------------------------------------------!
         call DGEMM('N','N',number_of_open_channels,number_of_open_channels,   &
            number_of_open_channels,-1.0_dp,s_tmp_matrix,                      &
            number_of_open_channels,k_matrix,number_of_open_channels,0.0_dp,   &
            s_matrix_imag,number_of_open_channels)
         !---------------------------------------------------------------------!
         do open_channel_index_1_ = 1, number_of_open_channels
            do open_channel_index_2_ = 1, number_of_open_channels
               s_matrix_real(open_channel_index_1_, open_channel_index_2_) =   &
                  s_tmp_matrix(open_channel_index_1_, open_channel_index_2_)
            enddo
            s_matrix_real(open_channel_index_1_, open_channel_index_1_) =      &
               s_matrix_real(open_channel_index_1_, open_channel_index_1_) - 1.0_dp
         enddo
         !---------------------------------------------------------------------!
      end subroutine calculate_s_matrix