ESMF_XGridGetDefault Subroutine

private subroutine ESMF_XGridGetDefault(xgrid, keywordEnforcer, sideAGridCount, sideBGridCount, sideAMeshCount, sideBMeshCount, coordSys, dimCount, elementCount, sideAGrid, sideBGrid, sideAMesh, sideBMesh, mesh, area, centroid, distgridA, distgridB, distgridM, sparseMatA2X, sparseMatX2A, sparseMatB2X, sparseMatX2B, name, rc)

Arguments

Type IntentOptional Attributes Name
type(ESMF_XGrid), intent(in) :: xgrid
type(ESMF_KeywordEnforcer), optional :: keywordEnforcer
integer, intent(out), optional :: sideAGridCount
integer, intent(out), optional :: sideBGridCount
integer, intent(out), optional :: sideAMeshCount
integer, intent(out), optional :: sideBMeshCount
type(ESMF_CoordSys_Flag), intent(out), optional :: coordSys
integer, intent(out), optional :: dimCount
integer, intent(out), optional :: elementCount
type(ESMF_Grid), intent(out), optional :: sideAGrid(:)
type(ESMF_Grid), intent(out), optional :: sideBGrid(:)
type(ESMF_Mesh), intent(out), optional :: sideAMesh(:)
type(ESMF_Mesh), intent(out), optional :: sideBMesh(:)
type(ESMF_Mesh), intent(out), optional :: mesh
real(kind=ESMF_KIND_R8), intent(out), optional :: area(:)
real(kind=ESMF_KIND_R8), intent(out), optional :: centroid(:,:)
type(ESMF_DistGrid), intent(out), optional :: distgridA(:)
type(ESMF_DistGrid), intent(out), optional :: distgridB(:)
type(ESMF_DistGrid), intent(out), optional :: distgridM
type(ESMF_XGridSpec), intent(out), optional :: sparseMatA2X(:)
type(ESMF_XGridSpec), intent(out), optional :: sparseMatX2A(:)
type(ESMF_XGridSpec), intent(out), optional :: sparseMatB2X(:)
type(ESMF_XGridSpec), intent(out), optional :: sparseMatX2B(:)
character(len=*), intent(out), optional :: name
integer, intent(out), optional :: rc

Source Code

subroutine ESMF_XGridGetDefault(xgrid, keywordEnforcer, &
    sideAGridCount, sideBGridCount, sideAMeshCount, sideBMeshCount, &
    coordSys, &
    dimCount, elementCount, &
    sideAGrid, sideBGrid, sideAMesh, sideBMesh, &
    mesh, &
    area, centroid, &
    distgridA, distgridB, distgridM, &
    sparseMatA2X, sparseMatX2A, sparseMatB2X, sparseMatX2B, &
    name, &
    rc) 

!
! !ARGUMENTS:
type(ESMF_XGrid),     intent(in)            :: xgrid
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
integer,              intent(out), optional :: sideAGridCount, sideBGridCount
integer,              intent(out), optional :: sideAMeshCount, sideBMeshCount
type(ESMF_CoordSys_Flag), intent(out), optional :: coordSys
integer,              intent(out), optional :: dimCount
integer,              intent(out), optional :: elementCount
type(ESMF_Grid),      intent(out), optional :: sideAGrid(:), sideBGrid(:)
type(ESMF_Mesh),      intent(out), optional :: sideAMesh(:), sideBMesh(:)
type(ESMF_Mesh),      intent(out), optional :: mesh
real(ESMF_KIND_R8),   intent(out), optional :: area(:)
real(ESMF_KIND_R8),   intent(out), optional :: centroid(:,:)
type(ESMF_DistGrid),  intent(out), optional :: distgridA(:)
type(ESMF_DistGrid),  intent(out), optional :: distgridB(:)
type(ESMF_DistGrid),  intent(out), optional :: distgridM
type(ESMF_XGridSpec), intent(out), optional :: sparseMatA2X(:)
type(ESMF_XGridSpec), intent(out), optional :: sparseMatX2A(:)
type(ESMF_XGridSpec), intent(out), optional :: sparseMatB2X(:)
type(ESMF_XGridSpec), intent(out), optional :: sparseMatX2B(:)
character (len=*),    intent(out), optional :: name
integer,              intent(out), optional :: rc 
!
! !DESCRIPTION:
!      Get information about XGrid
!
!     The arguments are:
!     \begin{description}
!     \item [xgrid]
!       The {\tt ESMF\_XGrid} object used to retrieve information from.
!     \item [{[sideAGridCount]}]
!           Total Number of Grids on the A side.
!     \item [{[sideBGridCount]}]
!           Total Number of Grids on the B side.
!     \item [{[sideAMeshCount]}]
!           Total Number of Meshes on the A side.
!     \item [{[sideBMeshCount]}]
!           Total Number of Meshes on the B side.
!     \item[{[coordSys]}]
!           The coordinate system of the XGrid's coordinate data.
!     \item [{[dimCount]}]
!           Number of dimension of the xgrid.
!     \item [{[elementCount]}]
!          Number of elements in exclusive region of the xgrid on this PET.
!     \item [{[sideAGrid]}]
!           List of 2D Grids on side A. Must enter with shape(sideAGrid)=(/sideAGridCount/).
!     \item [{[sideBGrid]}]
!           List of 2D Grids on side B. Must enter with shape(sideBGrid)=(/sideBGridCount/).
!     \item [{[sideAMesh]}]
!           List of 2D Meshes on side A. Must enter with shape(sideAMesh)=(/sideAMeshCount/).
!     \item [{[sideBMesh]}]
!           List of 2D Meshes on side B. Must enter with shape(sideBMesh)=(/sideBMeshCount/).
!     \item [{[mesh]}]
!           Super mesh stored in XGrid when storeOverlay is set true during XGrid creation.
!     \item [{[area]}]
!           Area of the xgrid cells on this PET. Must enter with shape(area)=(/elementCount/).
!     \item [{[centroid]}]
!           Coordinates at the area weighted center of the xgrid cells on this PET. Must enter with shape(centroid)=(/elementCount, dimCount/).
!     \item [{[distgridA]}]
!           List of distgrids whose sequence index list is an overlap between a Grid
!           on sideA and the xgrid object. Must enter with shape(distgridA)=(/sideAGridCount+sideAMeshCount/).
!     \item [{[distgridB]}]
!           List of distgrids whose sequence index list is an overlap between a Grid
!           on sideB and the xgrid object. Must enter with shape(distgridB)=(/sideBGridCount+sideBMeshCount/).
!     \item [{[distgridM]}]
!           The distgrid whose sequence index list fully describes the xgrid object.
!     \item [{[sparseMatA2X]}]
!           Indexlist from a Grid index space on side A to xgrid index space; 
!           indexFactorlist from a Grid index space on side A to xgrid index space. Must enter with shape(sparsematA2X)=(/sideAGridCount+sideAMeshCount/).
!     \item [{[sparseMatX2A]}]
!           Indexlist from xgrid index space to a Grid index space on side A; 
!           indexFactorlist from xgrid index space to a Grid index space on side A. Must enter with shape(sparsematX2A)=(/sideAGridCount+sideAMeshCount/).
!     \item [{[sparseMatB2X]}]
!           Indexlist from a Grid index space on side B to xgrid index space; 
!           indexFactorlist from a Grid index space on side B to xgrid index space. Must enter with shape(sparsematB2X)=(/sideBGridCount+sideBMeshCount/).
!     \item [{[sparseMatX2B]}]
!           Indexlist from xgrid index space to a Grid index space on side B; 
!           indexFactorlist from xgrid index space to a Grid index space on side B. Must enter with shape(sparsematX2B)=(/sideBGridCount+sideBMeshCount/).
!     \item [{[name]}]
!           Name of the xgrid object.
!     \item [{[rc]}]
!           Return code; equals {\tt ESMF\_SUCCESS} only if the {\tt ESMF\_XGrid} 
!           is created.
!     \end{description}
!
!EOP

    integer :: localrc, ngrid_a, ngrid_b, n_idx_a2x, n_idx_x2a, n_idx_b2x, n_idx_x2b
    integer :: n_wgts_a, n_wgts_b, ndim, ncells, i, count, pdim
    type(ESMF_XGridType), pointer :: xgtypep
    type(ESMF_DELayout)           :: delayout
    type(ESMF_XGridGeomType_Flag) :: xggt
    type(ESMF_CoordSys_Flag) :: coord_sys

    ! Initialize
    localrc = ESMF_RC_NOT_IMPL

    ! Initialize return code   
    if(present(rc)) rc = ESMF_RC_NOT_IMPL

    ! check init status of input XGrid
    ESMF_INIT_CHECK_DEEP(ESMF_XGridGetInit,xgrid,rc)

    xgtypep => xgrid%xgtypep

    if(present(sideAGridCount)) then
        count = 0
        do i = 1, size(xgtypep%sideA, 1)
          call ESMF_XGridGeomBaseGet(xgtypep%sideA(i), geomtype=xggt, rc=localrc)
          if (ESMF_LogFoundError(localrc, &
              ESMF_ERR_PASSTHRU, &
              ESMF_CONTEXT, rcToReturn=rc)) return
          if(xggt == ESMF_XGRIDGEOMTYPE_GRID) count = count + 1
        enddo
        sideAGridCount = count
    endif
    if(present(sideBGridCount)) then
        count = 0
        do i = 1, size(xgtypep%sideB, 1)
          call ESMF_XGridGeomBaseGet(xgtypep%sideB(i), geomtype=xggt, rc=localrc)
          if (ESMF_LogFoundError(localrc, &
              ESMF_ERR_PASSTHRU, &
              ESMF_CONTEXT, rcToReturn=rc)) return
          if(xggt == ESMF_XGRIDGEOMTYPE_GRID) count = count + 1
        enddo
        sideBGridCount = count
    endif
    if(present(sideAMeshCount)) then
        count = 0
        do i = 1, size(xgtypep%sideA, 1)
          call ESMF_XGridGeomBaseGet(xgtypep%sideA(i), geomtype=xggt, rc=localrc)
          if (ESMF_LogFoundError(localrc, &
              ESMF_ERR_PASSTHRU, &
              ESMF_CONTEXT, rcToReturn=rc)) return
          if(xggt == ESMF_XGRIDGEOMTYPE_MESH) count = count + 1
        enddo
        sideAMeshCount = count
    endif
    if(present(sideBMeshCount)) then
        count = 0
        do i = 1, size(xgtypep%sideB, 1)
          call ESMF_XGridGeomBaseGet(xgtypep%sideB(i), geomtype=xggt, rc=localrc)
          if (ESMF_LogFoundError(localrc, &
              ESMF_ERR_PASSTHRU, &
              ESMF_CONTEXT, rcToReturn=rc)) return
          if(xggt == ESMF_XGRIDGEOMTYPE_MESH) count = count + 1
        enddo
        sideBMeshCount = count
    endif

    if(present(sideAGrid)) then
        ngrid_a = size(sideAGrid, 1)
        count = 0
        do i = 1, size(xgtypep%sideA, 1)
          call ESMF_XGridGeomBaseGet(xgtypep%sideA(i), geomtype=xggt, rc=localrc)
          if (ESMF_LogFoundError(localrc, &
              ESMF_ERR_PASSTHRU, &
              ESMF_CONTEXT, rcToReturn=rc)) return
          if(xggt == ESMF_XGRIDGEOMTYPE_GRID) count = count + 1
        enddo
        if(ngrid_a /= count) then
            call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, & 
               msg="- size of sideAGrid doesn't match the number Grids on sideA in the XGrid", &
               ESMF_CONTEXT, rcToReturn=rc) 
            return
        endif
        count = 0
        do i = 1, size(xgtypep%sideA, 1)
          call ESMF_XGridGeomBaseGet(xgtypep%sideA(i), geomtype=xggt, rc=localrc)
          if (ESMF_LogFoundError(localrc, &
              ESMF_ERR_PASSTHRU, &
              ESMF_CONTEXT, rcToReturn=rc)) return
          if(xggt == ESMF_XGRIDGEOMTYPE_GRID) then
            count = count + 1
            sideAGrid(count) = xgtypep%sideA(i)%gbcp%grid
          endif
        enddo 
    endif
    if(present(sideAMesh)) then
        ngrid_a = size(sideAMesh, 1)
        count = 0
        do i = 1, size(xgtypep%sideA, 1)
          call ESMF_XGridGeomBaseGet(xgtypep%sideA(i), geomtype=xggt, rc=localrc)
          if (ESMF_LogFoundError(localrc, &
              ESMF_ERR_PASSTHRU, &
              ESMF_CONTEXT, rcToReturn=rc)) return
          if(xggt == ESMF_XGRIDGEOMTYPE_MESH) count = count + 1
        enddo
        if(ngrid_a /= count) then
            call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, & 
               msg="- size of sideAMesh doesn't match the number Meshes on sideA in the XGrid", &
               ESMF_CONTEXT, rcToReturn=rc) 
            return
        endif
        count = 0
        do i = 1, size(xgtypep%sideA, 1)
          call ESMF_XGridGeomBaseGet(xgtypep%sideA(i), geomtype=xggt, rc=localrc)
          if (ESMF_LogFoundError(localrc, &
              ESMF_ERR_PASSTHRU, &
              ESMF_CONTEXT, rcToReturn=rc)) return
          if(xggt == ESMF_XGRIDGEOMTYPE_MESH) then
            count = count + 1
            sideAMesh(count) = xgtypep%sideA(i)%gbcp%mesh
          endif
        enddo 
    endif

    if(present(sideBGrid)) then
        ngrid_b = size(sideBGrid, 1)
        count = 0
        do i = 1, size(xgtypep%sideB, 1)
          call ESMF_XGridGeomBaseGet(xgtypep%sideB(i), geomtype=xggt, rc=localrc)
          if (ESMF_LogFoundError(localrc, &
              ESMF_ERR_PASSTHRU, &
              ESMF_CONTEXT, rcToReturn=rc)) return
          if(xggt == ESMF_XGRIDGEOMTYPE_GRID) count = count + 1
        enddo
        if(ngrid_b /= count) then
            call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, & 
               msg="- size of sideBGrid doesn't match the number Grids on sideB in the XGrid", &
               ESMF_CONTEXT, rcToReturn=rc) 
            return
        endif
        count = 0
        do i = 1, size(xgtypep%sideB, 1)
          call ESMF_XGridGeomBaseGet(xgtypep%sideB(i), geomtype=xggt, rc=localrc)
          if (ESMF_LogFoundError(localrc, &
              ESMF_ERR_PASSTHRU, &
              ESMF_CONTEXT, rcToReturn=rc)) return
          if(xggt == ESMF_XGRIDGEOMTYPE_GRID) then
            count = count + 1
            sideBGrid(count) = xgtypep%sideB(i)%gbcp%grid
          endif
        enddo 
    endif
    if(present(sideBMesh)) then
        ngrid_b = size(sideBMesh, 1)
        count = 0
        do i = 1, size(xgtypep%sideB, 1)
          call ESMF_XGridGeomBaseGet(xgtypep%sideB(i), geomtype=xggt, rc=localrc)
          if (ESMF_LogFoundError(localrc, &
              ESMF_ERR_PASSTHRU, &
              ESMF_CONTEXT, rcToReturn=rc)) return
          if(xggt == ESMF_XGRIDGEOMTYPE_MESH) count = count + 1
        enddo
        if(ngrid_b /= count) then
            call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, & 
               msg="- size of sideBMesh doesn't match the number Meshes on sideB in the XGrid", &
               ESMF_CONTEXT, rcToReturn=rc) 
            return
        endif
        count = 0
        do i = 1, size(xgtypep%sideB, 1)
          call ESMF_XGridGeomBaseGet(xgtypep%sideB(i), geomtype=xggt, rc=localrc)
          if (ESMF_LogFoundError(localrc, &
              ESMF_ERR_PASSTHRU, &
              ESMF_CONTEXT, rcToReturn=rc)) return
          if(xggt == ESMF_XGRIDGEOMTYPE_MESH) then
            count = count + 1
            sideBMesh(count) = xgtypep%sideB(i)%gbcp%mesh
          endif
        enddo 
    endif

    if(present(mesh)) then
      if(.not. xgtypep%storeOverlay) then
          call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, &
             msg="- Cannot retrieve super mesh when storeOverylay is false.", &
             ESMF_CONTEXT, rcToReturn=rc)
          return
      endif    
      mesh = xgtypep%mesh
    endif

    if(present(area)) then
        ncells = size(area,1)
        if(.not. associated(xgtypep%area)) then
            call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, &
               msg="- uninitialized area in the XGrid", &
               ESMF_CONTEXT, rcToReturn=rc)
            return
        endif    
        if(ncells /= size(xgtypep%area, 1)) then
            call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, &
               msg=" size of area array doesn't match the size of area in the XGrid", &
               ESMF_CONTEXT, rcToReturn=rc)
            return
        endif
        area = xgtypep%area
    endif 
    if(present(centroid)) then
        ndim = size(centroid, 1)
        ncells = size(centroid, 2)
        if(.not. associated(xgtypep%centroid)) then
            call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, &
               msg="- uninitialized centroid in the XGrid", &
               ESMF_CONTEXT, rcToReturn=rc)
            return
        endif    
        if(ncells /= size(xgtypep%centroid, 2) .or. &
           ndim  /= size(xgtypep%centroid, 1)) then
            call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, &
               msg="- size of centroid doesn't match the size of centroid in the XGrid", &
               ESMF_CONTEXT, rcToReturn=rc)
            return
        endif
        centroid = xgtypep%centroid
    endif 

    if(present(distgridA)) then
        ngrid_a = size(distgridA)
        if(ngrid_a /= size(xgtypep%distgridA, 1)) then
            call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, & 
               msg="- size of distgridA doesn't match the size of distgridA in the XGrid", &
               ESMF_CONTEXT, rcToReturn=rc) 
            return
        endif
        do i = 1, ngrid_a
            distgridA(i) = xgtypep%distgridA(i)
        enddo 
    endif

    if(present(distgridB)) then
        ngrid_b = size(distgridB)
        if(ngrid_b /= size(xgtypep%distgridB, 1)) then
            call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, & 
               msg="- size of distgridB doesn't match the size of distgridB in the XGrid", &
               ESMF_CONTEXT, rcToReturn=rc) 
            return
        endif
        do i = 1, ngrid_b
            distgridB(i) = xgtypep%distgridB(i)
        enddo 
    endif

    if(present(distgridM)) then
        distgridM = xgtypep%distgridM
    endif

    if(present(sparseMatA2X)) then
        ngrid_a = size(sparseMatA2X, 1)
        if(.not. associated(xgtypep%sparseMatA2X)) then
            call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, &
               msg="- uninitialized sparseMatA2X in the XGrid", &
               ESMF_CONTEXT, rcToReturn=rc)
            return
        endif    
        if(ngrid_a /= size(xgtypep%sparseMatA2X, 1)) then
            call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, & 
               msg="- size of sparseMatA2X doesn't match the size of sparseMatA2X in the XGrid", &
               ESMF_CONTEXT, rcToReturn=rc) 
            return
        endif
        do i = 1, ngrid_a
            sparseMatA2X(i) = xgtypep%sparseMatA2X(i)
        enddo 
    endif

    if(present(sparseMatX2A)) then
        ngrid_a = size(sparseMatX2A, 1)
        if(.not. associated(xgtypep%sparseMatX2A)) then
            call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, &
               msg="- uninitialized sparseMatX2A in the XGrid", &
               ESMF_CONTEXT, rcToReturn=rc)
            return
        endif    
        if(ngrid_a /= size(xgtypep%sparseMatX2A, 1)) then
            call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, & 
               msg="- size of sparseMatX2A doesn't match the size of sparseMatX2A in the XGrid", &
               ESMF_CONTEXT, rcToReturn=rc) 
            return
        endif
        do i = 1, ngrid_a
            sparseMatX2A(i) = xgtypep%sparseMatX2A(i)
        enddo 
    endif

    if(present(sparseMatB2X)) then
        ngrid_a = size(sparseMatB2X, 1)
        if(.not. associated(xgtypep%sparseMatB2X)) then
            call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, &
               msg="- uninitialized sparseMatB2X in the XGrid", &
               ESMF_CONTEXT, rcToReturn=rc)
            return
        endif    
        if(ngrid_a /= size(xgtypep%sparseMatB2X, 1)) then
            call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, & 
               msg="- size of sparseMatB2X doesn't match the size of sparseMatB2X in the XGrid", &
               ESMF_CONTEXT, rcToReturn=rc) 
            return
        endif
        do i = 1, ngrid_a
            sparseMatB2X(i) = xgtypep%sparseMatB2X(i)
        enddo 
    endif

    if(present(sparseMatX2B)) then
        ngrid_a = size(sparseMatX2B, 1)
        if(.not. associated(xgtypep%sparseMatX2B)) then
            call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, &
               msg="- uninitialized sparseMatX2B in the XGrid", &
               ESMF_CONTEXT, rcToReturn=rc)
            return
        endif    
        if(ngrid_a /= size(xgtypep%sparseMatX2B, 1)) then
            call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, & 
               msg="- size of sparseMatX2B doesn't match the size of sparseMatX2B in the XGrid", &
               ESMF_CONTEXT, rcToReturn=rc) 
            return
        endif
        do i = 1, ngrid_a
            sparseMatX2B(i) = xgtypep%sparseMatX2B(i)
        enddo 
    endif
    
    if (present(name)) then
        call ESMF_GetName(xgtypep%base, name, localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
    endif

    ! Get coordSys
    if(present(coordSys)) then
       coordSys=xgtypep%coordSys
    endif

    if(present(dimCount)) then
      if(.not. xgtypep%storeOverlay) then
          call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, &
             msg="- Cannot retrieve dimCount of the super mesh when storeOverylay is false.", &
             ESMF_CONTEXT, rcToReturn=rc)
          return
      endif
      call C_ESMC_MeshGetDimensions(xgtypep%mesh%this, dimCount, pdim, coord_sys, localrc);
      if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
        ESMF_CONTEXT, rcToReturn=rc)) return
    endif
    if(present(elementCount)) then
        call ESMF_DistGridGet(xgtypep%distgridM, localDe=0, elementCount=elementCount, &
            rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
    endif

    ! success
    if(present(rc)) rc = ESMF_SUCCESS

end subroutine ESMF_XGridGetDefault