save_s_matrix_block_info Subroutine

public subroutine save_s_matrix_block_info(total_angular_momentum, parity_exponent, number_of_open_channels, channel_indices, channel_l_values, block_wavevectors, s_matrix_real, s_matrix_imag)

save information about current block -- total angular momentum, parity exponent, number of open channels in the current block -- array of indices pointing to the basis arrays, array holding \(l\) values, wavenumbers -- real part of the S-matrix -- imaginary part of the S-matrix

Arguments

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

total angular momentum of the current block

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

parity exponent of the current block

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

number of open channels in the block

integer(kind=int32), intent(in) :: channel_indices(number_of_open_channels)

holds the indices pointing to the basis arrays

integer(kind=int32), intent(in) :: channel_l_values(number_of_open_channels)

holds all values of \(l\)

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

holds all values of wavevectors in the block

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

real part of the S-matrix

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

imaginary part of the S-matrix


Called by

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

Contents


Variables

Type Visibility Attributes Name Initial
integer(kind=int32), private :: channel_index_
integer(kind=int32), private :: channel_index_2_

Source Code

      subroutine save_s_matrix_block_info(total_angular_momentum,              &
         parity_exponent, number_of_open_channels, channel_indices,            &
         channel_l_values, block_wavevectors, s_matrix_real, s_matrix_imag)
         !! save information about current block
         !! -- total angular momentum, parity exponent, number of open channels
         !!    in the current block
         !! -- array of indices pointing to the basis arrays, array holding
         !!    \\(l\\) values, wavenumbers
         !! -- real part of the S-matrix
         !! -- imaginary part of the S-matrix
         !---------------------------------------------------------------------!
         integer(int32), intent(in) :: total_angular_momentum
            !! total angular momentum of the current block
         integer(int32), intent(in) :: parity_exponent
            !! parity exponent of the current block
         integer(int32), intent(in) :: number_of_open_channels
            !! number of open channels in the block
         integer(int32), intent(in) :: channel_indices(number_of_open_channels)
            !! holds the indices pointing to the basis arrays
         integer(int32), intent(in) :: channel_l_values(number_of_open_channels)
            !! holds all values of \\(l\\)
         real(dp), intent(in) :: block_wavevectors(number_of_open_channels)
            !! holds all values of wavevectors in the block
         real(dp), intent(in) :: s_matrix_real(number_of_open_channels,        &
            number_of_open_channels)
            !! real part of the S-matrix
         real(dp), intent(in) :: s_matrix_imag(number_of_open_channels,        &
            number_of_open_channels)
            !! imaginary part of the S-matrix
         !---------------------------------------------------------------------!
         integer(int32) :: channel_index_, channel_index_2_
         !---------------------------------------------------------------------!
         write(s_matrix_unit) total_angular_momentum, parity_exponent,         &
            number_of_open_channels
         write(s_matrix_unit) (channel_indices(channel_index_),                &
            channel_l_values(channel_index_),block_wavevectors(channel_index_),&
            channel_index_ = 1, number_of_open_channels)
         write(s_matrix_unit)((s_matrix_real(channel_index_, channel_index_2_),&
            channel_index_2_ = 1, channel_index_), channel_index_=1,           &
            number_of_open_channels)
         write(s_matrix_unit) ((s_matrix_imag(channel_index_,channel_index_2_),&
            channel_index_2_ = 1, channel_index_), channel_index_=1,           &
            number_of_open_channels)
         !---------------------------------------------------------------------!
      end subroutine save_s_matrix_block_info