ESMF_PointListCreateFrmLocStream Function

private function ESMF_PointListCreateFrmLocStream(locstream, maskValues, addOrigCoords, rc)

Arguments

Type IntentOptional Attributes Name
type(ESMF_LocStream), intent(in) :: locstream
integer(kind=ESMF_KIND_I4), intent(in), optional :: maskValues(:)
logical, optional :: addOrigCoords
integer, intent(out), optional :: rc

Return Value type(ESMF_PointList)


Source Code

  function ESMF_PointListCreateFrmLocStream(locstream, maskValues, addOrigCoords, rc)
!
! !RETURN VALUE:
    type(ESMF_PointList) :: ESMF_PointListCreateFrmLocStream
!
! !ARGUMENTS:
    type(ESMF_LocStream), intent(in)  :: locstream
    integer(ESMF_KIND_I4), intent(in), optional   :: maskValues(:)
    logical, optional                 :: addOrigCoords    
    integer, intent(out), optional    :: rc
!
! !DESCRIPTION:
!   Allocates memory for a new {\tt ESMF\_PointList} object and
!   constructs its internals from input LocStream.
!
!   The arguments are:
!   \begin{description}
!   \item[{locStream}]
!     The Location Stream to get the information from to create the PointList.
!   \item[{maskValues}]
!     Values to set as masked
!   \item[{[addOrigCoords]}]
!     Also put the original (spherical) coordinate values in the PointList.
!     Defaults to false.    
!   \item[{[rc]}]
!     Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!   \end{description}
!
!EOPI
!------------------------------------------------------------------------------
    integer             :: localrc      ! local return code
    integer             :: i,j,k,dimcount,regrid_dims, num_local_pts, num_maskValues
    integer             :: cl,cu,cc
    integer             :: localDECount,lDE
    logical             :: d3Present,maskPresent,masked_value
    character(len=ESMF_MAXSTR)    :: keystring
    character(len=ESMF_MAXSTR), pointer  :: myKeyNames(:)
    type(ESMF_Array),pointer           :: myKeyArrays(:)
    type(ESMF_Array),pointer           :: thisKeyArray
    integer,allocatable           :: seqInd(:)
    real(ESMF_KIND_R8), pointer :: farrayPtrX(:),farrayPtrY(:),farrayPtrZ(:)
    real(ESMF_KIND_R8), allocatable :: mycoords(:), cart_coords(:)
    integer(ESMF_KIND_I4), pointer :: maskarray(:)
    type(ESMF_DistGrid) :: distgridOut
    integer :: seqCount
    type(ESMF_Array)       :: XArr,YArr,ZArr,MArr
    integer :: origCoordDim
    type(ESMF_CoordSys_Flag) :: coordSysLocal, origCoordSys
    type(ESMF_PointList)  :: pointlist
    logical :: localAddOrigCoords
    
    
    ! initialize return code; assume routine not implemented
    localrc = ESMF_RC_NOT_IMPL
    if (present(rc)) rc = ESMF_RC_NOT_IMPL
    pointlist%this = ESMF_NULL_POINTER


    ! Set default of localAddOrigCoords
    localAddOrigCoords=.false.
    if (present(addOrigCoords)) localAddOrigCoords=addOrigCoords
    
    ! Get information from LocStream
    call ESMF_LocStreamGet(locstream, coordSys=coordSysLocal, localDECount=localDECount, &
                           distgrid=distgridOut, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
      ESMF_CONTEXT, rcToReturn=rc)) return

    dimcount = 0
    !dimension 1
    if (coordSysLocal .eq. ESMF_COORDSYS_CART) then
      keystring='ESMF:X'
    else
      keystring='ESMF:Lon'
    endif
    call ESMF_LocStreamGetKey(locstream, keyName=keystring, &
                              keyArray=XArr, rc=localrc)
    if (.not. ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
      ESMF_CONTEXT, rcToReturn=rc)) then
      dimcount = dimcount+1
    else
      return
    endif

    !dimension 2
    if (coordSysLocal .eq. ESMF_COORDSYS_CART) then
      keystring='ESMF:Y'
    else
      keystring='ESMF:Lat'
    endif
    call ESMF_LocStreamGetKey(locstream, keyName=keystring, &
                              keyArray=YArr, rc=localrc)
    if (.not. ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
      ESMF_CONTEXT, rcToReturn=rc)) then
      dimcount = dimcount+1
    else
      return
    endif

    !dimension 3?
    if (coordSysLocal .eq. ESMF_COORDSYS_CART) then
      keystring='ESMF:Z'
    else
      keystring='ESMF:Radius'
    endif
    !only deal with third dimension if it's present
    call ESMF_LocStreamGetKey(locstream, keyName=keystring, &
                              isPresent=d3Present, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
      ESMF_CONTEXT, rcToReturn=rc)) return
    if (d3Present) then
      call ESMF_LocStreamGetKey(locstream, keyName=keystring, &
                                keyArray=ZArr, rc=localrc)
      if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
      dimcount = dimcount+1
    endif

    allocate(mycoords(dimcount), stat=localrc)
    if (ESMF_LogFoundAllocError(localrc, ESMF_ERR_PASSTHRU, &
      ESMF_CONTEXT, rcToReturn=rc)) return
    call c_ESMC_PointListCalcCartDim(coordSysLocal, dimcount, regrid_dims, localrc)
    if (ESMF_LogFoundAllocError(localrc, ESMF_ERR_PASSTHRU, &
      ESMF_CONTEXT, rcToReturn=rc)) return
    allocate(cart_coords(regrid_dims), stat=localrc)
    if (ESMF_LogFoundAllocError(localrc, ESMF_ERR_PASSTHRU, &
      ESMF_CONTEXT, rcToReturn=rc)) return


    !only deal with mask info if it's present
    call ESMF_LocStreamGetKey(locstream, keyName='ESMF:Mask', &
                              isPresent=maskPresent, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
      ESMF_CONTEXT, rcToReturn=rc)) return
    if (maskPresent) then
      call ESMF_LocStreamGetKey(locstream, keyName='ESMF:Mask', &
                                keyArray=MArr, rc=localrc)
      if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
        ESMF_CONTEXT, rcToReturn=rc)) return
    endif

    if (present(maskValues)) then
      if (.not. maskPresent) then
        call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, &
             msg='- LocStream has no masking info for use with specified mask values', &
             ESMF_CONTEXT, rcToReturn=rc)
        return
      endif
      num_maskValues = size(maskValues)
    else
      num_maskValues = 0
    endif

    !must count the points first to create the pointlist with the proper size
    num_local_pts=0
    do lDE=0,localDECount-1
      if (maskPresent) then
        call ESMF_ArrayGet(MArr, localDE=lDE, farrayPtr=maskarray, rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
      endif

      call  ESMF_LocStreamGetBounds(locstream, localDE=lDE, &
                              computationalLBound=cl, computationalUBound=cu, &
                              rc=localrc)
      if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return

      do j=cl,cu
        masked_value=.false.

        !redundant if statement here is needed to avoid bug in Intel optimizer, which
        !occurs with versions 15.0.3, 16.0.0 and possibly others
        if (num_maskValues .gt. 0) then
        do k=1,num_maskValues
          if (maskArray(j) .eq. maskValues(k)) then
            masked_value=.true.
            exit
          endif
        enddo
        endif
        if (.not. masked_value) num_local_pts = num_local_pts + 1
      enddo
    enddo

    ! Set original coordinates info based on whether they are requested
    origCoordDim=0 ! 0 -> don't add original coords
    if (localAddOrigCoords) origCoordDim=dimcount

    origCoordSys=ESMF_COORDSYS_UNINIT
    if (localAddOrigCoords) origCoordSys=coordSysLocal

    
    ! Create pointlist 
    pointlist = ESMF_PointListCreate(maxpts=num_local_pts,numdims=regrid_dims, &
         origCoordDim=origCoordDim, origCoordSys=origCoordSys, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
        ESMF_CONTEXT, rcToReturn=rc)) return

    !now we add the points
    do lDE=0,localDECount-1
      call ESMF_ArrayGet(XArr, localDE=lDE, farrayPtr=farrayPtrX, rc=localrc)
      if (ESMF_LogFoundError(localrc, &
          msg="expecting coordinate keys to be REAL*8", &
          ESMF_CONTEXT, rcToReturn=rc)) return
      call ESMF_ArrayGet(YArr, localDE=lDE, farrayPtr=farrayPtrY, rc=localrc)
      if (ESMF_LogFoundError(localrc, &
          msg="expecting coordinate keys to be REAL*8", &
          ESMF_CONTEXT, rcToReturn=rc)) return
      if (dimcount .eq. 3) then
        call ESMF_ArrayGet(ZArr, localDE=lDE, farrayPtr=farrayPtrZ, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
          msg="expecting coordinate keys to be REAL*8", &
          ESMF_CONTEXT, rcToReturn=rc)) return
      endif

      !Allocate space for seqInd
      allocate(seqInd(size(farrayPtrX)), stat=localrc)
      if (ESMF_LogFoundAllocError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
      call ESMF_DistGridGet(distgridOut, localDe=lDE, &
                            seqIndexList=seqInd, rc=localrc)
      if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return

      call  ESMF_LocStreamGetBounds(locstream, localDE=lDE, &
                              computationalLBound=cl, computationalUBound=cu, &
                              rc=localrc)
      if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
      cc=cu-cl+1

      if (size(farrayPtrX) .ne. size(farrayPtrY) .or. size(farrayPtrX) .ne. cc) then
        call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, &
          msg='- coord arrays must be equal in size and greater than size 0', &
          ESMF_CONTEXT, rcToReturn=rc)
        return
      endif

      if (dimcount .eq. 3 ) then
        if (size(farrayPtrZ) .ne. size(farrayPtrX)) then
          call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, &
            msg='- coord arrays must be equal in size and greater than size 0', &
            ESMF_CONTEXT, rcToReturn=rc)
          return
        endif
      endif

      if (maskPresent) then
        call ESMF_ArrayGet(MArr, localDE=lDE, farrayPtr=maskarray, rc=localrc)
        if (ESMF_LogFoundError(localrc, &
            msg="expecting mask key to be INTEGER*4", &
            ESMF_CONTEXT, rcToReturn=rc)) return

        if (size(maskarray) .ne. size(farrayPtrX)) then
          call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_WRONG, &
            msg='- mask array must be equal in size to coordinate arrays', &
            ESMF_CONTEXT, rcToReturn=rc)
          return
        endif
      endif

      do j=cl,cu
        masked_value=.false.

        !redundant if statement here is needed to avoid bug in Intel optimizer, which
        !occurs with versions 15.0.3, 16.0.0 and possibly others
        if (num_maskValues .gt. 0) then
        do k=1,num_maskValues
          if (maskArray(j) .eq. maskValues(k)) then
            masked_value=.true.
            exit
          endif
        enddo
        endif
        if (.not. masked_value) then
          mycoords(1)=farrayPtrX(j)
          mycoords(2)=farrayPtrY(j)
          if (dimcount .eq. 3) mycoords(3)=farrayPtrZ(j)

          call c_ESMC_PointListSph2CartCoord(coordSysLocal, dimcount, mycoords, &
                                             cart_coords, localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
              ESMF_CONTEXT, rcToReturn=rc)) return

          if (localAddOrigCoords) then
             call ESMF_PointListAdd(pointlist=pointlist,id=seqInd(j-cl+1), &
                  loc_coords=cart_coords, loc_orig_coords=mycoords, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                  ESMF_CONTEXT, rcToReturn=rc)) return
          else
             call ESMF_PointListAdd(pointlist=pointlist,id=seqInd(j-cl+1), &
                  loc_coords=cart_coords, rc=localrc)
             if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                  ESMF_CONTEXT, rcToReturn=rc)) return            
          endif
       endif
      enddo
      deallocate(seqInd, stat=localrc)
      if (ESMF_LogFoundAllocError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
    enddo

    deallocate(mycoords, stat=localrc)
    if (ESMF_LogFoundAllocError(localrc, ESMF_ERR_PASSTHRU, &
      ESMF_CONTEXT, rcToReturn=rc)) return
    deallocate(cart_coords, stat=localrc)
    if (ESMF_LogFoundAllocError(localrc, ESMF_ERR_PASSTHRU, &
      ESMF_CONTEXT, rcToReturn=rc)) return

    ! Set return values
    ESMF_PointListCreateFrmLocStream = pointlist

    ESMF_INIT_SET_CREATED(ESMF_PointListCreateFrmLocStream)

    ! Return successfully
    if (present(rc)) rc = ESMF_SUCCESS

  end function ESMF_PointListCreateFrmLocStream