ESMF_LocStreamCreateByBkgMesh Function

private function ESMF_LocStreamCreateByBkgMesh(locstream, background, keywordEnforcer, unmappedaction, name, rc)

Arguments

Type IntentOptional Attributes Name
type(ESMF_LocStream), intent(in) :: locstream
type(ESMF_Mesh), intent(in) :: background
type(ESMF_KeywordEnforcer), optional :: keywordEnforcer
type(ESMF_UnmappedAction_Flag), intent(in), optional :: unmappedaction
character(len=*), intent(in), optional :: name
integer, intent(out), optional :: rc

Return Value type(ESMF_LocStream)


Source Code

      function ESMF_LocStreamCreateByBkgMesh(locstream, &
                 background, keywordEnforcer, unmappedaction, name, rc)

!
! !RETURN VALUE:
      type(ESMF_LocStream) :: ESMF_LocStreamCreateByBkgMesh

!
! !ARGUMENTS:
      type(ESMF_LocStream),           intent(in)           :: locstream
      type(ESMF_Mesh),                intent(in)           :: background
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
      type(ESMF_UnmappedAction_Flag), intent(in), optional :: unmappedaction
      character (len=*),              intent(in), optional :: name
      integer,                        intent(out),optional :: rc
!
! !DESCRIPTION:
!
!     Create an location stream from an existing one in accordance with 
!     the distribution of the background Mesh.  The entries
!     in the new location stream are redistributed, so that they lie on the same PET
!     as the piece of Mesh which contains the coordinates of the entries. The coordinates
!     of the entries are the data in the keys named ESMF:Lon, ESMF:Lat, ESMF:Radius in the 
!     case of a spherical system and ESMF:X, ESMF:Y, ESMF:Z for cartesian. To copy data in
!     Fields or FieldBundles built on {\tt locstream} to the new one simply use {\tt ESMF\_FieldRedist()}
!     or {\tt ESMF\_FieldBundleRedist()}.
!
!     The arguments are:
!     \begin{description}
!      \item[locstream]
!          Location stream from which the new location stream is to be created
!      \item[background]
!          Background Mesh which determines the distribution of entries in the new locatiion stream.
!     \item [{[unmappedaction]}]
!           Specifies what should happen if there are destination points that
!           can't be mapped to a source cell. Please see Section~\ref{const:unmappedaction} for a 
!           list of valid options. If not specified, {\tt unmappedaction} defaults to {\tt ESMF\_UNMAPPEDACTION\_ERROR}. [NOTE: the {\tt unmappedaction=ESMF\_UNMAPPEDACTION\_IGNORE} option is currently not implemented.]
!      \item[{[name]}]
!          Name of the resulting location stream
!      \item[{[rc]}]
!          Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!
!EOP
      type(ESMF_LocStreamType), pointer :: oldLStypep, newLStypep
      type(ESMF_UnmappedAction_Flag) :: localunmappedaction
      type(ESMF_DistGrid) :: newDistGrid
      type(ESMF_TypeKind_Flag) ::keyTypeKind
      character(len=ESMF_MAXSTR)    :: keytemp, string
      integer :: keyCount,i,regrid_dims,idx,idx_cart
      integer :: localrc
      integer :: pntDim, pntCount
      real(ESMF_KIND_R8), pointer  :: pntList(:)
      real(ESMF_KIND_R8), pointer  :: pntList_cart(:)
      integer, pointer :: petList(:)
      character (len=ESMF_MAXSTR)            :: coordKeyNames
      type(ESMF_CoordSys_Flag) :: coordSysLocal, coordSys_ofBkgMesh
      logical :: three_dims

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

      ! Check Variables
      ESMF_INIT_CHECK_DEEP(ESMF_LocStreamGetInit,locstream,rc)
      ESMF_INIT_CHECK_DEEP(ESMF_MeshGetInit,background,rc)


      ! Set default vale for unmappedaction
      if (present(unmappedaction)) then
         localunmappedaction=unmappedaction
      else
         localunmappedaction=ESMF_UNMAPPEDACTION_ERROR
      endif

      ! Currently ESMF_UNMAPPEDACTION_IGNORE not implemented here
      if (localunmappedaction .eq. ESMF_UNMAPPEDACTION_IGNORE) then
        if (ESMF_LogFoundError(ESMF_RC_NOT_IMPL, &
           msg=" - ESMF_UNMAPPEDACTION_IGNORE option currently not implemented ", &
            ESMF_CONTEXT, rcToReturn=rc)) return
      endif


      ! Get old locstream internal pointer
      oldLStypep=>locstream%lstypep

      call ESMF_MeshGet(background, coordSys=coordSys_ofBkgMesh, rc=localrc)
      if (ESMF_LogFoundError(localrc, &
          ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return     

       call ESMF_LocStreamGet(locstream, coordSys=coordSysLocal, rc=localrc)
      if (ESMF_LogFoundError(localrc, &
          ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return     

      if ((coordSysLocal .eq. ESMF_COORDSYS_CART .and. &
           coordSys_ofBkgMesh .ne. ESMF_COORDSYS_CART) .or. &
          (coordSysLocal .ne. ESMF_COORDSYS_CART .and. &
           coordSys_ofBkgMesh .eq. ESMF_COORDSYS_CART)) then
        if (ESMF_LogFoundError(ESMF_RC_OBJ_BAD, &
           msg=" - coordinate systems of LocStream and Mesh are not compatible ", &
            ESMF_CONTEXT, rcToReturn=rc)) return
      endif

      if (coordSysLocal .eq. ESMF_COORDSYS_CART) then
        call ESMF_LocStreamGetKey(locstream, keyName="ESMF:Z", &
                                  isPresent=three_dims, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
           ESMF_ERR_PASSTHRU, &
           ESMF_CONTEXT, rcToReturn=rc)) return
        if (three_dims) then
          coordKeyNames = "ESMF:X,ESMF:Y,ESMF:Z"
        else
          coordKeyNames = "ESMF:X,ESMF:Y"
        endif

      else
        call ESMF_LocStreamGetKey(locstream, keyName="ESMF:Radius", &
                                  isPresent=three_dims, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
           ESMF_ERR_PASSTHRU, &
           ESMF_CONTEXT, rcToReturn=rc)) return
        if (three_dims) then
          coordKeyNames = "ESMF:Lon,ESMF:Lat,ESMF:Radius"
        else
          coordKeyNames = "ESMF:Lon,ESMF:Lat"
        endif

      endif

      ! Calculate pntDim 
       pntDim = 0
       string = trim(coordKeyNames )
       do while ( string /= '' )
          call ESMF_StripKey( string, keytemp )
          pntDim = pntDim + 1
       enddo

      ! Calculate number of local points
      call ESMF_LocStreamGetNumLocal(locstream, localCount=pntCount, &
                                     rc=localrc)      
      if (ESMF_LogFoundError(localrc, &
         ESMF_ERR_PASSTHRU, &
         ESMF_CONTEXT, rcToReturn=rc)) return


      ! Allocate memory for points
      allocate(pntList(pntDim*pntCount), stat=localrc)
      if (ESMF_LogFoundAllocError(localrc, msg="Allocating pntList", &
        ESMF_CONTEXT, rcToReturn=rc)) return   

      ! Allocate memory for pets
      allocate(petList(pntCount), stat=localrc)
      if (ESMF_LogFoundAllocError(localrc, msg="Allocating petList", &
        ESMF_CONTEXT, rcToReturn=rc)) return   

      ! Get Points 
      call ESMF_LocStreamGetPntList(locstream, coordKeyNames, pntDim, &
               pntCount, pntList, rc=localrc)
      if (ESMF_LogFoundError(localrc, &
         ESMF_ERR_PASSTHRU, &
         ESMF_CONTEXT, rcToReturn=rc)) return


      call c_ESMC_PointListCalcCartDim(coordSysLocal, pntDim, regrid_dims, localrc)
      if (ESMF_LogFoundAllocError(localrc, ESMF_ERR_PASSTHRU, &
                                ESMF_CONTEXT, rcToReturn=rc)) return

      ! Allocate memory for points in cart
      allocate(pntList_cart(regrid_dims*pntCount), stat=localrc)
      if (ESMF_LogFoundAllocError(localrc, msg="Allocating pntList_cart", &
        ESMF_CONTEXT, rcToReturn=rc)) return   

      do i=1,pntCount
        idx = (i-1) * pntDim + 1
        idx_cart = (i-1) * regrid_dims + 1
        call c_ESMC_PointListSph2CartCoord(coordSysLocal, pntDim, &
                   pntList(idx),pntList_cart(idx_cart),localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
      enddo


      ! Find out where points lie on Mesh
      call ESMF_MeshFindPnt(background, localunmappedaction, &
                                regrid_dims, pntCount, pntList_cart, &
                                petList, rc=localrc)

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

      ! Can now get rid of pntList
      deallocate(pntList)
      deallocate(pntList_cart)

      ! Create a new location stream by shifting the entries between 
      ! the pets based on petList
      ESMF_LocStreamCreateByBkgMesh=ESMF_LocStreamCreatePetList(locstream, name, &
                                  petList, rc=localrc)
      if (ESMF_LogFoundError(localrc, &
         ESMF_ERR_PASSTHRU, &
         ESMF_CONTEXT, rcToReturn=rc)) return
    
     ! Can now get rid of petList
      deallocate(petList)

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

      end function ESMF_LocStreamCreateByBkgMesh