determine_largest_cross_sections Subroutine

public subroutine determine_largest_cross_sections(open_basis_levels_, cross_sections_, max_elastic_cross_section_, max_inelastic_cross_section_, max_elastic_index_, max_inelastic_index_1_, max_inelastic_index_2_)

Determine the largest partial elastic and inelastic cross-sections in a given set of cross-sections.

Arguments

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

array holding indices to open basis levels

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

array holding partial cross-sections

real(kind=dp), intent(out) :: max_elastic_cross_section_

largest elastic cross-section

real(kind=dp), intent(out) :: max_inelastic_cross_section_

largest inelastic cross-section

integer(kind=int32), intent(out) :: max_elastic_index_

index pointing to the largest elastic cross-section

integer(kind=int32), intent(out) :: max_inelastic_index_1_

first index pointing to the largest inelastic cross-section

integer(kind=int32), intent(out) :: max_inelastic_index_2_

second index pointing to the largest inelastic cross-section


Called by

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

Contents


Variables

Type Visibility Attributes Name Initial
integer(kind=int32), private :: cross_section_index_
integer(kind=int32), private :: index_1_
integer(kind=int32), private :: index_2_
integer(kind=int32), private :: number_of_open_basis_levels_

Source Code

      subroutine determine_largest_cross_sections(open_basis_levels_,          &
         cross_sections_, max_elastic_cross_section_,                          &
         max_inelastic_cross_section_, max_elastic_index_,                     &
         max_inelastic_index_1_, max_inelastic_index_2_)
         !! Determine the largest partial elastic and inelastic cross-sections
         !! in a given set of cross-sections.
         !---------------------------------------------------------------------!
         integer(int32), intent(in) :: open_basis_levels_(:)
            !! array holding indices to open basis levels
         real(dp), intent(in) :: cross_sections_(:)
            !! array holding partial cross-sections
         real(dp), intent(out) :: max_elastic_cross_section_
            !! largest elastic cross-section
         real(dp), intent(out) :: max_inelastic_cross_section_
            !! largest inelastic cross-section
         integer(int32), intent(out) :: max_elastic_index_
            !! index pointing to the largest elastic cross-section
         integer(int32), intent(out) :: max_inelastic_index_1_
            !! first index pointing to the largest inelastic cross-section
         integer(int32), intent(out) :: max_inelastic_index_2_
            !! second index pointing to the largest inelastic cross-section
         !---------------------------------------------------------------------!
         integer(int32) :: number_of_open_basis_levels_, index_1_, index_2_,   &
            cross_section_index_
         !---------------------------------------------------------------------!
         number_of_open_basis_levels_ = size(open_basis_levels_)
         !---------------------------------------------------------------------!
         ! Initialize output values
         !---------------------------------------------------------------------!
         max_elastic_cross_section_ = 0.0_dp
         max_inelastic_cross_section_ = 0.0_dp
         max_elastic_index_ = 0
         max_inelastic_index_1_ = 0
         max_inelastic_index_2_ = 0
         !---------------------------------------------------------------------!
         ! Iterate over all combinations of open basis levels
         !---------------------------------------------------------------------!
         do index_1_ = 1, number_of_open_basis_levels_
            do index_2_ = 1, number_of_open_basis_levels_
               cross_section_index_ = (index_1_-1)                             &
                  * number_of_open_basis_levels_ + index_2_
               if (open_basis_levels_(index_2_) == open_basis_levels_(index_1_)) then
                  !------------------------------------------------------------!
                  ! Elastic cross-section
                  !------------------------------------------------------------!
                  if (cross_sections_(cross_section_index_)                    &
                     > max_elastic_cross_section_) then
                     max_elastic_cross_section_                                &
                        = cross_sections_(cross_section_index_)
                     max_elastic_index_ = index_1_
                  endif
               else
                  !------------------------------------------------------------!
                  ! Inelastic cross-section
                  !------------------------------------------------------------!
                  if (cross_sections_(cross_section_index_)                    &
                     > max_inelastic_cross_section_) then
                     max_inelastic_cross_section_                              &
                        = cross_sections_(cross_section_index_)
                     max_inelastic_index_1_ = index_1_
                     max_inelastic_index_2_ = index_2_
                  endif
               endif
            enddo
         enddo
         !---------------------------------------------------------------------!
      end subroutine determine_largest_cross_sections