ESMF_FieldRegridStoreX Subroutine

private subroutine ESMF_FieldRegridStoreX(xgrid, srcField, dstField, keywordEnforcer, regridmethod, srcTermProcessing, pipeLineDepth, routehandle, srcFracField, dstFracField, srcMergeFracField, dstMergeFracField, rc)

Arguments

Type IntentOptional Attributes Name
type(ESMF_XGrid), intent(in) :: xgrid
type(ESMF_Field), intent(in) :: srcField
type(ESMF_Field), intent(inout) :: dstField
type(ESMF_KeywordEnforcer), optional :: keywordEnforcer
type(ESMF_RegridMethod_Flag), intent(in), optional :: regridmethod
integer, intent(inout), optional :: srcTermProcessing
integer, intent(inout), optional :: pipeLineDepth
type(ESMF_RouteHandle), intent(inout), optional :: routehandle
type(ESMF_Field), intent(inout), optional :: srcFracField
type(ESMF_Field), intent(inout), optional :: dstFracField
type(ESMF_Field), intent(inout), optional :: srcMergeFracField
type(ESMF_Field), intent(inout), optional :: dstMergeFracField
integer, intent(out), optional :: rc

Source Code

    subroutine ESMF_FieldRegridStoreX(xgrid, srcField, dstField, keywordEnforcer, &
                    regridmethod, &
                    srcTermProcessing, pipeLineDepth, &
                    routehandle, &
                    srcFracField, dstFracField, &
                    srcMergeFracField, dstMergeFracField, rc)
!      
! !ARGUMENTS:
      type(ESMF_XGrid),       intent(in)                      :: xgrid
      type(ESMF_Field),       intent(in)                      :: srcField
      type(ESMF_Field),       intent(inout)                   :: dstField
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
      type(ESMF_RegridMethod_Flag),   intent(in),    optional :: regridmethod
      integer,                        intent(inout), optional :: srcTermProcessing
      integer,                        intent(inout), optional :: pipeLineDepth
      type(ESMF_RouteHandle), intent(inout), optional         :: routehandle
      type(ESMF_Field),       intent(inout), optional         :: srcFracField
      type(ESMF_Field),       intent(inout), optional         :: dstFracField
      type(ESMF_Field),       intent(inout), optional         :: srcMergeFracField
      type(ESMF_Field),       intent(inout), optional         :: dstMergeFracField
      integer,                intent(out),   optional         :: rc 
!
! !STATUS:
! \begin{itemize}
! \item\apiStatusCompatibleVersion{5.2.0r}
! \item\apiStatusModifiedSinceVersion{5.2.0r}
! \begin{description}
! \item[5.3.0] Added arguments {\tt srcFracField}, {\tt dstFracField}, {\tt srcMergeFracField}, and {\tt dstMergeFracField}.
! These fraction Fields allow a user to calculate correct flux regridded through {\tt ESMF\_XGrid}.
! \item[7.1.0r] Added argument {\tt regridmethod}. This new argument allows the user to choose the regrid method
!               to apply when computing the routehandle. 
! \item[8.5.0] Added arguments {\tt srcTermProcessing} and {\tt pipelineDepth} to
!              provide access to the tuning parameters affecting the sparse matrix
!              execution. See the text for details on the impact
!              {\tt srcTermProcessing} can have on bit-for-bit reproducibility.
! \end{description}
! \end{itemize}
!
! !DESCRIPTION:
!       \begin{sloppypar}
!       This method creates a RouteHandle to do conservative interpolation specifically between a
!       Field built on an XGrid and a Field build on one of the Grids or Meshes used to create that XGrid. 
!       (To do more general interpolation use the {\tt ESMF\_FieldRegridStore()} method
!       in section~\ref{api:esmf_fieldregridstorenx}.) The RouteHandle produced by this method can then be used in the call
!       {\tt ESMF\_FieldRegrid()} to interpolate from the {\tt srcField} to the {\tt dstField}. 
!       \end{sloppypar}
!       
!       The RouteHandle generated by this call is based just on the 
!       coordinates in the Grids, XGrids, or Meshes contained in the Fields.  If those
!       coordinates don't change the RouteHandle can
!       be used repeatedly to interpolate from the source Field to the 
!       destination Field.  This is true even if the data in the Fields 
!       changes. The RouteHandle may also be used to interpolate between any
!       source and destination Field which are created 
!       on the same Grid, XGrid, or Mesh as the original Fields.        
!
!       When it's no longer needed the RouteHandle should be destroyed by 
!       using {\tt ESMF\_FieldRegridRelease()} to free the memory it's using. 
!
!       Note, as a side effect, that this call may change the data in {\tt dstField}. If
!       this is undesirable, then an easy work around is to create a second temporary Field
!       with the same structure as {\tt dstField} and pass that in instead. 
!
!     The arguments are:
!     \begin{description}
!     \item [xgrid]
!           Exchange Grid.
!     \item [srcField]
!           Source Field built on either {\tt xgrid} or one of the Grids or Meshes used to create {\tt xgrid}.
!     \item [dstField]
!           Destination Field built on either {\tt xgrid} or one of the Grids or Meshes used to create {\tt xgrid}. 
!           The data in this Field may be overwritten by this call. 
!     \item [{[regridmethod]}]
!           The type of interpolation. For this method only 
!           {\tt ESMF\_REGRIDMETHOD\_CONSERVE} and {\tt ESMF\_REGRIDMETHOD\_CONSERVE\_2ND} are
!           supported. If not specified, defaults to {\tt ESMF\_REGRIDMETHOD\_CONSERVE}.
!     \item [{[srcTermProcessing]}]
!           The {\tt srcTermProcessing} parameter controls how many source terms,
!           located on the same PET and summing into the same destination element,
!           are summed into partial sums on the source PET before being transferred
!           to the destination PET. A value of 0 indicates that the entire arithmetic
!           is done on the destination PET; source elements are neither multiplied 
!           by their factors nor added into partial sums before being sent off by the
!           source PET. A value of 1 indicates that source elements are multiplied
!           by their factors on the source side before being sent to the destination
!           PET. Larger values of {\tt srcTermProcessing} indicate the maximum number
!           of terms in the partial sums on the source side.
!
!     Note that partial sums may lead to bit-for-bit differences in the results.
!     See section \ref{RH:bfb} for an in-depth discussion of {\em all}
!     bit-for-bit reproducibility aspects related to route-based communication
!     methods.
!
!     \begin{sloppypar}
!     The {\tt ESMF\_FieldRegridStore()} method implements an auto-tuning scheme
!     for the {\tt srcTermProcessing} parameter. The intent on the 
!     {\tt srcTermProcessing} argument is "{\tt inout}" in order to 
!     support both overriding and accessing the auto-tuning parameter.
!     If an argument $>= 0$ is specified, it is used for the 
!     {\tt srcTermProcessing} parameter, and the auto-tuning phase is skipped.
!     In this case the {\tt srcTermProcessing} argument is not modified on
!     return. If the provided argument is $< 0$, the {\tt srcTermProcessing}
!     parameter is determined internally using the auto-tuning scheme. In this
!     case the {\tt srcTermProcessing} argument is re-set to the internally
!     determined value on return. Auto-tuning is also used if the optional 
!     {\tt srcTermProcessing} argument is omitted.
!     \end{sloppypar}
!     
!   \item [{[pipelineDepth]}]
!     The {\tt pipelineDepth} parameter controls how many messages a PET
!     may have outstanding during a sparse matrix exchange. Larger values
!     of {\tt pipelineDepth} typically lead to better performance. However,
!     on some systems too large a value may lead to performance degradation,
!     or runtime errors.
!
!     Note that the pipeline depth has no effect on the bit-for-bit
!     reproducibility of the results. However, it may affect the performance
!     reproducibility of the exchange.
!
!     The {\tt ESMF\_FieldRegridStore()} method implements an auto-tuning scheme
!     for the {\tt pipelineDepth} parameter. The intent on the 
!     {\tt pipelineDepth} argument is "{\tt inout}" in order to 
!     support both overriding and accessing the auto-tuning parameter.
!     If an argument $>= 0$ is specified, it is used for the 
!     {\tt pipelineDepth} parameter, and the auto-tuning phase is skipped.
!     In this case the {\tt pipelineDepth} argument is not modified on
!     return. If the provided argument is $< 0$, the {\tt pipelineDepth}
!     parameter is determined internally using the auto-tuning scheme. In this
!     case the {\tt pipelineDepth} argument is re-set to the internally
!     determined value on return. Auto-tuning is also used if the optional 
!     {\tt pipelineDepth} argument is omitted.
!     \item [{[routehandle]}]
!           The handle that implements the regrid and that can be used in later 
!           {\tt ESMF\_FieldRegrid}.
!     \item [{[srcFracField]}] 
!           The fraction of each source cell participating in the regridding returned from this call. 
!           This Field needs to be created on the same Grid and location (e.g staggerloc) 
!           as the srcField.
!     \item [{[dstFracField]}] 
!           The fraction of each destination cell participating in the regridding returned from this call. 
!           This Field needs to be created on the same Grid and location (e.g staggerloc) 
!           as the dstField.
!     \item [{[srcMergeFracField]}] 
!           The fraction of each source cell as a result of Grid merge returned from this call.
!           This Field needs to be created on the same Grid and location (e.g staggerloc) 
!           as the srcField.
!     \item [{[dstMergeFracField]}] 
!           The fraction of each destination cell as a result of Grid merge returned from this call.
!           This Field needs to be created on the same Grid and location (e.g staggerloc) 
!           as the dstField.
!     \item [{[rc]}]
!           Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!
!EOP
        integer :: localrc, i

        type(ESMF_GeomType_Flag)  :: geomtype, srcgeomtype, dstgeomtype
        type(ESMF_XGrid)     :: srcXGrid, dstXGrid        
        type(ESMF_Mesh)      :: srcMesh, dstMesh

        integer :: srcIdx, dstIdx, ngrid_a, ngrid_b
        integer :: sideAGC, sideAMC, sideBGC, sideBMC
        type(ESMF_XGridSide_Flag) :: srcSide, dstSide
        type(ESMF_XGridGeomBase), allocatable :: gridA(:), gridB(:)
        type(ESMF_Grid)      :: srcGrid
        type(ESMF_Grid)      :: dstGrid
        type(ESMF_XGridSpec) :: sparseMat
        logical :: found, match
        type(ESMF_Array)     :: srcFracArray
        type(ESMF_Array)     :: dstFracArray
        type(ESMF_STAGGERLOC):: interpFieldStaggerloc, fracFieldStaggerloc
        type(ESMF_MESHLOC)   :: interpFieldMeshloc, fracFieldMeshloc
        type(ESMF_RegridMethod_Flag) :: lregridmethod
        type(ESMF_Mesh)      :: superMesh
        type(ESMF_Field)     :: tmpSrcField, tmpDstField
        type(ESMF_Typekind_Flag) :: fieldTypeKind

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

        ! Set optional method argument
        if (present(regridmethod)) then
           lregridmethod=regridmethod
        else     
           lregridmethod=ESMF_REGRIDMETHOD_CONSERVE
        endif

        ! Only conservative methods supported for now
        if ((lregridmethod .ne. ESMF_REGRIDMETHOD_CONSERVE) .and. &
             (lregridmethod .ne. ESMF_REGRIDMETHOD_CONSERVE_2ND)) then
           call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, & 
                msg="- Only conservative regrid methods supported through XGrid", & 
                ESMF_CONTEXT, rcToReturn=rc) 
           return
        endif

        ! look for the correct Grid to use
        ! first Get necessary information from XGrid and Fields
        call ESMF_XGridGet(xgrid, &
            sideAGridCount=sideAGC, sideAMeshCount=sideAMC, &
            sideBGridCount=sideBGC, sideBMeshCount=sideBMC, &
            rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
        ngrid_a = sideAGC + sideAMC
        ngrid_b = sideBGC + sideBMC
        allocate(gridA(ngrid_a), gridB(ngrid_b))

        call ESMF_XGridGet(xgrid, gridA, gridB, rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return

        call ESMF_FieldGet(srcField, geomtype=geomtype, &
               rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
        srcgeomtype = geomtype

        ! locate the Grid or XGrid contained in srcField
        if(geomtype == ESMF_GEOMTYPE_GRID) then
            call ESMF_FieldGet(srcField, grid=srcGrid, &
                   rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
              ESMF_CONTEXT, rcToReturn=rc)) return

            found = .false.
            do i = 1, ngrid_a
                !if(ESMF_GridMatch(srcGrid, gridA(i)%gbcp%grid, &
                !     globalflag=.true.) >=ESMF_GRIDMATCH_EXACT) then
                if(srcGrid == gridA(i)%gbcp%grid) then
                    srcIdx = i
                    srcSide = ESMF_XGRIDSIDE_A
                    found = .true.
                    exit
                endif
            enddo 
            do i = 1, ngrid_b
                !if(ESMF_GridMatch(srcGrid, gridB(i)%gbcp%grid, &
                !     globalflag=.true.) >=ESMF_GRIDMATCH_EXACT) then
                if(srcGrid == gridB(i)%gbcp%grid) then
                    if(found) then
                      ! TODO: maybe we should attach standard attibute
                      ! to differentiate src and dst side for regridding
                      call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, & 
                        msg="- duplication of Grid found in XGrid", &
                        ESMF_CONTEXT, rcToReturn=rc) 
                      return
                    endif
                    srcIdx = i
                    srcSide = ESMF_XGRIDSIDE_B
                    found = .true.
                    exit
                endif
            enddo 

            if(.not. found) then
                call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, & 
                   msg="- cannot Locate src Field Grid in XGrid", &
                   ESMF_CONTEXT, rcToReturn=rc) 
                return
            endif

        else if(geomtype == ESMF_GEOMTYPE_MESH) then
            call ESMF_FieldGet(srcField, mesh=srcMesh, &
                   rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
              ESMF_CONTEXT, rcToReturn=rc)) return

            found = .false.
            do i = 1, ngrid_a
                !if(ESMF_MeshMatch(srcMesh, gridA(i)%gbcp%mesh)) then
                if(srcMesh == gridA(i)%gbcp%mesh) then
                    srcIdx = i
                    srcSide = ESMF_XGRIDSIDE_A
                    found = .true.
                    exit
                endif
            enddo 
            do i = 1, ngrid_b
                !if(ESMF_MeshMatch(srcMesh, gridB(i)%gbcp%mesh)) then
                if(srcMesh == gridB(i)%gbcp%mesh) then
                    if(found) then
                      call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, & 
                        msg="- duplication of Mesh found in XGrid", &
                        ESMF_CONTEXT, rcToReturn=rc) 
                      return
                    endif
                    srcIdx = i
                    srcSide = ESMF_XGRIDSIDE_B
                    found = .true.
                    exit
                endif
            enddo 

            if(.not. found) then
                call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, & 
                   msg="- cannot Locate src Field Mesh in XGrid", &
                   ESMF_CONTEXT, rcToReturn=rc) 
                return
            endif
                
        else if(geomtype == ESMF_GEOMTYPE_XGRID) then
            call ESMF_FieldGet(srcField, xgrid=srcXGrid, &
                   rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
              ESMF_CONTEXT, rcToReturn=rc)) return
            
            match = ESMF_XGridMatch(xgrid, srcXGrid, rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
              ESMF_CONTEXT, rcToReturn=rc)) return

            if(match) then
                srcSide = ESMF_XGRIDSIDE_BALANCED
                srcIdx = 1
            else
                call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, &
                  msg="- XGrid in srcField doesn't match the input XGrid", &
                  ESMF_CONTEXT, rcToReturn=rc) 
                return
            endif
        else
            call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, &
              msg="- src Field is not built on Grid, Mesh, or XGrid", &
              ESMF_CONTEXT, rcToReturn=rc) 
            return
        endif

        ! locate the Grid or XGrid contained in dstField
        call ESMF_FieldGet(dstField, geomtype=geomtype, &
               rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
        dstgeomtype = geomtype

        if(geomtype == ESMF_GEOMTYPE_GRID) then
            call ESMF_FieldGet(dstField, grid=dstGrid, &
                   rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
              ESMF_CONTEXT, rcToReturn=rc)) return

            found = .false.
            do i = 1, ngrid_a
                !if(ESMF_GridMatch(dstGrid, gridA(i)%gbcp%grid, &
                !   globalflag=.true.) >=ESMF_GRIDMATCH_EXACT) then
                if(dstGrid == gridA(i)%gbcp%grid) then
                    dstIdx = i
                    dstSide = ESMF_XGRIDSIDE_A
                    found = .true.
                    exit
                endif
            enddo 
            do i = 1, ngrid_b
                !if(ESMF_GridMatch(dstGrid, gridB(i)%gbcp%grid, &
                !   globalflag=.true.) >=ESMF_GRIDMATCH_EXACT) then
                if(dstGrid == gridB(i)%gbcp%grid) then
                    if(found) then
                      call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, & 
                        msg="- duplication of Grid found in XGrid", &
                        ESMF_CONTEXT, rcToReturn=rc) 
                      return
                    endif
                    dstIdx = i
                    dstSide = ESMF_XGRIDSIDE_B
                    found = .true.
                    exit
                endif
            enddo 

            if(.not. found) then
                call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, & 
                   msg="- cannot Locate dst Field Grid in XGrid", &
                   ESMF_CONTEXT, rcToReturn=rc) 
                return
            endif

        else if(geomtype == ESMF_GEOMTYPE_MESH) then
            call ESMF_FieldGet(dstField, mesh=dstMesh, &
                   rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
              ESMF_CONTEXT, rcToReturn=rc)) return

            found = .false.
            do i = 1, ngrid_a
                !if(ESMF_MeshMatch(dstMesh, gridA(i)%gbcp%mesh)) then
                if(dstMesh == gridA(i)%gbcp%mesh) then
                    dstIdx = i
                    dstSide = ESMF_XGRIDSIDE_A
                    found = .true.
                    exit
                endif
            enddo 
            do i = 1, ngrid_b
                !if(ESMF_MeshMatch(dstMesh, gridB(i)%gbcp%mesh)) then
                if(dstMesh == gridB(i)%gbcp%mesh) then
                    if(found) then
                      call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, & 
                        msg="- duplication of Mesh found in XGrid", &
                        ESMF_CONTEXT, rcToReturn=rc) 
                      return
                    endif
                    dstIdx = i
                    dstSide = ESMF_XGRIDSIDE_B
                    found = .true.
                    exit
                endif
            enddo 

            if(.not. found) then
                call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, & 
                   msg="- cannot Locate dst Field Mesh in XGrid", &
                   ESMF_CONTEXT, rcToReturn=rc) 
                return
            endif

        else if(geomtype == ESMF_GEOMTYPE_XGRID) then
            call ESMF_FieldGet(dstField, xgrid=dstXGrid, &
                   rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
              ESMF_CONTEXT, rcToReturn=rc)) return

            match = ESMF_XGridMatch(xgrid, dstXGrid, rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
              ESMF_CONTEXT, rcToReturn=rc)) return

            if(match) then
                dstSide = ESMF_XGRIDSIDE_BALANCED
                dstIdx = 1
            else
                call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, &
                  msg="- XGrid in dstField doesn't match the input XGrid", &
                  ESMF_CONTEXT, rcToReturn=rc) 
                return
            endif
        else
            call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, &
              msg="- src Field is not built on Grid or XGrid", &
              ESMF_CONTEXT, rcToReturn=rc) 
            return
        endif

        ! src and dst Fields should not be on the same side
        if ( srcSide == dstSide ) then
            call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, & 
               msg="- src and dst Fields should not be on same side of the XGrid", &
               ESMF_CONTEXT, rcToReturn=rc) 
            return
        endif

        ! retrieve regridding fraction Fields on demand
        if(present(srcFracField)) then
          if(srcgeomtype == ESMF_GEOMTYPE_GRID) then
            call ESMF_FieldGet(srcFracField, staggerloc=fracFieldStaggerloc, &
                 array=srcFracArray, rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, & 
              ESMF_CONTEXT, rcToReturn=rc)) return

            call ESMF_FieldGet(srcField, staggerloc=interpFieldStaggerloc, &
                 rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, & 
              ESMF_CONTEXT, rcToReturn=rc)) return
      
            ! Make sure the staggerlocs match
            if (interpFieldStaggerloc .ne. fracFieldStaggerloc) then
               call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, &
                 msg="- fracField Field staggerloc must match interpField staggerloc", &
                 ESMF_CONTEXT, rcToReturn=rc)
               return
            endif
          else if(srcgeomtype == ESMF_GEOMTYPE_MESH) then
            call ESMF_FieldGet(srcFracField, meshloc=fracFieldMeshloc, &
                 array=srcFracArray, rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, & 
              ESMF_CONTEXT, rcToReturn=rc)) return

            call ESMF_FieldGet(srcField, meshloc=interpFieldMeshloc, &
                 rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, & 
              ESMF_CONTEXT, rcToReturn=rc)) return
      
            ! Make sure the staggerlocs match
            if (interpFieldMeshloc .ne. fracFieldMeshloc) then
               call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, &
                 msg="- fracField Field Meshloc must match interpField Meshloc", &
                 ESMF_CONTEXT, rcToReturn=rc)
               return
            endif
          else
            call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, &
              msg="- src Field has unrecognized geom type", &
              ESMF_CONTEXT, rcToReturn=rc)
            return
          endif

          call ESMF_XGridGet(xgrid, srcSide, srcIdx, &
              dstSide, dstIdx, srcFracArray=srcFracArray, rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
        endif

        if(present(dstFracField)) then
          if(dstgeomtype == ESMF_GEOMTYPE_GRID) then
            call ESMF_FieldGet(dstFracField, staggerloc=fracFieldStaggerloc, &
                 array=dstFracArray, rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, & 
              ESMF_CONTEXT, rcToReturn=rc)) return

            call ESMF_FieldGet(dstField, staggerloc=interpFieldStaggerloc, &
                 rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, & 
              ESMF_CONTEXT, rcToReturn=rc)) return
      
            ! Make sure the staggerlocs match
            if (interpFieldStaggerloc .ne. fracFieldStaggerloc) then
               call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, &
                 msg="- fracField Field staggerloc must match interpField staggerloc", &
                 ESMF_CONTEXT, rcToReturn=rc)
               return
            endif
          else if(dstgeomtype == ESMF_GEOMTYPE_MESH) then
            call ESMF_FieldGet(dstFracField, meshloc=fracFieldMeshloc, &
                 array=dstFracArray, rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, & 
              ESMF_CONTEXT, rcToReturn=rc)) return

            call ESMF_FieldGet(dstField, meshloc=interpFieldMeshloc, &
                 rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, & 
              ESMF_CONTEXT, rcToReturn=rc)) return
      
            ! Make sure the staggerlocs match
            if (interpFieldMeshloc .ne. fracFieldMeshloc) then
               call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, &
                 msg="- fracField Field Meshloc must match interpField Meshloc", &
                 ESMF_CONTEXT, rcToReturn=rc)
               return
            endif
          else
            call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, &
              msg="- dst Field has unrecognized geom type", &
              ESMF_CONTEXT, rcToReturn=rc)
            return
          endif

          call ESMF_XGridGet(xgrid, srcSide, srcIdx, &
              dstSide, dstIdx, dstFracArray=dstFracArray, rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
        endif

        ! retrieve regridding fraction2 Fields on demand
        if(present(srcMergeFracField)) then
          if(srcgeomtype == ESMF_GEOMTYPE_GRID) then
            call ESMF_FieldGet(srcMergeFracField, staggerloc=fracFieldStaggerloc, &
                 array=srcFracArray, rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, & 
              ESMF_CONTEXT, rcToReturn=rc)) return

            call ESMF_FieldGet(srcField, staggerloc=interpFieldStaggerloc, &
                 rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, & 
              ESMF_CONTEXT, rcToReturn=rc)) return
      
            ! Make sure the staggerlocs match
            if (interpFieldStaggerloc .ne. fracFieldStaggerloc) then
               call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, &
                 msg="- fracField Field staggerloc must match interpField staggerloc", &
                 ESMF_CONTEXT, rcToReturn=rc)
               return
            endif
          else if(srcgeomtype == ESMF_GEOMTYPE_MESH) then
            call ESMF_FieldGet(srcFracField, meshloc=fracFieldMeshloc, &
                 array=srcFracArray, rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, & 
              ESMF_CONTEXT, rcToReturn=rc)) return

            call ESMF_FieldGet(srcField, meshloc=interpFieldMeshloc, &
                 rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, & 
              ESMF_CONTEXT, rcToReturn=rc)) return
      
            ! Make sure the staggerlocs match
            if (interpFieldMeshloc .ne. fracFieldMeshloc) then
               call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, &
                 msg="- fracField Field Meshloc must match interpField Meshloc", &
                 ESMF_CONTEXT, rcToReturn=rc)
               return
            endif
          else
            call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, &
              msg="- src Field has unrecognized geom type", &
              ESMF_CONTEXT, rcToReturn=rc)
            return
          endif

          call ESMF_XGridGet(xgrid, srcSide, srcIdx, &
              dstSide, dstIdx, srcFrac2Array=srcFracArray, rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
        endif

        if(present(dstMergeFracField)) then
          if(dstgeomtype == ESMF_GEOMTYPE_GRID) then
            call ESMF_FieldGet(dstMergeFracField, staggerloc=fracFieldStaggerloc, &
                 array=dstFracArray, rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, & 
              ESMF_CONTEXT, rcToReturn=rc)) return

            call ESMF_FieldGet(dstField, staggerloc=interpFieldStaggerloc, &
                 rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, & 
              ESMF_CONTEXT, rcToReturn=rc)) return
      
            ! Make sure the staggerlocs match
            if (interpFieldStaggerloc .ne. fracFieldStaggerloc) then
               call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, &
                 msg="- fracField Field staggerloc must match interpField staggerloc", &
                 ESMF_CONTEXT, rcToReturn=rc)
               return
            endif
          else if(dstgeomtype == ESMF_GEOMTYPE_MESH) then
            call ESMF_FieldGet(dstFracField, meshloc=fracFieldMeshloc, &
                 array=dstFracArray, rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, & 
              ESMF_CONTEXT, rcToReturn=rc)) return

            call ESMF_FieldGet(dstField, meshloc=interpFieldMeshloc, &
                 rc=localrc)
            if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, & 
              ESMF_CONTEXT, rcToReturn=rc)) return
      
            ! Make sure the staggerlocs match
            if (interpFieldMeshloc .ne. fracFieldMeshloc) then
               call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, &
                 msg="- fracField Field Meshloc must match interpField Meshloc", &
                 ESMF_CONTEXT, rcToReturn=rc)
               return
            endif
          else
            call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, &
              msg="- dst Field has unrecognized geom type", &
              ESMF_CONTEXT, rcToReturn=rc)
            return
          endif

          call ESMF_XGridGet(xgrid, srcSide, srcIdx, &
              dstSide, dstIdx, dstFrac2Array=dstFracArray, rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
        endif

        ! Create routehandle based on regrid method
        if (lregridmethod .eq. ESMF_REGRIDMETHOD_CONSERVE) then

           ! retrieve the correct sparseMat structure
           call ESMF_XGridGet(xgrid, srcSide, srcIdx, &
                dstSide, dstIdx, sparseMat=sparseMat, rc=localrc)
           if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                ESMF_CONTEXT, rcToReturn=rc)) return
           
           ! call FieldSMMStore
           call ESMF_FieldSMMStore(srcField, dstField, routehandle, &
                sparseMat%factorList, sparseMat%factorIndexList, &
                srcTermProcessing=srcTermProcessing, &
                pipeLineDepth=pipeLineDepth, &
                rc=localrc)
           if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                ESMF_CONTEXT, rcToReturn=rc)) return
        else

           ! Set temporary field for source
           if (srcSide == ESMF_XGRIDSIDE_BALANCED) then

              ! Get Field typekind
              call ESMF_FieldGet(srcField, typekind=fieldTypeKind, &
                   rc=localrc) 
              if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                   ESMF_CONTEXT, rcToReturn=rc)) return

              ! Get Super Mesh
              call ESMF_XGridGet(xgrid, mesh=superMesh, rc=localrc)
              if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                   ESMF_CONTEXT, rcToReturn=rc)) return

              ! Create temporary field
              tmpSrcField=ESMF_FieldCreate(superMesh, &
                   typekind=fieldTypeKind, &
                   meshloc=ESMF_MESHLOC_ELEMENT, &
                   rc=localrc)
              if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                   ESMF_CONTEXT, rcToReturn=rc)) return
           else 
              tmpSrcField=srcField
           endif

           ! Set temporary field for dst
           if (dstSide == ESMF_XGRIDSIDE_BALANCED) then

              ! Get Field typekind
              call ESMF_FieldGet(dstField, typekind=fieldTypeKind, &
                   rc=localrc) 
              if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                   ESMF_CONTEXT, rcToReturn=rc)) return

              ! Get Super Mesh
              call ESMF_XGridGet(xgrid, mesh=superMesh, rc=localrc)
              if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                   ESMF_CONTEXT, rcToReturn=rc)) return

              ! Create temporary field
              tmpDstField=ESMF_FieldCreate(superMesh, &
                   typekind=fieldTypeKind, &
                   meshloc=ESMF_MESHLOC_ELEMENT, &
                   rc=localrc)
              if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                   ESMF_CONTEXT, rcToReturn=rc)) return
           else 
              tmpDstField=dstField
           endif

           ! Generate routehandle other that 1st order conserve
           call ESMF_FieldRegridStoreNX(&
                srcField=tmpSrcField, &
                dstField=tmpDstField, &
! ??            srcMaskValues, dstMaskValues, &
                regridmethod=lregridmethod, &
                srcTermProcessing=srcTermProcessing, &
                pipeLineDepth=pipeLineDepth, &
                routehandle=routehandle, &
                srcFracField=srcFracField, &
                dstFracField=dstFracField, &
                rc=localrc)
           if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                ESMF_CONTEXT, rcToReturn=rc)) return

           ! Get rid of temporary source Field if necessary
           if (srcSide == ESMF_XGRIDSIDE_BALANCED) then

              call ESMF_FieldDestroy(tmpSrcField, rc=localrc)
              if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                   ESMF_CONTEXT, rcToReturn=rc)) return
           endif

           ! Get rid of temporary destination Field if necessary
           if (dstSide == ESMF_XGRIDSIDE_BALANCED) then

              call ESMF_FieldDestroy(tmpDstField, rc=localrc)
              if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                   ESMF_CONTEXT, rcToReturn=rc)) return
           endif

        endif


        if(present(rc)) rc = ESMF_SUCCESS

    end subroutine ESMF_FieldRegridStoreX