print_largest_partial_cross_sections Subroutine

public subroutine print_largest_partial_cross_sections(total_angular_momentum_, largest_elastic_xs_, largest_inelastic_xs_, elastic_index_, inelastic_index_1_, inelastic_index_2_, open_basis_levels_)

Print the largest partial elastic and inelastic state-to-state cross-sections in a given block.

Arguments

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

total angular momentum

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

the largest partial elastic state-to-state XS in the block

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

the largest partial inelastic state-to-state XS in the block

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

index pointing indirectly to quantum numbers associated with the largest partial elastic state-to-state XS in the block

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

indices pointing indirectly to quantum numbers associated with the largest partial inelastic state-to-state XS in the block

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

indices pointing indirectly to quantum numbers associated with the largest partial inelastic state-to-state XS in the block

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

array holding indices to open basis levels


Calls

proc~~print_largest_partial_cross_sections~~CallsGraph proc~print_largest_partial_cross_sections print_largest_partial_cross_sections proc~print_basic_cross_section_info print_basic_cross_section_info proc~print_largest_partial_cross_sections->proc~print_basic_cross_section_info proc~print_detailed_cross_section_info print_detailed_cross_section_info proc~print_largest_partial_cross_sections->proc~print_detailed_cross_section_info proc~integer_to_character integer_to_character proc~print_basic_cross_section_info->proc~integer_to_character proc~float_to_character float_to_character proc~print_basic_cross_section_info->proc~float_to_character proc~write_message write_message proc~print_basic_cross_section_info->proc~write_message proc~print_detailed_cross_section_info->proc~integer_to_character proc~print_detailed_cross_section_info->proc~write_message

Called by

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

Contents


Source Code

      subroutine print_largest_partial_cross_sections(total_angular_momentum_, &
         largest_elastic_xs_, largest_inelastic_xs_, elastic_index_,           &
         inelastic_index_1_, inelastic_index_2_, open_basis_levels_)
         !! Print the largest partial elastic and inelastic state-to-state
         !! cross-sections in a given block.
         !---------------------------------------------------------------------!
         integer(int32), intent(in) :: total_angular_momentum_
            !! total angular momentum
         real(dp), intent(in) :: largest_elastic_xs_
            !! the largest partial elastic state-to-state XS in the block
         real(dp), intent(in) :: largest_inelastic_xs_
            !! the largest partial inelastic state-to-state XS in the block
         integer(int32), intent(in) :: elastic_index_
            !! index pointing indirectly to quantum numbers associated with
            !! the largest partial elastic state-to-state XS in the block
         integer(int32), intent(in) :: inelastic_index_1_, inelastic_index_2_
            !! indices pointing indirectly to quantum numbers associated with
            !! the largest partial inelastic state-to-state XS in the block
         integer(int32), intent(in) :: open_basis_levels_(:)
            !! array holding indices to open basis levels 
         !---------------------------------------------------------------------!
         if (print_level <= 2) then
            call print_basic_cross_section_info(total_angular_momentum_,       &
               largest_elastic_xs_, "elastic")
            call print_basic_cross_section_info(total_angular_momentum_,       &
               largest_inelastic_xs_, "inelastic")
         else if (print_level >= 3) then
            call print_detailed_cross_section_info(total_angular_momentum_,    &
               largest_elastic_xs_, elastic_index_, elastic_index_,            &
               open_basis_levels_, "elastic")
            call print_detailed_cross_section_info(total_angular_momentum_,    &
               largest_inelastic_xs_, inelastic_index_1_, inelastic_index_2_,  &
               open_basis_levels_, "inelastic")
         endif
         !---------------------------------------------------------------------!
      end subroutine print_largest_partial_cross_sections