test_regrid_store_from_file Subroutine

public subroutine test_regrid_store_from_file(rc)

Arguments

Type IntentOptional Attributes Name
integer, intent(out) :: rc

Source Code

  subroutine test_regrid_store_from_file(rc)
  integer, intent(out)  :: rc
  logical :: correct
  integer :: localrc
  type(ESMF_VM) :: vm
  type(ESMF_Grid) :: srcGrid, dstGrid
  type(ESMF_Field) :: srcField, dstField
  type(ESMF_RouteHandle) :: routeHandle
  real(ESMF_KIND_R8), pointer :: src(:,:), dst(:,:)
  real(ESMF_KIND_R8) :: lon, lat, theta, phi
  real(ESMF_KIND_R8), pointer :: s_x(:,:), s_y(:,:), d_x(:,:), d_y(:,:)
  integer(ESMF_KIND_I4) :: lbnd(2), ubnd(2)
  integer :: localPet, petCount
  integer :: i, j, m, n
  real(ESMF_KIND_R8), pointer :: factorList(:) 
  integer(ESMF_KIND_I4), pointer :: factorIndexList(:,:) 

  ! init success flag
  correct=.true.
  rc=ESMF_FAILURE

  ! get pet info
  call ESMF_VMGetGlobal(vm, rc=localrc)
  if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, &
    line=__LINE__, file=FILENAME, rcToReturn=rc)) return

  call ESMF_VMGet(vm, petCount=petCount, localPet=localpet, rc=localrc)
  if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, &
    line=__LINE__, file=FILENAME, rcToReturn=rc)) return

  m = 8
  n = 9

  srcGrid = ESMF_GridCreateNoPeriDim(maxIndex=(/n,n/), &
                                     coordSys=ESMF_COORDSYS_CART, &
                                     indexflag=ESMF_INDEX_GLOBAL, &
                                     rc=localrc)
  if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, &
    line=__LINE__, file=FILENAME, rcToReturn=rc)) return
  dstGrid = ESMF_GridCreateNoPeriDim(maxIndex=(/m,m/), &
                                     coordSys=ESMF_COORDSYS_CART, &
                                     indexflag=ESMF_INDEX_GLOBAL,  &
                                     rc=localrc)
  if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, &
    line=__LINE__, file=FILENAME, rcToReturn=rc)) return

  call ESMF_GridAddCoord(srcGrid, staggerloc=ESMF_STAGGERLOC_CENTER, rc=localrc)
  if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, &
    line=__LINE__, file=FILENAME, rcToReturn=rc)) return
  call ESMF_GridAddCoord(dstGrid, staggerloc=ESMF_STAGGERLOC_CENTER, rc=localrc)
  if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, &
    line=__LINE__, file=FILENAME, rcToReturn=rc)) return

  call ESMF_GridGetCoord(srcGrid, 1, staggerloc=ESMF_STAGGERLOC_CENTER, &
                         computationalLBound=lbnd, computationalUBound=ubnd, &
                         farrayPtr=s_x, rc=localrc)
  if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, &
    line=__LINE__, file=FILENAME, rcToReturn=rc)) return

  call ESMF_GridGetCoord(srcGrid, 2, staggerloc=ESMF_STAGGERLOC_CENTER, &
                         computationalLBound=lbnd, computationalUBound=ubnd, &
                         farrayPtr=s_y, rc=localrc)
  if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, &
    line=__LINE__, file=FILENAME, rcToReturn=rc)) return

  do i = lbnd(1), ubnd(1)
    do j = lbnd(2), ubnd(2)
        s_x(i, j) = i - 0.5
        s_y(i, j) = j - 0.5
    enddo
  enddo


  call ESMF_GridGetCoord(dstGrid, 1, staggerloc=ESMF_STAGGERLOC_CENTER, &
                         computationalLBound=lbnd, computationalUBound=ubnd, &
                         farrayPtr=d_x, rc=localrc)
  if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, &
    line=__LINE__, file=FILENAME, rcToReturn=rc)) return

  call ESMF_GridGetCoord(dstGrid, 2, staggerloc=ESMF_STAGGERLOC_CENTER, &
                         computationalLBound=lbnd, computationalUBound=ubnd, &
                         farrayPtr=d_y, rc=localrc)
  if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, &
    line=__LINE__, file=FILENAME, rcToReturn=rc)) return

    do i = lbnd(1), ubnd(1)
      do j = lbnd(2), ubnd(2)
          d_x(i, j) = i
          d_y(i, j) = j
      enddo
    enddo

  ! Create source/destination fields
   srcField = ESMF_FieldCreate(dstGrid, typekind=ESMF_TYPEKIND_R8, &
                               staggerloc=ESMF_STAGGERLOC_CENTER, &
                               name="source", rc=localrc)
  if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, &
    line=__LINE__, file=FILENAME, rcToReturn=rc)) return

   dstField = ESMF_FieldCreate(dstGrid, ESMF_TYPEKIND_R8, &
                               staggerloc=ESMF_STAGGERLOC_CENTER, &
                               name="dest", rc=localrc)
  if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, &
    line=__LINE__, file=FILENAME, rcToReturn=rc)) return


  ! Get arrays
  ! dstArray
  call ESMF_FieldGet(dstField, farrayPtr=dst, rc=localrc)
  if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, &
    line=__LINE__, file=FILENAME, rcToReturn=rc)) return
  dst = 0.

  ! srcArray
  call ESMF_FieldGet(srcField, farrayPtr=src, rc=localrc)
  if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, &
    line=__LINE__, file=FILENAME, rcToReturn=rc)) return
  src = 42.

  ! Do regrid
  call ESMF_FieldRegridStore(srcField, dstField, routehandle=routeHandle, &
    factorList=factorList, factorIndexList=factorIndexList, rc=localrc)
  if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, &
    line=__LINE__, file=FILENAME, rcToReturn=rc)) return

  call ESMF_SparseMatrixWrite(factorList, factorIndexList, &
                              "weights_esmf_smmsff.nc", rc=localrc)
  if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, &
    line=__LINE__, file=FILENAME, rcToReturn=rc)) return

  ! SMM store
  call ESMF_FieldSMMStore(srcField, dstField, "weights_esmf_smmsff.nc", &
                          routeHandle, rc=localrc)
  if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, &
    line=__LINE__, file=FILENAME, rcToReturn=rc)) return

  call ESMF_FieldRegrid(srcField, dstField, routeHandle, rc=localrc)
  if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, &
    line=__LINE__, file=FILENAME, rcToReturn=rc)) return

  do i = lbnd(1), ubnd(1)
    do j = lbnd(2), ubnd(2)
        if (dst(i, j) /= 42.) then 
          correct = .false.
        endif
    enddo
  enddo

  deallocate(factorList)
  deallocate(factorIndexList)

  call ESMF_FieldRegridRelease(routeHandle, rc=localrc)
  if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, &
    line=__LINE__, file=FILENAME, rcToReturn=rc)) return

  ! Destroy the Fields
   call ESMF_FieldDestroy(srcField, rc=localrc)
  if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, &
    line=__LINE__, file=FILENAME, rcToReturn=rc)) return

   call ESMF_FieldDestroy(dstField, rc=localrc)
  if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, &
    line=__LINE__, file=FILENAME, rcToReturn=rc)) return

  ! Free the grids
  call ESMF_GridDestroy(srcGrid, rc=localrc)
  if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, &
    line=__LINE__, file=FILENAME, rcToReturn=rc)) return

  call ESMF_GridDestroy(dstGrid, rc=localrc)
  if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, &
    line=__LINE__, file=FILENAME, rcToReturn=rc)) return

  ! return answer based on correct flag
  if (correct) then
    rc=ESMF_SUCCESS
  else
    rc=ESMF_FAILURE
  endif

  end subroutine test_regrid_store_from_file