subroutine test_sph_vec_blnr_identical(rc)
integer, intent(out) :: rc
logical :: correct
integer :: localrc
type(ESMF_Grid) :: srcGrid
type(ESMF_Grid) :: dstGrid
type(ESMF_Field) :: srcField
type(ESMF_Field) :: dstField
type(ESMF_Field) :: tmpField
type(ESMF_Field) :: xdstField
type(ESMF_Array) :: dstArray
type(ESMF_Array) :: srcArray
type(ESMF_Array) :: tmpArray
type(ESMF_RouteHandle) :: routeHandle
type(ESMF_ArraySpec) :: arrayspec
type(ESMF_VM) :: vm
real(ESMF_KIND_R8), pointer :: farrayPtrXC(:,:)
real(ESMF_KIND_R8), pointer :: farrayPtrYC(:,:)
real(ESMF_KIND_R8), pointer :: farrayPtr1DXC(:)
real(ESMF_KIND_R8), pointer :: farrayPtr1DYC(:)
real(ESMF_KIND_R8), pointer :: farrayPtr(:,:,:), farrayPtr2(:,:)
real(ESMF_KIND_R8), pointer :: xfarrayPtr(:,:,:)
real(ESMF_KIND_R8), pointer :: tmpfarrayPtr(:,:)
integer :: clbnd(2),cubnd(2)
integer :: fclbnd(3),fcubnd(3)
integer :: i1,i2,i3, index(2)
integer :: lDE, srclocalDECount, dstlocalDECount
real(ESMF_KIND_R8) :: coord(2)
character(len=ESMF_MAXSTR) :: string
integer src_nx, src_ny, dst_nx, dst_ny
integer num_arrays
real(ESMF_KIND_R8) :: dx,dy
real(ESMF_KIND_R8) :: src_dx, src_dy
real(ESMF_KIND_R8) :: dst_dx, dst_dy
real(ESMF_KIND_R8) :: lon, lat, theta, phi
real(ESMF_KIND_R8), parameter :: DEG2RAD = 3.141592653589793/180.0_ESMF_KIND_R8
integer :: localPet, petCount
! result code
integer :: finalrc
! init success flag
correct=.true.
rc=ESMF_SUCCESS
! get pet info
call ESMF_VMGetGlobal(vm, rc=localrc)
if (ESMF_LogFoundError(localrc, &
ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
call ESMF_VMGet(vm, petCount=petCount, localPet=localpet, rc=localrc)
if (ESMF_LogFoundError(localrc, &
ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
! Src Grid
srcGrid=ESMF_GridCreate1PeriDimUfrm(maxIndex=(/100,50/),&
minCornerCoord=(/0.0_ESMF_KIND_R8,-90.0_ESMF_KIND_R8/), &
maxCornerCoord=(/360.0_ESMF_KIND_R8,90.0_ESMF_KIND_R8/), &
regDecomp=(/1,petCount/), &
staggerLocList=(/ESMF_STAGGERLOC_CENTER/), &
rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
! Src Field
srcField = ESMF_FieldCreate(srcGrid, typekind=ESMF_TYPEKIND_R8, &
ungriddedLBound=(/1/), ungriddedUBound=(/2/), & ! 2D vector
staggerloc=ESMF_STAGGERLOC_CENTER, name="source", rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
! Get srcArray from Field
call ESMF_FieldGet(srcField, array=srcArray, rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
! Get number of local DEs
call ESMF_GridGet(srcGrid, localDECount=srclocalDECount, rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
! Construct Src Grid
! (Get memory and set coords for src)
do lDE=0,srclocalDECount-1
!! get coord 1
call ESMF_GridGetCoord(srcGrid, localDE=lDE, staggerLoc=ESMF_STAGGERLOC_CENTER, coordDim=1, &
farrayPtr=farrayPtr1DXC, rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
call ESMF_GridGetCoord(srcGrid, localDE=lDE, staggerLoc=ESMF_STAGGERLOC_CENTER, coordDim=2, &
farrayPtr=farrayPtr1DYC, rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
! get src pointer
call ESMF_FieldGet(srcField, lDE, farrayPtr, &
computationalLBound=fclbnd, computationalUBound=fcubnd, &
rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
!! Set Field value
do i1=fclbnd(1),fcubnd(1)
! Get X coord from Grid
lon = farrayPtr1DXC(i1)
theta = DEG2RAD*(lon)
do i2=fclbnd(2),fcubnd(2)
! Get Y coord from Grid
lat = farrayPtr1DYC(i2)
phi = DEG2RAD*(90.-lat)
! Set exact dst data and init dst field to 0.0
do i3=fclbnd(3),fcubnd(3)
! initialize source field
farrayPtr(i1,i2,i3)=1.0
enddo
enddo
enddo
enddo ! lDE
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! Destination grid
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! setup dest. grid
dstGrid=ESMF_GridCreate1PeriDimUfrm(maxIndex=(/100,50/),&
minCornerCoord=(/0.0_ESMF_KIND_R8,-90.0_ESMF_KIND_R8/), &
maxCornerCoord=(/360.0_ESMF_KIND_R8,90.0_ESMF_KIND_R8/), &
regDecomp=(/1,petCount/), &
staggerLocList=(/ESMF_STAGGERLOC_CENTER/), &
rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
! Create Fields
dstField = ESMF_FieldCreate(dstGrid, typekind=ESMF_TYPEKIND_R8, &
ungriddedLBound=(/1/), ungriddedUBound=(/2/), & ! 2D vector
staggerloc=ESMF_STAGGERLOC_CENTER, name="dest", rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
xdstField = ESMF_FieldCreate(dstGrid, typekind=ESMF_TYPEKIND_R8, &
ungriddedLBound=(/1/), ungriddedUBound=(/2/), & ! 2D vector
staggerloc=ESMF_STAGGERLOC_CENTER, name="xdest", rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
! Get dstArray from Field
call ESMF_FieldGet(dstField, array=dstArray, rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
! Get number of local DEs
call ESMF_GridGet(dstGrid, localDECount=dstlocalDECount, rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
! Get memory and set coords for dst
do lDE=0,dstlocalDECount-1
!! get coords
call ESMF_GridGetCoord(dstGrid, localDE=lDE, staggerLoc=ESMF_STAGGERLOC_CENTER, coordDim=1, &
farrayPtr=farrayPtr1DXC,rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
call ESMF_GridGetCoord(dstGrid, localDE=lDE, staggerLoc=ESMF_STAGGERLOC_CENTER, coordDim=2, &
farrayPtr=farrayPtr1DYC, rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
call ESMF_FieldGet(dstField, lDE, farrayPtr, &
computationalLBound=fclbnd, computationalUBound=fcubnd, &
rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
call ESMF_FieldGet(xdstField, lDE, xfarrayPtr, rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
!! Set Field value
do i1=fclbnd(1),fcubnd(1)
! Get X coord from Grid
lon = farrayPtr1DXC(i1)
theta = DEG2RAD*(lon)
do i2=fclbnd(2),fcubnd(2)
! Get Y coord from Grid
lat = farrayPtr1DYC(i2)
phi = DEG2RAD*(90.-lat)
! Set exact dst data and init dst field to 0.0
do i3=fclbnd(3),fcubnd(3)
!xfarrayPtr(i1,i2) = 2. + cos(theta)**2.*cos(2.*phi)
xfarrayPtr(i1,i2,i3) = 1.0 ! for now all 1.0, eventually do something better..
! initialize destination field
farrayPtr(i1,i2,i3)=0.0
enddo
enddo
enddo
enddo ! lDE
#if 0
call ESMF_GridWriteVTK(srcGrid,staggerloc=ESMF_STAGGERLOC_CENTER, &
filename="srcGrid", rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
#endif
!!! Regrid forward from the A grid to the B grid
! Regrid store
call ESMF_FieldRegridStore( &
srcField, &
dstField=dstField, &
vectorRegrid=.true., &
routeHandle=routeHandle, &
regridmethod=ESMF_REGRIDMETHOD_BILINEAR, &
rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
! Do regrid
call ESMF_FieldRegrid(srcField, dstField, routeHandle, rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
call ESMF_FieldRegridRelease(routeHandle, rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
! Check results
do lDE=0,dstlocalDECount-1
call ESMF_FieldGet(dstField, lDE, farrayPtr, &
computationalLBound=fclbnd, computationalUBound=fcubnd, &
rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
!! make sure we're not using any bad points
do i1=fclbnd(1),fcubnd(1)
do i2=fclbnd(2),fcubnd(2)
do i3=fclbnd(3),fcubnd(3)
! if working everything should be close to exact answer
if (abs(farrayPtr(i1,i2,i3)-xfarrayPtr(i1,i2,i3)) .gt. 0.001) then
!write(*,*) i1,i2,i3," ",farrayPtr(i1,i2,i3),xfarrayPtr(i1,i2,i3)
correct=.false.
endif
enddo
enddo
enddo
enddo ! lDE
#if 0
call ESMF_GridWriteVTK(srcGrid,staggerloc=ESMF_STAGGERLOC_CENTER, &
filename="srcGrid", array1=srcArray, &
rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
call ESMF_GridWriteVTK(dstGrid,staggerloc=ESMF_STAGGERLOC_CENTER, &
filename="dstGrid", array1=dstArray, &
rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
#endif
! Destroy the Fields
call ESMF_FieldDestroy(srcField, rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
call ESMF_FieldDestroy(dstField, rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
! Free the grids
call ESMF_GridDestroy(srcGrid, rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
call ESMF_GridDestroy(dstGrid, rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
! return answer based on correct flag
if (correct) then
rc=ESMF_SUCCESS
else
rc=ESMF_FAILURE
endif
end subroutine test_sph_vec_blnr_identical