set_space_fixed_channels Subroutine

public subroutine set_space_fixed_channels(total_angular_momentum_, parity_exponent_, channel_l_values)

Prepares the channel_l_values array which holds values of orbital angular momentum, \(l\), a space-fixed-frame quantum number.

Arguments

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

total angular momentum

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

parity exponent of the block (0 if p = +1, 1 if p = -1)

integer(kind=int32), intent(inout) :: channel_l_values(:)

holds all values of l


Calls

proc~~set_space_fixed_channels~~CallsGraph proc~set_space_fixed_channels set_space_fixed_channels proc~write_error write_error proc~set_space_fixed_channels->proc~write_error proc~write_message write_message proc~write_error->proc~write_message

Called by

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

Contents


Variables

Type Visibility Attributes Name Initial
integer, private :: channel_index_
integer, private :: l_
integer, private :: l_max_
integer, private :: l_min_
integer, private :: level_index_

Source Code

      subroutine set_space_fixed_channels(total_angular_momentum_,             &
         parity_exponent_, channel_l_values)
         !! Prepares the channel_l_values array which holds values of orbital
         !! angular momentum, \\(l\\), a space-fixed-frame quantum number.
         !---------------------------------------------------------------------!
         integer(int32), intent(in) :: total_angular_momentum_
            !! total angular momentum
         integer(int32), intent(in) :: parity_exponent_
            !! parity exponent of the block (0 if p = +1, 1 if p = -1)
         integer(int32), intent(inout) :: channel_l_values(:)
            !! holds all values of l
         !---------------------------------------------------------------------!
         integer :: level_index_, l_min_, l_max_, l_, channel_index_
         !---------------------------------------------------------------------!
         channel_index_ = 0

         do level_index_ = 1, number_of_basis_levels
            l_min_ = abs(total_angular_momentum_ - rot_levels(level_index_))
            l_max_ = total_angular_momentum_ + rot_levels(level_index_)

            do l_ = l_min_, l_max_
               if (mod(l_ + rot_levels(level_index_), 2)==parity_exponent_) then
                  channel_index_ = channel_index_ + 1
                  if (channel_index_ > size(channel_l_values)) then
                     call write_error("channel_index_ out of bounds of " //    &
                        "channel_l_values in set_space_fixed_channels.")
                  end if
                  channel_l_values(channel_index_) = l_
               endif
            enddo
         enddo
         !---------------------------------------------------------------------!
      end subroutine set_space_fixed_channels