unitarity_check Subroutine

public subroutine unitarity_check(number_of_open_channels, s_matrix_real, s_matrix_imag, is_unitary)

checks the unitarity of the S-matrix (Eq. (13) 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) :: 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

logical, intent(inout) :: is_unitary

(output) if .true. unitarity is fulfilled, .false. otherwise


Calls

proc~~unitarity_check~~CallsGraph proc~unitarity_check unitarity_check proc~write_message write_message proc~unitarity_check->proc~write_message proc~handle_unitarity_output_message handle_unitarity_output_message proc~unitarity_check->proc~handle_unitarity_output_message proc~check_unitarity_for_each_channel check_unitarity_for_each_channel proc~unitarity_check->proc~check_unitarity_for_each_channel proc~calculate_sum_of_squares_for_each_channel calculate_sum_of_squares_for_each_channel proc~unitarity_check->proc~calculate_sum_of_squares_for_each_channel proc~handle_unitarity_output_message->proc~write_message proc~print_sum_of_squares print_sum_of_squares proc~handle_unitarity_output_message->proc~print_sum_of_squares proc~write_warning write_warning proc~handle_unitarity_output_message->proc~write_warning proc~print_sum_of_squares->proc~write_message proc~integer_to_character integer_to_character proc~print_sum_of_squares->proc~integer_to_character proc~float_to_character float_to_character proc~print_sum_of_squares->proc~float_to_character proc~write_warning->proc~write_message

Called by

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

Contents

Source Code


Variables

Type Visibility Attributes Name Initial
integer(kind=int32), private :: channel_index
real(kind=dp), private :: sum_of_squares(number_of_open_channels)

Source Code

      subroutine unitarity_check(number_of_open_channels, s_matrix_real,       &
         s_matrix_imag,is_unitary)
         !! checks the unitarity of the S-matrix
         !! (Eq. (13) in "Solution of coupled equations")
         !---------------------------------------------------------------------!
         integer(int32), intent(in) :: number_of_open_channels
            !! number of open channels
         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
         logical, intent(inout) :: is_unitary
            !! (output) if .true. unitarity is fulfilled, .false. otherwise
         !---------------------------------------------------------------------!
         integer(int32) :: channel_index
         real(dp) :: sum_of_squares(number_of_open_channels)
         !---------------------------------------------------------------------!
         is_unitary  = .true.
         !---------------------------------------------------------------------!
         call write_message("-- Check of the unitarity of the S-matrix...")
         !---------------------------------------------------------------------!
         ! Calculating sum of squares for each channel
         !---------------------------------------------------------------------!
         call calculate_sum_of_squares_for_each_channel(s_matrix_real,         &
            s_matrix_imag, sum_of_squares)
         !---------------------------------------------------------------------!
         ! Checking unitarity for each channel
         !---------------------------------------------------------------------!
         is_unitary = check_unitarity_for_each_channel(sum_of_squares)
         !---------------------------------------------------------------------!
         ! Handling the output message based on unitarity check
         !---------------------------------------------------------------------!
         call handle_unitarity_output_message(is_unitary, sum_of_squares)
         !---------------------------------------------------------------------!
      end subroutine unitarity_check