subroutine test_regridMeshSph3x3ToGrid(rc)
integer, intent(out) :: rc
logical :: correct
integer :: localrc
type(ESMF_Mesh) :: srcMesh
type(ESMF_Grid) :: dstGrid
type(ESMF_Field) :: srcField
type(ESMF_Field) :: dstField
type(ESMF_Field) :: xdstField
type(ESMF_Array) :: dstArray
type(ESMF_Array) :: lonArrayA
type(ESMF_Array) :: srcArrayA
type(ESMF_RouteHandle) :: routeHandle
type(ESMF_ArraySpec) :: arrayspec
type(ESMF_VM) :: vm
real(ESMF_KIND_R8), pointer :: farrayPtrXC(:,:), farrayPtr1D(:)
real(ESMF_KIND_R8), pointer :: farrayPtrYC(:,:)
real(ESMF_KIND_R8), pointer :: farrayPtr(:,:),farrayPtr2(:,:)
real(ESMF_KIND_R8), pointer :: xfarrayPtr(:,:)
integer :: clbnd(2),cubnd(2)
integer :: fclbnd(2),fcubnd(2)
integer :: i1,i2,i3, index(2)
integer :: lDE, localDECount
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) :: dst_minx,dst_miny
real(ESMF_KIND_R8) :: dst_maxx,dst_maxy
real(ESMF_KIND_R8) :: x,y
real(ESMF_KIND_R8) :: src_dx, src_dy
real(ESMF_KIND_R8) :: dst_dx, dst_dy
real(ESMF_KIND_R8) :: lon, lat, theta, phi, DEG2RAD
real(ESMF_KIND_R8) :: err, relErr
integer :: spherical_grid
integer, pointer :: larrayList(:)
integer :: localPet, petCount
integer, pointer :: nodeIds(:),nodeOwners(:)
real(ESMF_KIND_R8), pointer :: nodeCoords(:)
integer, pointer :: elemIds(:),elemTypes(:),elemConn(:)
integer :: numNodes, numElems
integer :: numElemConns,numQuadElems,numTriElems, numTotElems
real(ESMF_KIND_R8), pointer :: elemCoords(:)
! 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
! If we don't have 1 or 4 PETS then exit successfully
if ((petCount .ne. 1) .and. (petCount .ne. 4)) then
print*,'ERROR: test must be run using exactly 1 or 4 PETS - detected ',petCount
rc=ESMF_FAILURE
return
endif
src_nx = 3
src_ny = 3
! Establish the resolution of the grids
dst_nx = 15
dst_ny = 15
! Establish the coordinates of the grids
dst_minx = 0.1
dst_miny = 0.1
dst_maxx = 2.9
dst_maxy = 2.9
dst_dx=360.0/dst_nx
dst_dy=180.0/dst_ny
src_dx=360.0/src_nx
src_dy=180.0/src_ny
! degree to rad conversion
DEG2RAD = 3.141592653589793_ESMF_KIND_R8/180.0_ESMF_KIND_R8
! Setup mesh info depending on the
! number of PETs
if (petCount .eq. 1) then
! Fill in node data
numNodes=16
!! node ids
allocate(nodeIds(numNodes))
nodeIds=(/1,2,3,4,5,6,7,8, &
9,10,11,12,13,14,&
15,16/)
!! node Coords
allocate(nodeCoords(numNodes*2))
nodeCoords=(/0.0,0.0, & ! 1
1.0,0.0, & ! 2
2.0,0.0, & ! 3
3.0,0.0, & ! 4
0.0,1.0, & ! 5
1.0,1.0, & ! 6
2.0,1.0, & ! 7
3.0,1.0, & ! 8
0.0,2.0, & ! 9
1.0,2.0, & ! 10
2.0,2.0, & ! 11
3.0,2.0, & ! 12
0.0,3.0, & ! 13
1.0,3.0, & ! 14
2.0,3.0, & ! 15
3.0,3.0 /) ! 16
!! node owners
allocate(nodeOwners(numNodes))
nodeOwners=0 ! everything on proc 0
! Fill in elem data
numTriElems=2
numQuadElems=8
numElems=numTriElems+numQuadElems
numElemConns=3*numTriElems+4*numQuadElems
!! elem ids
allocate(elemIds(numElems))
elemIds=(/1,2,3,4,5,6,7,8,9,10/)
!! elem types
allocate(elemTypes(numElems))
elemTypes=(/ESMF_MESHELEMTYPE_QUAD, & ! 1
ESMF_MESHELEMTYPE_QUAD, & ! 2
ESMF_MESHELEMTYPE_QUAD, & ! 3
ESMF_MESHELEMTYPE_QUAD, & ! 4
ESMF_MESHELEMTYPE_QUAD, & ! 5
ESMF_MESHELEMTYPE_QUAD, & ! 6
ESMF_MESHELEMTYPE_QUAD, & ! 7
ESMF_MESHELEMTYPE_QUAD, & ! 8
ESMF_MESHELEMTYPE_TRI, & ! 9
ESMF_MESHELEMTYPE_TRI/) ! 10
!! elem coords
allocate(elemCoords(2*numElems))
elemCoords=(/0.5,0.5, & ! 1
1.5,0.5, & ! 2
2.5,0.5, & ! 3
0.5,1.5, & ! 4
1.5,1.5, & ! 5
2.5,1.5, & ! 6
0.5,2.5, & ! 7
1.5,2.5, & ! 8
2.75,2.25,& ! 9
2.25,2.75/) ! 10
!! elem conn
allocate(elemConn(numElemConns))
elemConn=(/1,2,6,5, & ! 1
2,3,7,6, & ! 2
3,4,8,7, & ! 3
5,6,10,9, & ! 4
6,7,11,10, & ! 5
7,8,12,11, & ! 6
9,10,14,13, & ! 7
10,11,15,14, & ! 8
11,12,16, & ! 9
11,16,15/) ! 10
else if (petCount .eq. 4) then
! Setup mesh data depending on PET
if (localPet .eq. 0) then
! Fill in node data
numNodes=4
!! node ids
allocate(nodeIds(numNodes))
nodeIds=(/1,2,5,6/)
!! node Coords
allocate(nodeCoords(numNodes*2))
nodeCoords=(/0.0,0.0, & ! 1
1.0,0.0, & ! 2
0.0,1.0, & ! 5
1.0,1.0 /) ! 6
!! node owners
allocate(nodeOwners(numNodes))
nodeOwners=0 ! everything on proc 0
! Fill in elem data
numTriElems=0
numQuadElems=1
numElems=numTriElems+numQuadElems
numElemConns=3*numTriElems+4*numQuadElems
!! elem ids
allocate(elemIds(numElems))
elemIds=(/1/)
!! elem types
allocate(elemTypes(numElems))
elemTypes=(/ESMF_MESHELEMTYPE_QUAD/) ! 1
!! elem coords
allocate(elemCoords(2*numElems))
elemCoords=(/0.5,0.5/) ! 1
!! elem conn
allocate(elemConn(numElemConns))
elemConn=(/1,2,4,3/) ! 1
else if (localPet .eq. 1) then
! Fill in node data
numNodes=6
!! node ids
allocate(nodeIds(numNodes))
nodeIds=(/2,3,4,6,7,8/)
!! node Coords
allocate(nodeCoords(numNodes*2))
nodeCoords=(/1.0,0.0, & ! 2
2.0,0.0, & ! 3
3.0,0.0, & ! 4
1.0,1.0, & ! 6
2.0,1.0, & ! 7
3.0,1.0 /) ! 8
!! node owners
allocate(nodeOwners(numNodes))
nodeOwners=(/0, & ! 2
1, & ! 3
1, & ! 4
0, & ! 6
1, & ! 7
1/) ! 8
! Fill in elem data
numTriElems=0
numQuadElems=2
numElems=numTriElems+numQuadElems
numElemConns=3*numTriElems+4*numQuadElems
!! elem ids
allocate(elemIds(numElems))
elemIds=(/2,3/)
!! elem types
allocate(elemTypes(numElems))
elemTypes=(/ESMF_MESHELEMTYPE_QUAD, & ! 2
ESMF_MESHELEMTYPE_QUAD/) ! 3
!! elem coords
allocate(elemCoords(2*numElems))
elemCoords=(/1.5,0.5, & ! 2
2.5,0.5/) ! 3
!! elem conn
allocate(elemConn(numElemConns))
elemConn=(/1,2,5,4, & ! 2
2,3,6,5/) ! 3
else if (localPet .eq. 2) then
! Fill in node data
numNodes=9
!! node ids
allocate(nodeIds(numNodes))
nodeIds=(/5,6,7, &
9,10,11, &
13,14,15/)
!! node Coords
allocate(nodeCoords(numNodes*2))
nodeCoords=(/0.0,1.0, & ! 5
1.0,1.0, & ! 6
2.0,1.0, & ! 7
0.0,2.0, & ! 9
1.0,2.0, & ! 10
2.0,2.0, & ! 11
0.0,3.0, & ! 13
1.0,3.0, & ! 14
2.0,3.0/) ! 15
!! node owners
allocate(nodeOwners(numNodes))
nodeOwners=(/0, & ! 5
0, & ! 6
1, & ! 7
2, & ! 9
2, & ! 10
3, & ! 11
2, & ! 13
2, & ! 14
3/) ! 15
! Fill in elem data
numTriElems=0
numQuadElems=4
numElems=numTriElems+numQuadElems
numElemConns=3*numTriElems+4*numQuadElems
!! elem ids
allocate(elemIds(numElems))
elemIds=(/4,5,7,8/)
!! elem types
allocate(elemTypes(numElems))
elemTypes=(/ESMF_MESHELEMTYPE_QUAD, & ! 4
ESMF_MESHELEMTYPE_QUAD, & ! 5
ESMF_MESHELEMTYPE_QUAD, & ! 7
ESMF_MESHELEMTYPE_QUAD/) ! 8
!! elem coords
allocate(elemCoords(2*numElems))
elemCoords=(/0.5,1.5, & ! 4
1.5,1.5, & ! 5
0.5,2.5, & ! 7
1.5,2.5/) ! 8
!! elem conn
allocate(elemConn(numElemConns))
elemConn=(/1,2,5,4, & ! 4
2,3,6,5, & ! 5
4,5,8,7, & ! 7
5,6,9,8/) ! 8
else if (localPet .eq. 3) then
! Fill in node data
numNodes=6
!! node ids
allocate(nodeIds(numNodes))
nodeIds=(/7,8,11,12,15,16/)
!! node Coords
allocate(nodeCoords(numNodes*2))
nodeCoords=(/2.0,1.0, & ! 7
3.0,1.0, & ! 8
2.0,2.0, & ! 11
3.0,2.0, & ! 12
2.0,3.0, & ! 15
3.0,3.0 /) ! 16
!! node owners
allocate(nodeOwners(numNodes))
nodeOwners=(/1, & ! 7
1, & ! 8
3, & ! 11
3, & ! 12
3, & ! 15
3/) ! 16
! Fill in elem data
numTriElems=2
numQuadElems=1
numElems=numTriElems+numQuadElems
numElemConns=3*numTriElems+4*numQuadElems
!! elem ids
allocate(elemIds(numElems))
elemIds=(/6,9,10/)
!! elem types
allocate(elemTypes(numElems))
elemTypes=(/ESMF_MESHELEMTYPE_QUAD, & ! 6
ESMF_MESHELEMTYPE_TRI, & ! 9
ESMF_MESHELEMTYPE_TRI/) ! 10
!! elem coords
allocate(elemCoords(2*numElems))
elemCoords=(/2.5,1.5, & ! 6
2.75,2.25,& ! 9
2.25,2.75/) ! 10
!! elem conn
allocate(elemConn(numElemConns))
elemConn=(/1,2,4,3, & ! 6
3,4,6, & ! 9
3,6,5/) ! 10
endif
endif
! Create Mesh structure in 1 step
srcMesh=ESMF_MeshCreate(parametricDim=2,spatialDim=2, &
nodeIds=nodeIds, nodeCoords=nodeCoords, &
nodeOwners=nodeOwners, elementIds=elemIds,&
elementTypes=elemTypes, elementConn=elemConn, &
elementCoords=elemCoords, &
coordSys=ESMF_COORDSYS_SPH_DEG, &
rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
! Create source field
call ESMF_ArraySpecSet(arrayspec, 1, ESMF_TYPEKIND_R8, rc=rc)
srcField = ESMF_FieldCreate(srcMesh, arrayspec, &
name="source", rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
! Load test data into the source Field
! Should only be 1 localDE
call ESMF_FieldGet(srcField, 0, farrayPtr1D, rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
! set interpolated function
i2=1
do i1=1,numNodes
if (nodeOwners(i1) .eq. localPet) then
! Get coordinates
lon=nodeCoords(2*i1-1)
lat=nodeCoords(2*i1)
theta = DEG2RAD*(lon)
phi = DEG2RAD*(90.-lat)
! set exact src data
farrayPtr1D(i2) = 2. + cos(theta)**2.*cos(2.*phi)
! Advance to next owner
i2=i2+1
endif
enddo
! deallocate node data
deallocate(nodeIds)
deallocate(nodeCoords)
deallocate(nodeOwners)
! deallocate elem data
deallocate(elemIds)
deallocate(elemTypes)
deallocate(elemConn)
! setup dest. grid
dstGrid=ESMF_GridCreateNoPeriDim(minIndex=(/1,1/),maxIndex=(/dst_nx,dst_ny/),regDecomp=(/2,2/), &
coordSys=ESMF_COORDSYS_SPH_DEG,indexflag=ESMF_INDEX_GLOBAL, &
rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
! Create source/destination fields
call ESMF_ArraySpecSet(arrayspec, 2, ESMF_TYPEKIND_R8, rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
dstField = ESMF_FieldCreate(dstGrid, arrayspec, &
staggerloc=ESMF_STAGGERLOC_CENTER, name="dest", rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
xdstField = ESMF_FieldCreate(dstGrid, arrayspec, &
staggerloc=ESMF_STAGGERLOC_CENTER, name="xdest", rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
call ESMF_GridAddCoord(dstGrid, staggerloc=ESMF_STAGGERLOC_CENTER, rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
! Get number of local DEs
call ESMF_GridGet(dstGrid, localDECount=localDECount, rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
! Get arrays
! dstArray
call ESMF_FieldGet(dstField, array=dstArray, rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
! srcArrayA
call ESMF_FieldGet(srcField, array=srcArrayA, rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! Destination grid
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! Get memory and set coords for dst
do lDE=0,localDECount-1
!! get coords
call ESMF_GridGetCoord(dstGrid, localDE=lDE, staggerLoc=ESMF_STAGGERLOC_CENTER, coordDim=1, &
computationalLBound=clbnd, computationalUBound=cubnd, farrayPtr=farrayPtrXC, rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
call ESMF_GridGetCoord(dstGrid, localDE=lDE, staggerLoc=ESMF_STAGGERLOC_CENTER, coordDim=2, &
computationalLBound=clbnd, computationalUBound=cubnd, farrayPtr=farrayPtrYC, 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 coords
do i1=clbnd(1),cubnd(1)
do i2=clbnd(2),cubnd(2)
! Set destination coordinates
farrayPtrXC(i1,i2) = ((dst_maxx-dst_minx)*REAL(i1-1)/REAL(dst_nx-1))+dst_minx
farrayPtrYC(i1,i2) = ((dst_maxy-dst_miny)*REAL(i2-1)/REAL(dst_ny-1))+dst_miny
lon = farrayPtrXC(i1,i2)
lat = farrayPtrYc(i1,i2)
theta = DEG2RAD*(lon)
phi = DEG2RAD*(90.-lat)
! set exact dst data
xfarrayPtr(i1,i2) = 2. + cos(theta)**2.*cos(2.*phi)
! initialize destination field
farrayPtr(i1,i2)=0.0
enddo
enddo
enddo ! lDE
!!! Regrid forward from the A grid to the B grid
! Regrid store
call ESMF_FieldRegridStore( &
srcField, &
dstField=dstField, &
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 error
do lDE=0,localDECount-1
!! get coords
call ESMF_GridGetCoord(dstGrid, localDE=lDE, staggerLoc=ESMF_STAGGERLOC_CENTER, coordDim=1, &
farrayPtr=farrayPtrXC, 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=farrayPtrYC, rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
call ESMF_FieldGet(dstField, lDE, farrayPtr, computationalLBound=clbnd, &
computationalUBound=cubnd, 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
!! check error
do i1=clbnd(1),cubnd(1)
do i2=clbnd(2),cubnd(2)
! compute relative error
err=abs(farrayPtr(i1,i2)-xfarrayPtr(i1,i2))
if (xfarrayPtr(i1,i2) .ne. 0.0) then
relErr=err/xfarrayPtr(i1,i2)
else
relErr=err
endif
! Return error if relative error too big
if (relErr > .001) then
correct=.false.
endif
enddo
enddo
enddo ! lDE
! Uncomment these calls to see some actual regrid results
#if 0
spherical_grid = 0
call ESMF_MeshIO(vm, srcMesh, ESMF_STAGGERLOC_EDGE1, &
"srcmesh", srcArrayA, rc=localrc, &
spherical=spherical_grid)
call ESMF_MeshIO(vm, dstGrid, ESMF_STAGGERLOC_CENTER, &
"dstmesh", dstArray, rc=localrc, &
spherical=spherical_grid)
#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
call ESMF_FieldDestroy(xdstField, rc=localrc)
if (localrc /=ESMF_SUCCESS) then
rc=ESMF_FAILURE
return
endif
! Free the grids
call ESMF_MeshDestroy(srcMesh, 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_regridMeshSph3x3ToGrid