ESMF_GridCreateFrmGridspec Function

private function ESMF_GridCreateFrmGridspec(grid_filename, regDecomp, indexflag, keywordEnforcer, decompflag, addMask, varname, coordNames, isSphere, polekindflag, addCornerStagger, coordTypeKind, rc)

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: grid_filename
integer, intent(in) :: regDecomp(:)
type(ESMF_Index_Flag), intent(in) :: indexflag
type(ESMF_KeywordEnforcer), optional :: keywordEnforcer
type(ESMF_Decomp_Flag), intent(in), optional :: decompflag(:)
logical, intent(in), optional :: addMask
character(len=*), intent(in), optional :: varname
character(len=*), intent(in), optional :: coordNames(:)
logical, intent(in), optional :: isSphere
type(ESMF_PoleKind_Flag), intent(in), optional :: polekindflag(2)
logical, intent(in), optional :: addCornerStagger
type(ESMF_TypeKind_Flag), intent(in), optional :: coordTypeKind
integer, intent(out), optional :: rc

Return Value type(ESMF_Grid)


Source Code

  function ESMF_GridCreateFrmGridspec(grid_filename, &
                                        regDecomp, indexflag, keywordEnforcer, decompflag, &
                                        addMask, varname, coordNames, &
                                        isSphere, polekindflag, &
                                        addCornerStagger, coordTypeKind, rc)

! !RETURN VALUE:
    type(ESMF_Grid) :: ESMF_GridCreateFrmGridspec
!
! !ARGUMENTS:
    character(len=*),      intent(in)             :: grid_filename
    integer,               intent(in)             :: regDecomp(:)
    type(ESMF_Index_Flag), intent(in)      :: indexflag
    type(ESMF_KeywordEnforcer), optional :: keywordEnforcer ! must use keywords below
    type(ESMF_Decomp_Flag), intent(in),   optional:: decompflag(:)
    logical,                intent(in),  optional  :: addMask
    character(len=*),       intent(in),  optional  :: varname
    character(len=*),       intent(in),  optional  :: coordNames(:)
    type(ESMF_PoleKind_Flag),  intent(in),  optional :: polekindflag(2)
    logical,               intent(in),  optional   :: isSphere
    logical,               intent(in),  optional   :: addCornerStagger
    type(ESMF_TypeKind_Flag),intent(in), optional  :: coordTypeKind
    integer,               intent(out), optional   :: rc

! !DESCRIPTION:
! This function creates a {\tt ESMF\_Grid} object using the grid definition from
! a GridSpec grid file.
! To specify the distribution, the user passes in an array
! ({\tt regDecomp}) specifying the number of DEs to divide each
! dimension into. The array {\tt decompflag} indicates how the division into DEs is to
! occur.  The default is to divide the range as evenly as possible.
! The grid defined in the file has to be a *** GENERALIZE  2D logically rectangular
! grid (i.e. {\tt grid\_rank} in the file needs to be 2 ***).
!
! This call is {\em collective} across the current VM.
!
! The arguments are:
! \begin{description}
! \item[{[grid_filename]}]
!      The GridSpec grid tile filename.
! \item[{[regDecomp]}]
!      A ndims-element array specifying how the grid is decomposed.
!      Each entry is the number of decounts for that dimension.
! \item[{[decompflag]}]
!      List of decomposition flags indicating how each dimension of the
!      tile is to be divided between the DEs. The default setting
!      is {\tt ESMF\_DECOMP\_BALANCED} in all dimensions. Please see
!      Section~\ref{opt:decompflag} for a full description of the
!      possible options. Note that currently the option
!      {\tt ESMF\_DECOMP\_CYCLIC} isn't supported in Grid creation.
! \item[{[addMask]}]
!      If .true., generate the mask using the missing value defined for varname
! \item[{[varname]}]
!      If addMask is true, provide a variable name stored in the grid file and
!      the mask will be generated using the missing value of the data value of
!      this variable.  The first two dimensions of the variable has to be the
!      the longitude and the latitude dimension and the mask is derived from the
!      first 2D values of this variable even if this data is 3D, or 4D array.
!\item[{[coordNames]}]
!      a two-element array containing the longitude and latitude variable names in a
!      GRIDSPEC file if there are multiple coordinates defined in the file
! \item[{[isSphere]}]
!      If .true., create a periodic Grid. If .false., create a regional Grid. Defaults to .true.
! \item[{[polekindflag]}]
!      Two item array which specifies the type of connection which occurs at the pole. The value in polekindflag(1)
!      specifies the connection that occurs at the minimum end of the pole dimension. The value in polekindflag(2)
!      specifies the connection that occurs at the maximum end of the pole dimension. Please see
!      Section~\ref{const:polekind} for a full list of options. If not specified,
!      the default is {\tt ESMF\_POLEKIND\_MONOPOLE} for both.
! \item[{[addCornerStagger]}]
!      Uses the information in the GridSpec file to add the Corner stagger to
!      the Grid. If not specified, defaults to true (since GridSpec defaults to
!      vertex-centered grids).
! \item[{[coordTypeKind]}]
!          The type/kind of the grid coordinate data. Only ESMF\_TYPEKIND\_R4
!          and ESMF\_TYPEKIND\_R8 are allowed.  
!          If not specified then defaults to ESMF\_TYPEKIND\_R8.
! \item[{[rc]}]
!      Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
! \end{description}
!
!EOP

#ifdef ESMF_NETCDF
    integer :: gridid, mosaicid, DimId, VarId, localrc, PetNo, PetCnt, ncStatus
    Integer :: dimids(2), coordids(2), gridims(2), datadims(2)
    integer :: ndims
    real(ESMF_KIND_R8),  allocatable :: loncoord1D(:), latcoord1D(:)
    real(ESMF_KIND_R8),  allocatable :: loncoord2D(:,:), latcoord2D(:,:)
    real(ESMF_KIND_R8),  allocatable :: cornerlon2D(:,:), cornerlat2D(:,:)
    real(ESMF_KIND_R8),  allocatable :: cornerlon3D(:,:,:), cornerlat3D(:,:,:)
    real(ESMF_KIND_R8),  allocatable :: corner1D(:), corner2D(:,:)
    real(ESMF_KIND_R4),  allocatable :: loncoord1DR4(:), latcoord1DR4(:)
    real(ESMF_KIND_R4),  allocatable :: loncoord2DR4(:,:), latcoord2DR4(:,:)
    real(ESMF_KIND_R4),  allocatable :: cornerlon2DR4(:,:), cornerlat2DR4(:,:)
    real(ESMF_KIND_R4),  allocatable :: cornerlon3DR4(:,:,:), cornerlat3DR4(:,:,:)
    real(ESMF_KIND_R4),  allocatable :: corner1DR4(:), corner2DR4(:,:)
    integer :: msgbuf(8)
    type(ESMF_CommHandle) :: commHandle
    integer :: localMinIndex(2), gridEdgeLWidth(2), gridEdgeUWidth(2)
    type(ESMF_Grid)  :: grid
    type(ESMF_Array) :: array
    type(ESMF_VM) :: vm
    type(ESMF_DistGrid) :: distgrid
    type(ESMF_Decomp_Flag):: decompflagLocal(2)
    type(ESMF_StaggerLoc) :: localStaggerLoc
    logical :: localAddCornerStagger, localIsSphere, localAddMask
    integer, allocatable :: dims(:), dims1(:)
    integer :: total(2), lbnd(2), ubnd(2), recv(1)
    real(ESMF_KIND_R8), allocatable :: varBuffer(:,:), recvbuf(:)
    real(ESMF_KIND_R8), pointer :: fptrlat(:,:), fptrlon(:,:)
    real(ESMF_KIND_R4), allocatable :: recvbufR4(:)
    real(ESMF_KIND_R4), pointer :: fptrlatR4(:,:), fptrlonR4(:,:)
    integer, pointer :: fptrmask(:,:), maskbuf(:)
    integer, allocatable :: mask2D(:,:)
    real(ESMF_KIND_R8) :: missing_value
    integer :: i,j,k,localroot
    integer :: maxIndex2D(2)
    integer, pointer :: minind(:,:)
    type(ESMF_CoordSys_Flag) :: coordsys
    character (len=256) :: units
    integer :: decnt
    integer, pointer                           :: minIndexPDe(:,:)
    integer, pointer                           :: maxIndexPDe(:,:)
    integer                                    :: start(2), count(2)
    logical                                    :: isGlobal
    logical                                    :: isSupergrid
    real(kind=ESMF_KIND_R8),  pointer          :: lonPtr(:,:), latPtr(:,:)
    real(kind=ESMF_KIND_R4),  pointer          :: lonPtrR4(:,:), latPtrR4(:,:)
    integer                                    :: localDe, deCount, s
    integer                                    :: sizex, sizey
    type(ESMF_StaggerLoc), allocatable         :: staggerLocList(:)
    type(ESMF_DELayout)                        :: delayout
    integer, allocatable                       :: demap(:)
    type(ESMF_TypeKind_Flag)                   :: localCoordTypeKind
    
   
    ! Initialize return code; assume failure until success is certain
    localrc = ESMF_RC_NOT_IMPL
    if (present(rc)) rc = ESMF_RC_NOT_IMPL

    ! get global vm information
    call ESMF_VMGetCurrent(vm, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return

    ! set up local pet info
    call ESMF_VMGet(vm, localPet=PetNo, petCount=PetCnt, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return

    if (present(decompFlag)) then
        decompFlagLocal(:)=decompFlag(:)
    else
        decompFlagLocal(:)=ESMF_DECOMP_BALANCED
    endif

    if (present(addCornerStagger)) then
        localAddCornerStagger=addCornerStagger
    else
        localAddCornerStagger=.false.
    endif

    if (present(isSphere)) then
        localIsSphere=isSphere
    else
        localIsSphere=.true.
    endif

    localAddMask = .false.
    if (present(addMask)) then
        if (addMask) then
          if (.not. present(varname)) then
             call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, &
                  msg="- need varname argument to create mask", &
                  ESMF_CONTEXT, rcToReturn=rc)
             return
          end if
          localAddMask = .true.
        endif
    endif
     
    if (present(coordTypeKind)) then
       if (coordTypeKind .ne. ESMF_TYPEKIND_R4 .and. &
           coordTypeKind .ne. ESMF_TYPEKIND_R8) then
          call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, &
             msg="- only ESMF_TYPEKIND_R4 and ESMF_TYPEKIND_R8 are allowed", &
             ESMF_CONTEXT, rcToReturn=rc)
          return
       endif      
       localCoordTypeKind=coordTypeKind
    else
       localCoordTypeKind=ESMF_TYPEKIND_R8
    endif

    call ESMF_GridspecQueryTileFile(grid_filename, isSupergrid, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
        ESMF_CONTEXT, rcToReturn=rc)) return

    if (isSupergrid) then
      call ESMF_GridspecQueryTileSize(grid_filename, sizex, sizey, rc=localrc)
      if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
        
      call ESMF_GridspecQueryTileGlobal(trim(grid_filename), isGlobal, rc=localrc)
      if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return

      if (isGlobal) then
        grid = ESMF_GridCreate1PeriDim(regDecomp, decompFlagLocal, &
          minIndex=(/1,1/), maxIndex=(/sizex,sizey/), &
          indexflag=indexflag, &
          coordSys=ESMF_COORDSYS_SPH_DEG, &
          coordTypeKind = localCoordTypeKind, &
          rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
      else
        grid = ESMF_GridCreateNoPeriDim(regDecomp, decompFlagLocal, &
          minIndex=(/1,1/), maxIndex=(/sizex,sizey/), &
          indexflag=indexflag, &
          coordSys=ESMF_COORDSYS_SPH_DEG, &
          coordTypeKind = localCoordTypeKind, &
          rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
      endif

      call ESMF_GridGet(grid, distgrid=distgrid, localDECount=decnt, rc=localrc)
      if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
         ESMF_CONTEXT, rcToReturn=rc)) return

      deCount = regDecomp(1)*regDecomp(2)
      allocate(minIndexPDe(2,deCount), maxIndexPDe(2,deCount))
      call ESMF_DistgridGet(distgrid, minIndexPDe=minIndexPDe, maxIndexPDe = maxIndexPDe, &
                          delayout=delayout, rc=localrc)
      if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
           ESMF_CONTEXT, rcToReturn=rc)) return
    
      allocate(demap(0:decnt-1))
      call ESMF_DELayoutGet(delayout, localDeToDeMap=demap, rc=localrc)
      if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
    
      if (localAddCornerStagger) then
         allocate(staggerLocList(2))
         staggerLocList(1) = ESMF_STAGGERLOC_CENTER
         staggerLocList(2) = ESMF_STAGGERLOC_CORNER
      else
         allocate(staggerLocList(1))
         staggerLocList(1) = ESMF_STAGGERLOC_CENTER
      endif      
      do s=1, size(staggerLocList)
         call ESMF_GridAddCoord(grid, staggerloc=staggerLocList(s), rc=localrc)
         if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
              ESMF_CONTEXT, rcToReturn=rc)) return
         if (localCoordTypeKind == ESMF_TYPEKIND_R8) then
           do localDe = 0,decnt-1
             call ESMF_GridGetCoord(grid, coordDim=1, localDe=localDe, &
                  staggerloc=staggerLocList(s), farrayPtr=lonPtr, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                  ESMF_CONTEXT, rcToReturn=rc)) return

             start(1)=minIndexPDe(1,demap(localDe)+1)
             start(2)=minIndexPDe(2,demap(localDe)+1)
             count=ubound(lonPtr)-lbound(lonPtr)+1
             call ESMF_GridGetCoord(grid, coordDim=2, localDe=localDe, &
                  staggerloc=staggerLocList(s), farrayPtr=latPtr, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                  ESMF_CONTEXT, rcToReturn=rc)) return
             ! Generate glocal edge coordinates and local center coordinates
             ! need to adjust the count???
             call ESMF_GridSpecReadStagger(trim(grid_filename),sizex, sizey, lonPtr, latPtr, &
                  staggerLoc=staggerLocList(s), &
                  start=start, count=count, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                  ESMF_CONTEXT, rcToReturn=rc)) return
           enddo
         else ! localCoordTypeKind == ESMF_TYPEKIND_R4
           do localDe = 0,decnt-1
             call ESMF_GridGetCoord(grid, coordDim=1, localDe=localDe, &
                  staggerloc=staggerLocList(s), farrayPtr=lonPtrR4, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                  ESMF_CONTEXT, rcToReturn=rc)) return

             start(1)=minIndexPDe(1,demap(localDe)+1)
             start(2)=minIndexPDe(2,demap(localDe)+1)
             count=ubound(lonPtrR4)-lbound(lonPtrR4)+1
             call ESMF_GridGetCoord(grid, coordDim=2, localDe=localDe, &
                  staggerloc=staggerLocList(s), farrayPtr=latPtrR4, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                  ESMF_CONTEXT, rcToReturn=rc)) return
             ! Generate glocal edge coordinates and local center coordinates
             ! need to adjust the count???
             call ESMF_GridSpecReadStagger(trim(grid_filename),sizex, sizey, lonPtrR4, latPtrR4, &
                  staggerLoc=staggerLocList(s), &
                  start=start, count=count, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                  ESMF_CONTEXT, rcToReturn=rc)) return
           enddo
         endif ! localCoordTypeKind == ESMF_TYPEKIND_RR
      enddo
      deallocate(minIndexPDe, maxIndexPDe, demap, staggerLocList)
    else  ! a regular CF Grid file containing center stagger 
   
      ! Get the grid rank and dimensions from the GridSpec file on PET 0, broadcast the
      ! data to all the PETs
      if (PetNo == 0) then
        call ESMF_GridspecInq(grid_filename, ndims, gridims, coord_names=coordNames, &
                dimids=dimids, coordids = coordids, units=units, rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                ESMF_CONTEXT, rcToReturn=rc)) return
        ! broadcast the values to other PETs (generalized)
        msgbuf(1)=ndims
        msgbuf(2:3) = gridims(:)
        msgbuf(4:5) = coordids(:)
        msgbuf(6:7) = dimids(:)
        if (trim(units) .eq. 'degrees') then
          msgbuf(8) = 0
          coordsys = ESMF_COORDSYS_SPH_DEG
        elseif (units(1:1) .eq. 'm') then
          msgbuf(8) = 1
          coordsys = ESMF_COORDSYS_CART
        elseif (units(1:1) .eq. 'k') then
          msgbuf(8) = 2
          coordsys = ESMF_COORDSYS_CART
        endif
        call ESMF_VMBroadcast(vm, msgbuf, 8, 0, rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
           ESMF_CONTEXT, rcToReturn=rc)) return
      else
        call ESMF_VMBroadcast(vm, msgbuf, 8, 0, rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
        ndims = msgbuf(1)
        gridims = msgbuf(2:3)
        coordids = msgbuf(4:5)
        dimids = msgbuf(6:7)
        if (msgbuf(8) == 0) then
          units = "degrees"
          coordsys = ESMF_COORDSYS_SPH_DEG
        elseif (msgbuf(8) == 1) then
          units = "meters"
          coordsys = ESMF_COORDSYS_CART
        elseif (msgbuf(8) == 2) then
          units = "kilometers"
          coordsys = ESMF_COORDSYS_CART
        endif
      endif

      if (localIsSphere .and. coordsys == ESMF_COORDSYS_CART) then
        call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, &
            msg="- A global grid cannot use Cartisian coordinates", &
            ESMF_CONTEXT, rcToReturn=rc)
        return
      endif      
      ! Create Grid based on the input distgrid
      gridEdgeLWidth=(/0,0/)
      if (localIsSphere) then
        gridEdgeUWidth=(/0,1/)
      else
        gridEdgeUWidth=(/1,1/)
      endif

      ! only parallelize the code if ndims == 2, so separate the code based on ndums
      if (ndims == 1) then
        if (localIsSphere) then
           grid = ESMF_GridCreate1PeriDim(minIndex=(/1,1/), maxIndex=gridims, &
                regDecomp=regDecomp, &
                gridEdgeLWidth=gridEdgeLWidth, gridEdgeUWidth=gridEdgeUWidth, &
                polekindflag=polekindflag, &
                coordDep1=(/1/), coordDep2=(/2/), &
                coordSys=ESMF_COORDSYS_SPH_DEG, &
                coordTypeKind = localCoordTypeKind, &
                indexflag=indexflag, rc=localrc)
           if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                ESMF_CONTEXT, rcToReturn=rc)) return
        else
            grid = ESMF_GridCreateNoPeriDim(minIndex=(/1,1/), maxIndex=gridims, &
                regDecomp=regDecomp, &
                gridEdgeLWidth=gridEdgeLWidth, gridEdgeUWidth=gridEdgeUWidth, &
                coordDep1=(/1/), coordDep2=(/2/), &
                coordSys=coordsys, &
                coordTypeKind = localCoordTypeKind, &
                indexflag=indexflag, rc=localrc)
           if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                ESMF_CONTEXT, rcToReturn=rc)) return
        endif
        
        if (localCoordTypeKind == ESMF_TYPEKIND_R8) then
          if (PetNo == 0) then
            ! Get coordinate info from the GridSpec file, if in radians, convert to degrees
            if (localAddCornerStagger) then
               allocate(loncoord1D(gridims(1)), latcoord1D(gridims(2)))
               allocate(cornerlon2D(2,gridims(1)), cornerlat2D(2, gridims(2)))
               call ESMF_GridspecGetVar1D(grid_filename, coordids, loncoord1D, latcoord1D,&
                                    cornerlon=cornerlon2D, cornerlat=cornerlat2D, rc=localrc)
               if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                  ESMF_CONTEXT, rcToReturn=rc)) return
            else
               allocate(loncoord1D(gridims(1)), latcoord1D(gridims(2)))
               call ESMF_GridspecGetVar1D(grid_filename, coordids, loncoord1D, latcoord1D,&
                                    rc=localrc)
               if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                 ESMF_CONTEXT, rcToReturn=rc)) return
            endif
            ! convert to kilometer if the units is "meters"
            if (units(1:1) .eq. 'm') then
               loncoord1D(:) = loncoord1D(:) * 1.d-3
               latcoord1D(:) = latcoord1D(:) * 1.d-3
               if (localAddCornerStagger) then
                  cornerlon2D(:,:) = cornerlon2D(:,:) * 1.d-3
                  cornerlat2D(:,:) = cornerlat2D(:,:) * 1.d-3
               endif
            endif
          endif
  
          ! Set coordinate tables -  Put Corners into coordinates
          localStaggerLoc = ESMF_STAGGERLOC_CENTER

          call ESMF_GridAddCoord(grid, staggerloc=localStaggerLoc, rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
              ESMF_CONTEXT, rcToReturn=rc)) return

          ! Set longitude coordinate
          call ESMF_GridGetCoord(grid, staggerloc=localStaggerLoc, coordDim=1, &
                     array=array, rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                     ESMF_CONTEXT, rcToReturn=rc)) return
          call ESMF_ArrayScatter(array, loncoord1D, rootPet=0, rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                     ESMF_CONTEXT, rcToReturn=rc)) return

          ! Set longitude coordinate
          call ESMF_GridGetCoord(grid, staggerloc=localStaggerLoc, coordDim=2, &
                     array=array, rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                     ESMF_CONTEXT, rcToReturn=rc)) return
          call ESMF_ArrayScatter(array, latcoord1D, rootPet=0, rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                     ESMF_CONTEXT, rcToReturn=rc)) return

          if (PetNo == 0) then
            deallocate(loncoord1D, latcoord1D)
          endif

          ! Add coordinates at the corner stagger location
          if (localAddCornerStagger) then
             call ESMF_GridAddCoord(grid, staggerloc=ESMF_STAGGERLOC_CORNER, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                ESMF_CONTEXT, rcToReturn=rc)) return

             ! Longitude
             call ESMF_GridGetCoord(grid, staggerloc=ESMF_STAGGERLOC_CORNER, coordDim=1, &
                  array = array, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                  ESMF_CONTEXT, rcToReturn=rc)) return

             call ESMF_GridGet(grid,tile=1,staggerloc=ESMF_STAGGERLOC_CORNER, maxIndex=maxIndex2D,&
                       rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                       ESMF_CONTEXT, rcToReturn=rc)) return
             if (PetNo == 0) then
                allocate(corner1D(maxIndex2D(1)))
                corner1D(1:gridims(1)) = cornerlon2D(1,:)
                if (maxIndex2D(1) > gridims(1)) then
                   corner1D(maxIndex2D(1)) = cornerlon2D(2,gridims(1))
                endif
             endif

             call ESMF_ArrayScatter(array, corner1D, rootPet=0, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                 ESMF_CONTEXT, rcToReturn=rc)) return

             if (PetNo == 0) deallocate(corner1D)

             ! Latitude
             call ESMF_GridGetCoord(grid, staggerloc=ESMF_STAGGERLOC_CORNER, coordDim=2, &
                 array = array, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                 ESMF_CONTEXT, rcToReturn=rc)) return

             if (PetNo == 0) then
                allocate(corner1D(maxIndex2D(2)))
                corner1D(1:gridims(2)) = cornerlat2D(1,:)
                if (maxIndex2D(2) > gridims(2)) then
                   corner1D(maxIndex2D(2)) = cornerlat2D(2,gridims(2))
                endif
             endif

             call ESMF_ArrayScatter(array, corner1D, rootPet=0, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                   ESMF_CONTEXT, rcToReturn=rc)) return
          endif
        else ! localCoordTypeKind == ESMF_TYPEKIND_R4
          if (PetNo == 0) then
            ! Get coordinate info from the GridSpec file, if in radians, convert to degrees
            if (localAddCornerStagger) then
              allocate(loncoord1DR4(gridims(1)), latcoord1DR4(gridims(2)))
              allocate(cornerlon2DR4(2,gridims(1)), cornerlat2DR4(2, gridims(2)))
              call ESMF_GridspecGetVar1DR4(grid_filename, coordids, loncoord1DR4, latcoord1DR4,&
                                    cornerlon=cornerlon2DR4, cornerlat=cornerlat2DR4, rc=localrc)
              if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                  ESMF_CONTEXT, rcToReturn=rc)) return
            else
              allocate(loncoord1DR4(gridims(1)), latcoord1DR4(gridims(2)))
              call ESMF_GridspecGetVar1DR4(grid_filename, coordids, loncoord1DR4, latcoord1DR4,&
                                    rc=localrc)
              if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                ESMF_CONTEXT, rcToReturn=rc)) return
            endif
            ! convert to kilometer if the units is "meters"
            if (units(1:1) .eq. 'm') then
              loncoord1DR4(:) = loncoord1DR4(:) * 1.d-3
              latcoord1DR4(:) = latcoord1DR4(:) * 1.d-3
              if (localAddCornerStagger) then
                 cornerlon2DR4(:,:) = cornerlon2DR4(:,:) * 1.d-3
                 cornerlat2DR4(:,:) = cornerlat2DR4(:,:) * 1.d-3
              endif
            endif
          endif
  
          ! Set coordinate tables -  Put Corners into coordinates
          localStaggerLoc = ESMF_STAGGERLOC_CENTER

          call ESMF_GridAddCoord(grid, staggerloc=localStaggerLoc, rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
              ESMF_CONTEXT, rcToReturn=rc)) return

          ! Set longitude coordinate
          call ESMF_GridGetCoord(grid, staggerloc=localStaggerLoc, coordDim=1, &
                     array=array, rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                     ESMF_CONTEXT, rcToReturn=rc)) return
          call ESMF_ArrayScatter(array, loncoord1DR4, rootPet=0, rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                     ESMF_CONTEXT, rcToReturn=rc)) return

          ! Set longitude coordinate
          call ESMF_GridGetCoord(grid, staggerloc=localStaggerLoc, coordDim=2, &
                     array=array, rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                     ESMF_CONTEXT, rcToReturn=rc)) return
          call ESMF_ArrayScatter(array, latcoord1DR4, rootPet=0, rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                     ESMF_CONTEXT, rcToReturn=rc)) return

          if (PetNo == 0) then
             deallocate(loncoord1DR4, latcoord1DR4)
          endif

          ! Add coordinates at the corner stagger location
          if (localAddCornerStagger) then
             call ESMF_GridAddCoord(grid, staggerloc=ESMF_STAGGERLOC_CORNER, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                ESMF_CONTEXT, rcToReturn=rc)) return

             ! Longitude
             call ESMF_GridGetCoord(grid, staggerloc=ESMF_STAGGERLOC_CORNER, coordDim=1, &
                  array = array, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                  ESMF_CONTEXT, rcToReturn=rc)) return

             call ESMF_GridGet(grid,tile=1,staggerloc=ESMF_STAGGERLOC_CORNER, maxIndex=maxIndex2D,&
                       rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                       ESMF_CONTEXT, rcToReturn=rc)) return
             if (PetNo == 0) then
                allocate(corner1DR4(maxIndex2D(1)))
                corner1DR4(1:gridims(1)) = cornerlon2DR4(1,:)
                if (maxIndex2D(1) > gridims(1)) then
                    corner1DR4(maxIndex2D(1)) = cornerlon2DR4(2,gridims(1))
                endif
             endif

             call ESMF_ArrayScatter(array, corner1DR4, rootPet=0, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                 ESMF_CONTEXT, rcToReturn=rc)) return

             if (PetNo == 0) deallocate(corner1DR4)

             ! Latitude
             call ESMF_GridGetCoord(grid, staggerloc=ESMF_STAGGERLOC_CORNER, coordDim=2, &
                  array = array, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                  ESMF_CONTEXT, rcToReturn=rc)) return

             if (PetNo == 0) then
                allocate(corner1DR4(maxIndex2D(2)))
                corner1DR4(1:gridims(2)) = cornerlat2DR4(1,:)
                if (maxIndex2D(2) > gridims(2)) then
                   corner1DR4(maxIndex2D(2)) = cornerlat2DR4(2,gridims(2))
                endif
             endif

             call ESMF_ArrayScatter(array, corner1DR4, rootPet=0, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                   ESMF_CONTEXT, rcToReturn=rc)) return
          endif !localAddCornerStagger
        endif  ! localCoordTypeKind == ESMF_TYPEKIND_R8
      elseif (ndims==2) then
        if (localIsSphere) then
           grid = ESMF_GridCreate1PeriDim(minIndex=(/1,1/), maxIndex=gridims, &
                  regDecomp=regDecomp, &
                  gridEdgeLWidth=gridEdgeLWidth, gridEdgeUWidth=gridEdgeUWidth, &
                  polekindflag=polekindflag, &
                  coordSys=ESMF_COORDSYS_SPH_DEG, &
                  coordTypeKind=localCoordTypeKind, &
                  indexflag=indexflag, rc=localrc)
           if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                  ESMF_CONTEXT, rcToReturn=rc)) return
        else
           grid = ESMF_GridCreateNoPeriDim(minIndex=(/1,1/), maxIndex=gridims, &
                  regDecomp=regDecomp, &
                  gridEdgeLWidth=gridEdgeLWidth, gridEdgeUWidth=gridEdgeUWidth, &
                  coordSys=coordsys, &
                  coordTypeKind=localCoordTypeKind, &
                  indexflag=indexflag, rc=localrc)
           if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                  ESMF_CONTEXT, rcToReturn=rc)) return
        endif

        ! Add center stagger coordinates
        localStaggerLoc = ESMF_STAGGERLOC_CENTER

        call ESMF_GridAddCoord(grid, staggerloc=localStaggerLoc, rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

        call ESMF_GridGet(grid, localDECount=decnt, rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
           ESMF_CONTEXT, rcToReturn=rc)) return
       
        if (localCoordTypeKind == ESMF_TYPEKIND_R8) then
          if (decnt > 0) then
            if ( mod(PetNo, regDecomp(1)) == 0) then
              call ESMF_GridGet(grid, distgrid=distgrid, rc=localrc)
              if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                 ESMF_CONTEXT, rcToReturn=rc)) return
              allocate(minind(2,PetCnt))
              call ESMF_DistGridGet(distgrid, minIndexPDe=minind, rc=localrc)
              if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                 ESMF_CONTEXT, rcToReturn=rc)) return
              call ESMF_GridGet(grid, ESMF_STAGGERLOC_CENTER, 0, exclusiveLBound=lbnd, &
                 exclusiveUBound=ubnd, exclusiveCount=total, rc=localrc)
              if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                 ESMF_CONTEXT, rcToReturn=rc)) return

              total(1)=gridims(1)
              allocate(loncoord2D(total(1),total(2)), latcoord2D(total(1),total(2)))
              if (localAddCornerStagger) then
                 allocate(cornerlon3D(4,total(1),total(2)), cornerlat3D(4, total(1), total(2)))
                 call ESMF_GridspecGetVar2D(grid_filename, coordids, &
                                    loncoord=loncoord2D, latcoord=latcoord2D, &
                                    cornerlon=cornerlon3D, cornerlat=cornerlat3D, &
                                    start=minind(:,PetNo+1), count=total, &
                                    rc=localrc)
                 if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                      ESMF_CONTEXT, rcToReturn=rc)) return
              else
                 call ESMF_GridspecGetVar2D(grid_filename, coordids,  &
                                    loncoord=loncoord2D, latcoord=latcoord2D, &
                                    start=minind(:,PetNo+1), count=total, &
                                    rc=localrc)
                 if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                      ESMF_CONTEXT, rcToReturn=rc)) return
              endif !localAddCornerStagger
              ! convert to kilometer if the units is "meters"
              if (units(1:1) .eq. 'm') then
                 loncoord2D(:,:) = loncoord2D(:,:) * 1.d-3
                 latcoord2D(:,:) = latcoord2D(:,:) * 1.d-3
                 if (localAddCornerStagger) then
                   cornerlon3D(:,:,:) = cornerlon3D(:,:,:) * 1.d-3
                   cornerlat3D(:,:,:) = cornerlat3D(:,:,:) * 1.d-3
                 endif
              endif
            endif ! mod(PetNo, regDecomp(1))==0

            ! Set longitude coordinate
            call ESMF_GridGetCoord(grid, staggerloc=localStaggerLoc, coordDim=1, &
                   farrayptr=fptrlon, rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                   ESMF_CONTEXT, rcToReturn=rc)) return

            ! Set latitude coordinate
            call ESMF_GridGetCoord(grid, staggerloc=localStaggerLoc, coordDim=2, &
                   farrayptr=fptrlat, rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                   ESMF_CONTEXT, rcToReturn=rc)) return

            if (mod(PetNo, regDecomp(1)) == 0) then
               allocate(dims(regdecomp(1)-1))
               do i=1, regDecomp(1)-1
                 call ESMF_VMRecv(vm, recv, 1, PetNo+i)
                 if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                     ESMF_CONTEXT, rcToReturn=rc)) return
                 dims(i)=recv(1)
              enddo
              call pack_and_send_float2D(vm, total, regDecomp(1), PetNo, lonCoord2D, fptrlon, dims)
              call pack_and_send_float2D(vm, total, regDecomp(1), PetNo, latCoord2D, fptrlat, dims)

              deallocate(loncoord2D, latcoord2D, dims)
            else
              localroot = (PetNo/regDecomp(1))*regDecomp(1)
              call ESMF_GridGet(grid, ESMF_STAGGERLOC_CENTER, 0, exclusiveLBound=lbnd, &
                 exclusiveUBound=ubnd, exclusiveCount=total, rc=localrc)
              if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                 ESMF_CONTEXT, rcToReturn=rc)) return
              allocate(recvbuf(total(1)*total(2)))
              ! First, send the xdim of the local array to the localroot
              call ESMF_VMSend(vm, total, 1, localroot, rc=localrc)
              if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                 ESMF_CONTEXT, rcToReturn=rc)) return

              ! Longitude coordinates
              call ESMF_VMRecv(vm, recvbuf, total(1)*total(2), localroot, rc=localrc)
              if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                 ESMF_CONTEXT, rcToReturn=rc)) return
              k=1
              do i=lbnd(2),ubnd(2)
                do j=lbnd(1),ubnd(1)
                  fptrlon(j,i) = recvbuf(k)
                  k=k+1
                enddo
              enddo
              ! Latitude coordinates
              call ESMF_VMRecv(vm, recvbuf, total(1)*total(2), localroot, rc=localrc)
              if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                  ESMF_CONTEXT, rcToReturn=rc)) return
              k=1
              do i=lbnd(2),ubnd(2)
                do j=lbnd(1),ubnd(1)
                  fptrlat(j,i) = recvbuf(k)
                  k=k+1
                enddo
              enddo
              deallocate(recvbuf)
            endif ! mod(PetNo, regDecomp(1))==0
          endif !decnt > 0

          ! Add coordinates at the corner stagger location
          if (localAddCornerStagger) then
             call ESMF_GridGet(grid,tile=1,staggerloc=ESMF_STAGGERLOC_CORNER, &
                    maxIndex=maxIndex2D, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                  ESMF_CONTEXT, rcToReturn=rc)) return

             call ESMF_GridAddCoord(grid, staggerloc=ESMF_STAGGERLOC_CORNER, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                 ESMF_CONTEXT, rcToReturn=rc)) return

             ! Longitude
             call ESMF_GridGetCoord(grid, staggerloc=ESMF_STAGGERLOC_CORNER, coordDim=1, &
                 farrayptr = fptrlon, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                ESMF_CONTEXT, rcToReturn=rc)) return

             ! Latitude
             call ESMF_GridGetCoord(grid, staggerloc=ESMF_STAGGERLOC_CORNER, coordDim=2, &
                  farrayptr = fptrlat, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                  ESMF_CONTEXT, rcToReturn=rc)) return

             if (decnt > 0) then
               if (mod(PetNo, regDecomp(1)) == 0) then
                 ! Get the x dimension of every local array
                 allocate(dims1(regDecomp(1)-1))
                 do i=1, regDecomp(1)-1
                   call ESMF_VMRecv(vm, recv, 1, PetNo+i)
                   if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                     ESMF_CONTEXT, rcToReturn=rc)) return
                   dims1(i)=recv(1)
                 enddo

                 call ESMF_GridGet(grid, ESMF_STAGGERLOC_CORNER, 0, exclusiveLBound=lbnd, &
                    exclusiveUBound=ubnd, exclusiveCount=total, rc=localrc)
                 if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                    ESMF_CONTEXT, rcToReturn=rc)) return
                 ! prepare the corner coordinates for the entire slab to send to the
                 ! member PETs, the horizontal size is the maxIndex for the corner stagger
                 total(1)=maxIndex2D(1)
                 allocate(corner2D(total(1),total(2)))
                 ! the input cornerlon3D may be smaller than the size of corner2D, so need to
                 ! fill the data at the right border or the bottom border (if it is the last
                 ! row of the grid
                 datadims(1)=size(cornerlon3D,2)
                 datadims(2)=size(cornerlon3D,3)
                 corner2D(1:datadims(1),1:datadims(2)) = cornerlon3D(1,:,:)
                 if (total(1) > datadims(1)) then
                    corner2D(total(1), 1:datadims(2)) = cornerlon3D(2,datadims(1),:)
                 end if
                 if (total(2) > datadims(2)) then
                    corner2D(1:datadims(1), total(2)) = cornerlon3D(4,:,datadims(2))
                 end if
                 if (total(1) > datadims(1) .and. total(2) > datadims(2)) then
                    corner2D(total(1),total(2))=cornerlon3D(3, datadims(1), datadims(2))
                 endif

                 call pack_and_send_float2D(vm, total, regDecomp(1), PetNo, corner2D, fptrlon, dims1)

                 corner2D(1:datadims(1),1:datadims(2)) = cornerlat3D(1,:,:)
                 if (total(1) > datadims(1)) then
                    corner2D(total(1), 1:datadims(2)) = cornerlat3D(2,datadims(1),:)
                 end if
                 if (total(2) > datadims(2)) then
                   corner2D(1:datadims(1), total(2)) = cornerlat3D(4,:,datadims(2))
                 end if
                 if (total(1) > datadims(1) .and. total(2) > datadims(2)) then
                   corner2D(total(1),total(2))=cornerlat3D(3, datadims(1), datadims(2))
                 endif
                 call pack_and_send_float2D(vm, total, regDecomp(1), PetNo, corner2D, fptrlat, dims1)
                 deallocate(dims1, minind)
               else ! mod(PetNo, regDecomp(1)) == 0
                 call ESMF_GridGet(grid, ESMF_STAGGERLOC_CORNER, 0, exclusiveLBound=lbnd, &
                    exclusiveUBound=ubnd, exclusiveCount=total, rc=localrc)
                 if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                   ESMF_CONTEXT, rcToReturn=rc)) return
                 allocate(recvbuf(total(1)*total(2)))

                 ! First, send the xdim of the local array to the localroot
                 call ESMF_VMSend(vm, total, 1, localroot, rc=localrc)
                 if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                    ESMF_CONTEXT, rcToReturn=rc)) return

                 ! Longitude coordinates
                 call ESMF_VMRecv(vm, recvbuf, total(1)*total(2), localroot, rc=localrc)
                 if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                    ESMF_CONTEXT, rcToReturn=rc)) return
                 k=1
                 do i=lbnd(2),ubnd(2)
                   do j=lbnd(1),ubnd(1)
                     fptrlon(j,i) = recvbuf(k)
                     k=k+1
                   enddo
                 enddo
                 ! Latitude coordinates
                 call ESMF_VMRecv(vm, recvbuf, total(1)*total(2), localroot, rc=localrc)
                 if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                     ESMF_CONTEXT, rcToReturn=rc)) return
                 k=1
                 do i=lbnd(2),ubnd(2)
                   do j=lbnd(1),ubnd(1)
                     fptrlat(j,i) = recvbuf(k)
                     k=k+1
                   enddo
                 enddo
                 deallocate(recvbuf)
               endif  ! end if (mod(PetNo, RegDecomp(1))==0)
             endif  ! decnt > 0
           endif  ! end if (AddCornerStagger)
         else  ! localCoordTypeKind == ESMF_TYPEKIND_R4
          if (decnt > 0) then
            if ( mod(PetNo, regDecomp(1)) == 0) then
              call ESMF_GridGet(grid, distgrid=distgrid, rc=localrc)
              if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                 ESMF_CONTEXT, rcToReturn=rc)) return
              allocate(minind(2,PetCnt))
              call ESMF_DistGridGet(distgrid, minIndexPDe=minind, rc=localrc)
              if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                 ESMF_CONTEXT, rcToReturn=rc)) return
              call ESMF_GridGet(grid, ESMF_STAGGERLOC_CENTER, 0, exclusiveLBound=lbnd, &
                 exclusiveUBound=ubnd, exclusiveCount=total, rc=localrc)
              if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                 ESMF_CONTEXT, rcToReturn=rc)) return

              total(1)=gridims(1)
              allocate(loncoord2DR4(total(1),total(2)), latcoord2DR4(total(1),total(2)))
              if (localAddCornerStagger) then
                allocate(cornerlon3DR4(4,total(1),total(2)), cornerlat3DR4(4, total(1), total(2)))
                call ESMF_GridspecGetVar2DR4(grid_filename, coordids, &
                                    loncoord=loncoord2DR4, latcoord=latcoord2DR4, &
                                    cornerlon=cornerlon3DR4, cornerlat=cornerlat3DR4, &
                                    start=minind(:,PetNo+1), count=total, &
                                    rc=localrc)
                if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                      ESMF_CONTEXT, rcToReturn=rc)) return
              else
                call ESMF_GridspecGetVar2DR4(grid_filename, coordids,  &
                                    loncoord=loncoord2DR4, latcoord=latcoord2DR4, &
                                    start=minind(:,PetNo+1), count=total, &
                                    rc=localrc)
                if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                      ESMF_CONTEXT, rcToReturn=rc)) return
              endif
              ! convert to kilometer if the units is "meters"
              if (units(1:1) .eq. 'm') then
                loncoord2DR4(:,:) = loncoord2DR4(:,:) * 1.d-3
                latcoord2DR4(:,:) = latcoord2DR4(:,:) * 1.d-3
                if (localAddCornerStagger) then
                  cornerlon3DR4(:,:,:) = cornerlon3DR4(:,:,:) * 1.d-3
                  cornerlat3DR4(:,:,:) = cornerlat3DR4(:,:,:) * 1.d-3
                endif
              endif
            endif ! mod(PetNo, regDecomp(1))==0
        
            ! Set longitude coordinate
            call ESMF_GridGetCoord(grid, staggerloc=localStaggerLoc, coordDim=1, &
                     farrayptr=fptrlonR4, rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                     ESMF_CONTEXT, rcToReturn=rc)) return

            ! Set latitude coordinate
            call ESMF_GridGetCoord(grid, staggerloc=localStaggerLoc, coordDim=2, &
                     farrayptr=fptrlatR4, rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                     ESMF_CONTEXT, rcToReturn=rc)) return

           if (mod(PetNo, regDecomp(1)) == 0) then
             allocate(dims(regdecomp(1)-1))
             do i=1, regDecomp(1)-1
               call ESMF_VMRecv(vm, recv, 1, PetNo+i)
               if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                   ESMF_CONTEXT, rcToReturn=rc)) return
               dims(i)=recv(1)
             enddo
             call pack_and_send_float2DR4(vm, total, regDecomp(1), PetNo, lonCoord2DR4, fptrlonR4, dims)
             call pack_and_send_float2DR4(vm, total, regDecomp(1), PetNo, latCoord2DR4, fptrlatR4, dims)

             deallocate(loncoord2DR4, latcoord2DR4, dims)
           else
             localroot = (PetNo/regDecomp(1))*regDecomp(1)
             call ESMF_GridGet(grid, ESMF_STAGGERLOC_CENTER, 0, exclusiveLBound=lbnd, &
                exclusiveUBound=ubnd, exclusiveCount=total, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                ESMF_CONTEXT, rcToReturn=rc)) return
             allocate(recvbufR4(total(1)*total(2)))
             ! First, send the xdim of the local array to the localroot
             call ESMF_VMSend(vm, total, 1, localroot, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                 ESMF_CONTEXT, rcToReturn=rc)) return

             ! Longitude coordinates
             call ESMF_VMRecv(vm, recvbufR4, total(1)*total(2), localroot, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                 ESMF_CONTEXT, rcToReturn=rc)) return
             k=1
             do i=lbnd(2),ubnd(2)
               do j=lbnd(1),ubnd(1)
                   fptrlonR4(j,i) = recvbufR4(k)
                   k=k+1
               enddo
             enddo
             ! Latitude coordinates
             call ESMF_VMRecv(vm, recvbufR4, total(1)*total(2), localroot, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                ESMF_CONTEXT, rcToReturn=rc)) return
             k=1
             do i=lbnd(2),ubnd(2)
               do j=lbnd(1),ubnd(1)
                 fptrlatR4(j,i) = recvbufR4(k)
                 k=k+1
               enddo
             enddo
             deallocate(recvbufR4)
            endif  ! end if (mod(PetNo, RegDecomp(1))==0)
          endif  ! decnt > 0

          ! Add coordinates at the corner stagger location
          if (localAddCornerStagger) then
             call ESMF_GridGet(grid,tile=1,staggerloc=ESMF_STAGGERLOC_CORNER, &
                    maxIndex=maxIndex2D, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                  ESMF_CONTEXT, rcToReturn=rc)) return

             call ESMF_GridAddCoord(grid, staggerloc=ESMF_STAGGERLOC_CORNER, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                  ESMF_CONTEXT, rcToReturn=rc)) return

             ! Longitude
             call ESMF_GridGetCoord(grid, staggerloc=ESMF_STAGGERLOC_CORNER, coordDim=1, &
                 farrayptr = fptrlonR4, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                ESMF_CONTEXT, rcToReturn=rc)) return

             ! Latitude
             call ESMF_GridGetCoord(grid, staggerloc=ESMF_STAGGERLOC_CORNER, coordDim=2, &
                  farrayptr = fptrlatR4, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                  ESMF_CONTEXT, rcToReturn=rc)) return

             if (decnt > 0) then
               if (mod(PetNo, regDecomp(1)) == 0) then
                 ! Get the x dimension of every local array
                 allocate(dims1(regDecomp(1)-1))
                 do i=1, regDecomp(1)-1
                   call ESMF_VMRecv(vm, recv, 1, PetNo+i)
                   if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                     ESMF_CONTEXT, rcToReturn=rc)) return
                   dims1(i)=recv(1)
                 enddo

               call ESMF_GridGet(grid, ESMF_STAGGERLOC_CORNER, 0, exclusiveLBound=lbnd, &
                   exclusiveUBound=ubnd, exclusiveCount=total, rc=localrc)
               if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                    ESMF_CONTEXT, rcToReturn=rc)) return
               ! prepare the corner coordinates for the entire slab to send to the
               ! member PETs, the horizontal size is the maxIndex for the corner stagger
               total(1)=maxIndex2D(1)
               allocate(corner2DR4(total(1),total(2)))
               ! the input cornerlon3D may be smaller than the size of corner2D, so need to
               ! fill the data at the right border or the bottom border (if it is the last
               ! row of the grid
               datadims(1)=size(cornerlon3DR4,2)
               datadims(2)=size(cornerlon3DR4,3)
               corner2DR4(1:datadims(1),1:datadims(2)) = cornerlon3DR4(1,:,:)
               if (total(1) > datadims(1)) then
                    corner2DR4(total(1), 1:datadims(2)) = cornerlon3DR4(2,datadims(1),:)
               end if
               if (total(2) > datadims(2)) then
                   corner2DR4(1:datadims(1), total(2)) = cornerlon3DR4(4,:,datadims(2))
               end if
               if (total(1) > datadims(1) .and. total(2) > datadims(2)) then
                   corner2DR4(total(1),total(2))=cornerlon3DR4(3, datadims(1), datadims(2))
               endif

               call pack_and_send_float2DR4(vm, total, regDecomp(1), PetNo, corner2DR4, fptrlonR4, dims1)

               corner2DR4(1:datadims(1),1:datadims(2)) = cornerlat3DR4(1,:,:)
               if (total(1) > datadims(1)) then
                    corner2DR4(total(1), 1:datadims(2)) = cornerlat3DR4(2,datadims(1),:)
               end if
               if (total(2) > datadims(2)) then
                   corner2DR4(1:datadims(1), total(2)) = cornerlat3DR4(4,:,datadims(2))
               end if
               if (total(1) > datadims(1) .and. total(2) > datadims(2)) then
                   corner2DR4(total(1),total(2))=cornerlat3DR4(3, datadims(1), datadims(2))
               endif
               call pack_and_send_float2DR4(vm, total, regDecomp(1), PetNo, corner2DR4, fptrlatR4, dims1)
               deallocate(dims1, minind)
             else ! mod(PetNo, regDecomp(1))==0
               call ESMF_GridGet(grid, ESMF_STAGGERLOC_CORNER, 0, exclusiveLBound=lbnd, &
                  exclusiveUBound=ubnd, exclusiveCount=total, rc=localrc)
               if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                   ESMF_CONTEXT, rcToReturn=rc)) return
               allocate(recvbufR4(total(1)*total(2)))

               ! First, send the xdim of the local array to the localroot
               call ESMF_VMSend(vm, total, 1, localroot, rc=localrc)
               if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                  ESMF_CONTEXT, rcToReturn=rc)) return

               ! Longitude coordinates
               call ESMF_VMRecv(vm, recvbufR4, total(1)*total(2), localroot, rc=localrc)
               if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                    ESMF_CONTEXT, rcToReturn=rc)) return
                 k=1
                 do i=lbnd(2),ubnd(2)
                   do j=lbnd(1),ubnd(1)
                     fptrlonR4(j,i) = recvbufR4(k)
                     k=k+1
                   enddo
                 enddo
                 ! Latitude coordinates
                 call ESMF_VMRecv(vm, recvbufR4, total(1)*total(2), localroot, rc=localrc)
                 if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                     ESMF_CONTEXT, rcToReturn=rc)) return
                 k=1
                 do i=lbnd(2),ubnd(2)
                   do j=lbnd(1),ubnd(1)
                     fptrlatR4(j,i) = recvbufR4(k)
                     k=k+1
                   enddo
                 enddo
                 deallocate(recvbufR4)
               endif  ! end if (mod(PetNo, RegDecomp(1))==0)
             endif    ! end if (decnt>0)
           endif      ! end if (localAddCornerStagger) then
        endif         ! end if (localCoordTypeKind == ESMF_TYPEKIND_R8)
      endif           ! end if ndims == 2

      ! Only add mask if localAddMask = .TRUE.
      ! This code is common whether it is ndims=1 or ndims=2
      if (localAddMask) then
        call ESMF_GridAddItem(grid, staggerloc=localStaggerLoc, &
                        itemflag = ESMF_GRIDITEM_MASK, rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                        ESMF_CONTEXT, rcToReturn=rc)) return
        call ESMF_GridGetItem(grid, staggerloc=localStaggerLoc,  &
                        itemflag=ESMF_GRIDITEM_MASK, farrayptr=fptrmask, rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                        ESMF_CONTEXT, rcToReturn=rc)) return
        call ESMF_GridGet(grid, localStaggerLoc, 0, exclusiveLBound=lbnd, &
                  exclusiveUBound=ubnd, exclusiveCount=total, rc=localrc)

        ! Check if we want to extract mask from a data variable
        if (mod(PetNo, regDecomp(1)) == 0) then
           call ESMF_GridGet(grid, distgrid=distgrid, rc=localrc)
           if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                ESMF_CONTEXT, rcToReturn=rc)) return
           allocate(minind(2,PetCnt))
           call ESMF_DistGridGet(distgrid, minIndexPDe=minind, rc=localrc)
           if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                ESMF_CONTEXT, rcToReturn=rc)) return
           allocate(dims(regdecomp(1)-1))
           do i=1, regDecomp(1)-1
              call ESMF_VMRecv(vm, recv, 1, PetNo+i)
              if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                  ESMF_CONTEXT, rcToReturn=rc)) return
              dims(i)=recv(1)
           enddo
           call ESMF_GridGet(grid, ESMF_STAGGERLOC_CENTER, 0, exclusiveLBound=lbnd, &
              exclusiveUBound=ubnd, exclusiveCount=total, rc=localrc)
           total(1)=gridims(1)
           allocate(varBuffer(total(1),total(2)),mask2D(total(1),total(2)))
           mask2D(:,:) = 1
           call ESMF_GridspecGetVarByName(grid_filename, varname, dimids, &
                                varBuffer, missing_value = missing_value, &
                                start=minind(:,PetNo+1), count=total, rc=localrc)
           if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                ESMF_CONTEXT, rcToReturn=rc)) return
           do i=1,size(varBuffer,2)
             do j=1,size(varBuffer,1)
               if (varBuffer(j,i) == missing_value) then
                 mask2D(j,i)=0
               endif
             enddo
           enddo
           call pack_and_send_int2D(vm, total, regDecomp(1), PetNo, mask2D, fptrmask, dims)
           deallocate(varBuffer)
           deallocate(mask2D)
           deallocate(dims)
           deallocate(minind)
        else
           call ESMF_GridGet(grid, ESMF_STAGGERLOC_CENTER, 0, exclusiveLBound=lbnd, &
              exclusiveUBound=ubnd, exclusiveCount=total, rc=localrc)

           localroot = (PetNo/regDecomp(1))*regDecomp(1)
           ! First, send the xdim of the local array to the localroot
           call ESMF_VMSend(vm, total, 1, localroot)
           if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
              ESMF_CONTEXT, rcToReturn=rc)) return
           allocate(maskbuf(total(1)*total(2)))
           call ESMF_VMRecv(vm, maskbuf, total(1)*total(2), localroot, rc=localrc)
           if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                ESMF_CONTEXT, rcToReturn=rc)) return
           k=1
           do i=lbnd(2),ubnd(2)
              do j=lbnd(1),ubnd(1)
                 fptrmask(j,i) = maskbuf(k)
                 k=k+1
              enddo
           enddo
           deallocate(maskbuf)
        endif
      endif
    endif ! if (isSuperGrid)

    ESMF_GridCreateFrmGridspec = grid

    if (present(rc)) rc=ESMF_SUCCESS
    return
#else
    call ESMF_LogSetError(rcToCheck=ESMF_RC_LIB_NOT_PRESENT, &
                 msg="- ESMF_NETCDF not defined when lib was compiled", &
                 ESMF_CONTEXT, rcToReturn=rc)
    return
#endif

end function ESMF_GridCreateFrmGridspec