function ESMF_DistGridCreateDBAI(arbSeqIndexList, arbDim, &
minIndexPTile, maxIndexPTile, keywordEnforcer, rc)
!
! !RETURN VALUE:
type(ESMF_DistGrid) :: ESMF_DistGridCreateDBAI
!
! !ARGUMENTS:
integer, intent(in) :: arbSeqIndexList(:)
integer, intent(in) :: arbDim
integer, intent(in) :: minIndexPTile(:)
integer, intent(in) :: maxIndexPTile(:)
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
integer, intent(out), optional :: rc
!
! !STATUS:
! \begin{itemize}
! \item\apiStatusCompatibleVersion{5.2.0r}
! \end{itemize}
!
! !DESCRIPTION:
! Create an {\tt ESMF\_DistGrid} of {\tt dimCount} $1+n$, where
! $n=$ {\tt size(minIndexPTile)} = {\tt size(maxIndexPTile)}.
!
! The resulting DistGrid will have a 1D distribution determined by the
! PET-local {\tt arbSeqIndexList}. The PET-local size of the
! {\tt arbSeqIndexList} argument determines the number of local elements
! along the arbitrarily distributed dimension in the created DistGrid. The
! sequence indices must be unique across all PETs. The associated,
! automatically created DELayout will have 1 DE per PET across all PETs of
! the current VM.
!
! In addition to the arbitrarily distributed dimension, regular DistGrid
! dimensions can be specified in {\tt minIndexPTile} and {\tt maxIndexPTile}. The
! $n$ dimensional subspace spanned by the regular dimensions is "multiplied"
! with the arbitrary dimension on each DE, to form a $1+n$ dimensional
! total index space described by the DistGrid object. The {\tt arbDim}
! argument allows to specify which dimension in the resulting DistGrid
! corresponds to the arbitrarily distributed one.
!
! This is a {\em collective} method across the current VM.
!
! The arguments are:
! \begin{description}
! \item[arbSeqIndexList]
! List of arbitrary sequence indices that reside on the local PET.
! \item[arbDim]
! Dimension of the arbitrary distribution.
! \item[minIndexPTile]
! Index space tuple of the lower corner of the tile. The
! arbitrary dimension is {\em not} included in this tile
! \item[maxIndexPTile]
! Index space tuple of the upper corner of the tile. The
! arbitrary dimension is {\em not} included in this tile
! \item[{[rc]}]
! Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
! \end{description}
!
!EOP
!------------------------------------------------------------------------------
integer :: localrc ! local return code
type(ESMF_DistGrid) :: distgrid ! opaque pointer to new C++ DistGrid
type(ESMF_VM) :: vm ! opaque pointer to VM object
type(ESMF_InterArray) :: indicesAux ! index helper
integer :: localSize(1) ! number of local indices
integer, allocatable :: globalSizes(:) ! array of all sizes
integer :: petCount ! num pets
integer, allocatable :: deblock(:,:,:) ! Array of sizes
integer :: i, j, jj, csum ! loop variable
integer, allocatable :: minC(:), maxC(:)! min/max corner
integer, allocatable :: collocationPDim(:)
integer :: dimCount ! number of dimension
! initialize return code; assume routine not implemented
localrc = ESMF_RC_NOT_IMPL
if (present(rc)) rc = ESMF_RC_NOT_IMPL
! invalidate return value
distgrid%this = ESMF_NULL_POINTER
ESMF_DistGridCreateDBAI = distgrid
! check input
dimCount = size(minIndexPTile)
if (dimCount /= size(maxIndexPTile)) then
call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_SIZE, &
msg="size(minIndexPTile) must match size(maxIndexPTile)", &
ESMF_CONTEXT, rcToReturn=rc)
return
endif
if (arbDim < 1 .or. arbDim > dimCount+1) then
call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_VALUE, &
msg="arbDim out of range", &
ESMF_CONTEXT, rcToReturn=rc)
return
endif
! get VM and related information
call ESMF_VMGetCurrent(vm, rc=localrc)
if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
call ESMF_VMGet(vm, petCount=petCount, rc=localrc)
if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
! number of local indices
localSize(1) = size(arbSeqIndexList)
! gather all sizes locally
allocate(globalSizes(petCount))
call ESMF_VMAllGather(vm, localSize, globalSizes, 1, rc=localrc)
if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
! set up the deblocks
allocate(deblock(dimCount+1,2,petCount))
csum = 0
do i=1,petCount
deblock(arbDim,1,i) = csum + 1 ! min
csum = csum + globalSizes(i)
deblock(arbDim,2,i) = csum ! max
do j=1,dimCount+1
if (j==arbDim) cycle
jj=j
if (j>arbDim) jj=jj-1
deblock(j,1,i) = minIndexPTile(jj) ! min
deblock(j,2,i) = maxIndexPTile(jj) ! max
enddo
enddo
! create fitting DistGrid
allocate(minC(dimCount+1))
allocate(maxC(dimCount+1))
do i=1, dimCount+1
minC(i) = deblock(i,1,1)
maxC(i) = deblock(i,2,petCount)
enddo
distgrid = ESMF_DistGridCreate(minC, maxC, deBlockList=deblock, &
indexflag=ESMF_INDEX_GLOBAL, rc=localrc)
if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
! garbage collection
deallocate(deblock)
deallocate(globalSizes)
deallocate(minC)
deallocate(maxC)
! set return value
ESMF_DistGridCreateDBAI = distgrid
! set collocations to separate arbDim from the reset
allocate(collocationPDim(dimCount+1))
collocationPDim = 2 ! initialize
collocationPDim(arbDim) = 1 ! arbDim singled out as collocation "1"
call ESMF_DistGridSet(distgrid, collocationPDim=collocationPDim, &
rc=localrc)
if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
deallocate(collocationPDim)
! prepare to set local arbitrary sequence indices
indicesAux = ESMF_InterArrayCreate(farray1D=arbSeqIndexList, rc=localrc)
if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
! set local arbitrary sequence indices in DistGrid object
! localDe=0, collocation=1, i.e. arbDim's collocation
call c_ESMC_DistGridSetArbSeqIndex(distgrid, indicesAux, 0, 1, localrc)
if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
! garbage collection
call ESMF_InterArrayDestroy(indicesAux, rc=localrc)
if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
! Set init code
ESMF_INIT_SET_CREATED(ESMF_DistGridCreateDBAI)
! return successfully
if (present(rc)) rc = ESMF_SUCCESS
end function ESMF_DistGridCreateDBAI