subroutine computeRedistAreaMesh(mesh, vm, petNo, petCnt, area, rc)
type(ESMF_Mesh) :: mesh
type(ESMF_VM) :: VM
integer :: petNo,petCnt
real (ESMF_KIND_R8), pointer :: area(:)
integer :: rc
real (ESMF_KIND_R8), pointer :: localArea(:)
integer :: localrc
integer :: localElemCount,i
integer (ESMF_KIND_I4) :: localCount(1), globalCount(1)
integer :: totalCount
type(ESMF_DistGrid) :: distgrid, justPet0Distgrid
type(ESMF_Array) :: areaArray, justPet0Array
type(ESMF_RouteHandle) :: rh
! Get local size of mesh areas
call ESMF_MeshGet(mesh, numOwnedElements=localElemCount, rc=localrc)
if (ESMF_LogFoundError(localrc, &
ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
! allocate space for areas
allocate(localArea(localElemCount))
! Get local Areas
call ESMF_MeshGetElemArea(mesh, areaList=localArea, rc=localrc)
if (ESMF_LogFoundError(localrc, &
ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
! The element disgrid is for the split elements, thus the following code
! doesn't work with split element
call ESMF_MeshGet(mesh, elementDistgrid=distgrid, rc=localrc)
if (ESMF_LogFoundError(localrc, &
ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
areaArray = ESMF_ArrayCreate(distgrid, localArea, rc=localrc)
if (ESMF_LogFoundError(localrc, &
ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
! Get total size
localCount(1)=localElemCount
globalCount(1)=0
call ESMF_VMAllReduce(vm,localCount,globalCount,count=1, &
reduceflag=ESMF_REDUCE_SUM, rc=localrc)
if (ESMF_LogFoundError(localrc, &
ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
totalCount=globalCount(1)
! Create distgrid with everything on PET 0
justPet0Distgrid = ESMF_DistGridCreate((/1/),(/totalCount/), regDecomp=(/1/),&
rc=localrc)
if (ESMF_LogFoundError(localrc, &
ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
if (PetNo == 0) then
! Allocate final area list
allocate(area(totalCount))
else
allocate(area(0))
endif
! Create array from distgrid
justPet0Array=ESMF_ArrayCreate(justPet0Distgrid, &
farrayPtr=area, rc=localrc)
if (ESMF_LogFoundError(localrc, &
ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
! Redist from one to the other
call ESMF_ArrayRedistStore(srcArray=areaArray, &
dstArray=justPet0Array, &
routehandle=rh, rc=localrc)
if (ESMF_LogFoundError(localrc, &
ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
call ESMF_ArrayRedist(srcArray=areaArray, &
dstArray=justPet0Array, &
routehandle=rh, rc=localrc)
if (ESMF_LogFoundError(localrc, &
ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
call ESMF_ArrayRedistRelease(routehandle=rh, rc=localrc)
if (ESMF_LogFoundError(localrc, &
ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
! Get rid of helper variables
call ESMF_ArrayDestroy(areaArray)
deallocate(localArea)
! call ESMF_ArrayDestroy(justPet0Array)
! call ESMF_DistGridDestroy(justPet0Distgrid)
end subroutine computeRedistAreaMesh