tem_identify_prevailDirections Subroutine

public subroutine tem_identify_prevailDirections(me, cxDir)

This subroutine fills the array of prevailing directions according to the given array of directions.

Arguments

Type IntentOptional Attributes Name
real(kind=rk), intent(out), allocatable :: me(:,:)

growing array of prevailing directions

integer, intent(in) :: cxDir(:,:)

array of given offsets


Contents


Source Code

  subroutine tem_identify_prevailDirections( me, cxDir )
    ! -------------------------------------------------------------------- !
    !> growing array of prevailing directions
    real( kind=rk ), allocatable, intent(out) :: me(:,:)
    !> array of given offsets
    integer, intent(in)  :: cxDir(:,:)
    ! -------------------------------------------------------------------- !
    real(kind=rk),allocatable :: tmpPrevailDir(:,:)
    real(kind=rk) :: length_rk
    integer :: length
    ! counters
    integer :: iOffset
    integer :: nEntries ! size of array
    integer :: counter
    ! -------------------------------------------------------------------- !

    nEntries = size(cxDir, 2)
    allocate( tmpPrevailDir(3, nEntries) )
    counter = 0

    do iOffset=1,nEntries

      length =   cxDir(1,iOffset)**2 &
        &      + cxDir(2,iOffset)**2 &
        &      + cxDir(3,iOffset)**2

      if (length /= 0) then
        counter = counter + 1
        length_rk = sqrt(real(length, kind=rk))
        tmpPrevailDir(1:3, counter) = real(cxDir(1:3, iOffset), kind=rk) &
          &                          / length_rk
      end if

    end do

    allocate( me( 3, counter ))

    do iOffset = 1, counter
      me(:, iOffset) = tmpPrevailDir( :, iOffset)
    end do

  end subroutine tem_identify_prevailDirections