ESMF_GridToMesh Function

public function ESMF_GridToMesh(grid, staggerLoc, isSphere, isLatLonDeg, maskValues, regridConserve, rc)

Arguments

Type IntentOptional Attributes Name
type(ESMF_Grid), intent(in) :: grid
type(ESMF_StaggerLoc), intent(in) :: staggerLoc
integer, intent(in) :: isSphere
logical, intent(in), optional :: isLatLonDeg
integer(kind=ESMF_KIND_I4), optional :: maskValues(:)
type(ESMF_RegridConserve), intent(in), optional :: regridConserve
integer, intent(out), optional :: rc

Return Value type(ESMF_Mesh)


Source Code

   function ESMF_GridToMesh(grid, staggerLoc, isSphere, isLatLonDeg, maskValues, regridConserve, rc)
!
!
! !RETURN VALUE:
    type(ESMF_Mesh)                               :: ESMF_GridToMesh
!
! !ARGUMENTS:
    type(ESMF_Grid), intent(in)                :: grid
    type(ESMF_StaggerLoc),  intent(in)            :: staggerLoc
    integer,                intent(in)            :: isSphere
    logical, intent(in),   optional               :: isLatLonDeg
    type(ESMF_RegridConserve), intent(in), optional :: regridConserve
    integer(ESMF_KIND_I4), optional               :: maskValues(:)
    integer, intent(out) , optional               :: rc
!
! !DESCRIPTION:
!   Create a mesh object with the same topology as the grid.
!
!   \begin{description}
!   \item [grid]
!         The grid to copy.
!   \item [staggerLoc]
!         Stagger location on grid to build.
!   \item [isSphhere]
!         1 = a spherical grid make peridoic
!   \item [isLatLonDeg]
!         true coords of grids are lat lon in deg, default depends on isSphere
!         if isSphere=1 then default is true, else is false.
!   \item [regridConserve]
!         ESMF\_REGRID\_CONSERVE\_ON turns on the conservative regridding
!   \item [{[rc]}]
!         Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!   \end{description}
!
!EOPI
!------------------------------------------------------------------------------
    integer :: localrc
    type(ESMF_Pointer) :: theMesh
    type(ESMF_InterArray) :: maskValuesArg
    type(ESMF_Index_Flag) :: indexflag
    type(ESMF_RegridConserve) :: lregridConserve
    integer :: localIsLatLonDeg
    type(ESMF_CoordSys_Flag) :: coordSys
    integer :: dimCount

    localrc = ESMF_SUCCESS

    ! Error check Grid
    call ESMF_GridGet(grid, indexflag=indexflag, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return

    ! Handle optional conserve argument
    if (present(regridConserve)) then
       lregridConserve=regridConserve
    else
       lregridConserve=ESMF_REGRID_CONSERVE_OFF
    endif

    ! Handle optional isLatLonDeg argument
    if (present(isLatLonDeg)) then
       if (isLatLonDeg) then
          localIsLatLonDeg=1
       else
          localIsLatLonDeg=0
       endif
    else
       if (isSphere .eq. 1) then
           localIsLatLonDeg=1
       else
           localIsLatLonDeg=0
       endif
    endif


    ! convert mask values
    maskValuesArg = ESMF_InterArrayCreate(maskValues, rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return

    call c_esmc_gridtomesh(grid, staggerLoc%staggerloc, isSphere, &
      localIsLatLonDeg, theMesh, maskValuesArg, &
      lregridConserve%regridconserve, localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

    call ESMF_InterArrayDestroy(maskValuesArg, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

    ESMF_GridToMesh = ESMF_MeshCreate(theMesh)

    ! Set these here, eventually this will happen automatically internally inside ESMF_MeshCreate()
    call ESMF_GridGet(grid,              &
                      coordSys=coordSys, &
                      dimCount=dimCount, &
                      rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return

    if (present(rc)) rc = ESMF_SUCCESS

    end function ESMF_GridToMesh