ESMF_GridConvertIndex Subroutine

public subroutine ESMF_GridConvertIndex(grid, gridindex, distgridindex, rc)

Arguments

Type IntentOptional Attributes Name
type(ESMF_Grid), intent(in) :: grid
integer, intent(in) :: gridindex(:)
integer, intent(out) :: distgridindex(:)
integer, intent(out), optional :: rc

Source Code

      subroutine ESMF_GridConvertIndex(grid,gridindex, distgridindex, rc)
!
! !ARGUMENTS:
         type(ESMF_Grid),       intent(in)            :: grid
         integer              , intent(in)            :: gridindex(:)
         integer              , intent(out)            :: distgridindex(:)
         integer              , intent(out), optional :: rc
!
! !DESCRIPTION:
!
!  Convert a multi-dimensional index of the arbitrarily distributed grid into the
!  index of the 1D DistGrid.  The associated DistGrid for an arbitrarily distributed
!  grid is 1D plus any undistributed dimension.  The function
!  calculates the index of the DistGrid for a given index from the original Grid.
!
! The arguments are:
! \begin{description}
! \item[{grid}]
!     The grid to get the information from to create the Array.
! \item[{[gridindex]}]
!     The Grid index to be converted.
! \item[{[distgridindex]}]
!     The DistGrid index to be returned.
! \item[{[rc]}]
!      Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!\end{description}
!
!EOPI

    integer ::  localrc
    integer ::  DimCount, distDimCount, undistDimCount
    integer, pointer ::  minIndex(:)
    integer, pointer ::  maxIndex(:)
    integer, pointer ::  distgridToGridMap(:)
    integer          :: i,j,k
    integer ::  index1D    ! the return value
    type(ESMF_InterArray)   :: gridIndexArg
    type(ESMF_GridDecompType) :: decompType
    type(ESMF_DistGrid) :: distGrid
    integer, allocatable :: undistdim(:)
    logical  :: found
    integer :: distGridDimCount, arbDim


    ! Initialize return code; assume failure until success is certain
    localrc = ESMF_RC_NOT_IMPL
    if (present(rc)) rc = ESMF_RC_NOT_IMPL

    ! Check init status of arguments
    ESMF_INIT_CHECK_DEEP_SHORT(ESMF_GridGetInit, grid, rc)

    ! Get Grid decomposition type
    call ESMF_GridGetDecompType(grid, decompType, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
      ESMF_CONTEXT, rcToReturn=rc)) return
        
    ! Check if the grid is arbitrary
    if (decompType /= ESMF_GRID_ARBITRARY) then
       call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, &
                 msg="- ESMF_GridConvertIndex only works for arbritrarily distributed grid", &
                 ESMF_CONTEXT, rcToReturn=rc)
       return
    endif

    ! Get info from Grid
    call ESMF_GridGet(grid, distgrid= distGrid, DimCount=DimCount, &
                      rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
      ESMF_CONTEXT, rcToReturn=rc)) return

    ! allocate minIndex and maxIndex
    allocate(minIndex(DimCount), maxIndex(DimCount), stat=localrc)
    if (ESMF_LogFoundAllocError(localrc, msg="Allocating minIndex and maxIndex", &
                                     ESMF_CONTEXT, rcToReturn=rc)) return

    ! Get minIndex and maxIndex from the grid
    call ESMF_GridGetIndex(grid, minIndex= minIndex, maxIndex=maxIndex, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
      ESMF_CONTEXT, rcToReturn=rc)) return

    ! find out how many dimensions are arbitrarily distributed
     call ESMF_DistGridGet(distGrid, dimcount = distGridDimCount, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
      ESMF_CONTEXT, rcToReturn=rc)) return
    if (distGridDimCount > dimCount) then
        call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, &
                   msg="- distgrid dimension has to be less than or equal to dimCount", &
                     ESMF_CONTEXT, rcToReturn=rc)
        return
     endif

    ! set distDimCount - number of dimensions arbitrarily distributed
    !     undistDimCount - number of dimensions not arbitrarily distributed
    if (distGridDimCount == 1) then
       ! all dimensions are arbitrarily distributed
       distDimCount = dimCount
       undistDimCount = 0
    else
       undistDimCount = distGridDimCount - 1
       distDimCount = dimCount - undistDimCount
    endif

    ! Check index dimension
    if (size(gridindex) /= dimCount) then
       call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, &
            msg="- gridindex dimension is different from the grid DimCount", &
            ESMF_CONTEXT, rcToReturn=rc)
       return
    endif

    ! Check index out of bound
    do i=1,dimCount
       if (gridindex(i) .lt. minIndex(i) .and. gridindex(i) > maxIndex(i)) then
           call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, &
               msg="- gridindex is out of bound", &
                      ESMF_CONTEXT, rcToReturn=rc)
           return
       endif
    enddo

    ! clean up memory allocation
    deallocate(minIndex)
    deallocate(maxIndex)

    ! Call the C function to get the index of the 1D distgrid
    !! index
    gridIndexArg = ESMF_InterArrayCreate(gridindex, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
      ESMF_CONTEXT, rcToReturn=rc)) return

    call c_ESMC_gridconvertindex(grid%this, gridIndexArg, index1D, localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
      ESMF_CONTEXT, rcToReturn=rc)) return

    if (undistDimCount /= 0) then
      allocate(distgridToGridMap(dimCount), stat=localrc)
      call ESMF_GridGet(grid, arbDim=arbDim, &
        distgridToGridMap=distgridToGridMap, rc=localrc)
      if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
        ESMF_CONTEXT, rcToReturn=rc)) return
      k=1
      allocate(undistdim(undistDimCount))
      do i=1,dimCount
        found = .false.
        do j=1,distDimCount
          if (i == distgridToGridMap(j)) found=.true.
        enddo
        if (.not. found) then
          undistdim(k)=i
          k=k+1
        endif
      enddo

      k=1
      do i=1,distGridDimCount
        if (i == arbDim) then
           distgridindex(i)=index1D
        else
           distgridindex(i)=gridindex(undistdim(k))
           k=k+1
        endif
      enddo
      deallocate(undistdim)
      deallocate(distgridToGridMap)
    else
      distgridindex(1)=index1D
    endif

    ! clean up memory allocation
    call ESMF_InterArrayDestroy(GridIndexArg, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
      ESMF_CONTEXT, rcToReturn=rc)) return

    ! Return successfully
    if (present(rc)) rc = ESMF_SUCCESS
    return

end subroutine ESMF_GridConvertIndex