compute_individual_cross_section Function

private function compute_individual_cross_section(initial_state_, final_state_, open_basis_levels_, basis_wavevectors_, s_matrix_real_, s_matrix_imag_, channel_indices_, channel_l_values_, total_angular_momentum_) result(cross_section_)

Calculates cross-section for a given initial and final state.

Arguments

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

index pointing to the intial state in basis arrays

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

index pointing to the final state in basis arrays

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

holds indices to the basis arrays which correspond to open channels

real(kind=dp), intent(in) :: basis_wavevectors_(:)

holds wavevectors k_{i}

real(kind=dp), intent(in) :: s_matrix_real_(:,:)

real and imaginary parts of the S-matrix

real(kind=dp), intent(in) :: s_matrix_imag_(:,:)

real and imaginary parts of the S-matrix

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

holds the indices pointing to the basis arrays

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

holds all values of l

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

total angular momentum

Return Value real(kind=dp)

(output) cross-section


Calls

proc~~compute_individual_cross_section~~CallsGraph proc~compute_individual_cross_section compute_individual_cross_section proc~get_block_indices get_block_indices proc~compute_individual_cross_section->proc~get_block_indices proc~sum_cross_section_contributions sum_cross_section_contributions proc~compute_individual_cross_section->proc~sum_cross_section_contributions proc~compute_real_component compute_real_component proc~sum_cross_section_contributions->proc~compute_real_component proc~compute_imag_component compute_imag_component proc~sum_cross_section_contributions->proc~compute_imag_component

Called by

proc~~compute_individual_cross_section~~CalledByGraph proc~compute_individual_cross_section compute_individual_cross_section proc~calculate_state_to_state_cross_section calculate_state_to_state_cross_section proc~calculate_state_to_state_cross_section->proc~compute_individual_cross_section program~scattering SCATTERING program~scattering->proc~calculate_state_to_state_cross_section

Contents


Variables

Type Visibility Attributes Name Initial
integer(kind=int32), private, allocatable :: final_block_indices_(:)
integer(kind=int32), private :: final_index_
integer(kind=int32), private, allocatable :: init_block_indices_(:)
integer(kind=int32), private :: initial_index_
integer(kind=int32), private :: j_final_
integer(kind=int32), private :: j_initial_
real(kind=dp), private :: sum_contributions_
integer(kind=int32), private :: v_final_
integer(kind=int32), private :: v_initial_
real(kind=dp), private :: wavevector_initial_

Source Code

      function compute_individual_cross_section(initial_state_, final_state_,  &
         open_basis_levels_, basis_wavevectors_, s_matrix_real_,          &
         s_matrix_imag_, channel_indices_, channel_l_values_,                  &
         total_angular_momentum_) result(cross_section_)
         !! Calculates cross-section for a given initial and final state.
         !---------------------------------------------------------------------!
         integer(int32), intent(in) :: initial_state_
            !! index pointing to the intial state in basis arrays
         integer(int32), intent(in) :: final_state_
            !! index pointing to the final state in basis arrays
         integer(int32), intent(in) :: open_basis_levels_(:)
            !! holds indices to the basis arrays which correspond to open channels
         real(dp), intent(in) :: basis_wavevectors_(:)
            !! holds wavevectors k_{i}
         real(dp), intent(in) :: s_matrix_real_(:,:), s_matrix_imag_(:,:)
            !! real and imaginary parts of the S-matrix
         integer(int32), intent(in) :: channel_indices_(:)
            !! holds the indices pointing to the basis arrays
         integer(int32), intent(in) :: channel_l_values_(:)
            !! holds all values of l
         integer(int32), intent(in) :: total_angular_momentum_
            !! total angular momentum
         real(dp) :: cross_section_
            !! (output) cross-section
         !---------------------------------------------------------------------!
         integer(int32) :: initial_index_, final_index_, v_initial_,           &
            j_initial_, v_final_, j_final_
         real(dp) :: wavevector_initial_, sum_contributions_
         integer(int32), allocatable :: init_block_indices_(:),                &
            final_block_indices_(:)
         !---------------------------------------------------------------------!
         initial_index_ = open_basis_levels_(initial_state_)
         v_initial_ = vib_levels(initial_index_)
         j_initial_ = rot_levels(initial_index_)
         wavevector_initial_ = basis_wavevectors_(initial_state_)
         !---------------------------------------------------------------------!
         final_index_ = open_basis_levels_(final_state_)
         v_final_ = vib_levels(final_index_)
         j_final_ = rot_levels(final_index_)
         !---------------------------------------------------------------------!
         init_block_indices_ = get_block_indices(v_initial_, j_initial_,       &
            channel_indices_)
         final_block_indices_ = get_block_indices(v_final_, j_final_,          &
            channel_indices_)
         !---------------------------------------------------------------------!
         sum_contributions_ = sum_cross_section_contributions(                 &
            init_block_indices_, final_block_indices_, s_matrix_real_,         &
            s_matrix_imag_, channel_l_values_)
         !---------------------------------------------------------------------!
         cross_section_ = sum_contributions_                                   &
            * real(2 * total_angular_momentum_ + 1, dp)                        &
            / (real(2 * j_initial_ + 1, dp) * wavevector_initial_**2.0_dp) * pi
         !---------------------------------------------------------------------!
      end function compute_individual_cross_section