subroutine ESMF_ArraySMMStoreInd4R8(srcArray, dstArray, routehandle, &
factorList, factorIndexList, keywordEnforcer, ignoreUnmatchedIndices, &
srcTermProcessing, pipelineDepth, rc)
!
! !ARGUMENTS:
type(ESMF_Array), intent(in) :: srcArray
type(ESMF_Array), intent(inout) :: dstArray
type(ESMF_RouteHandle), intent(inout) :: routehandle
real(ESMF_KIND_R8), target, intent(in) :: factorList(:)
integer, intent(in) :: factorIndexList(:,:)
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
logical, intent(in), optional :: ignoreUnmatchedIndices
integer, intent(inout), optional :: srcTermProcessing
integer, intent(inout), optional :: pipelineDepth
integer, intent(out), optional :: rc
!
!EOPI
!------------------------------------------------------------------------------
integer :: localrc ! local return code
real(ESMF_KIND_R8), pointer :: opt_factorList(:) ! helper variable
integer :: len_factorList ! helper variable
type(ESMF_InterArray) :: factorIndexListArg ! helper variable
type(ESMF_Logical) :: opt_ignoreUnmatched ! helper variable
! initialize return code; assume routine not implemented
localrc = ESMF_RC_NOT_IMPL
if (present(rc)) rc = ESMF_RC_NOT_IMPL
! Check init status of arguments
ESMF_INIT_CHECK_DEEP(ESMF_ArrayGetInit, srcArray, rc)
ESMF_INIT_CHECK_DEEP(ESMF_ArrayGetInit, dstArray, rc)
! Wrap factor arguments
len_factorList = size(factorList)
opt_factorList => factorList
factorIndexListArg = &
ESMF_InterArrayCreate(farray2D=factorIndexList, rc=localrc)
if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
! Set default flags
opt_ignoreUnmatched = ESMF_FALSE
if (present(ignoreUnmatchedIndices)) opt_ignoreUnmatched = ignoreUnmatchedIndices
! Call into the C++ interface, which will sort out optional arguments
call c_ESMC_ArraySMMStoreInd4(srcArray, dstArray, routehandle, &
ESMF_TYPEKIND_R8, opt_factorList, len_factorList, factorIndexListArg, &
opt_ignoreUnmatched, srcTermProcessing, pipelineDepth, localrc)
if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
! Garbage collection
call ESMF_InterArrayDestroy(factorIndexListArg, rc=localrc)
if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
! Mark routehandle object as being created
call ESMF_RouteHandleSetInitCreated(routehandle, rc=localrc)
if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
ESMF_CONTEXT, rcToReturn=rc)) return
! return successfully
if (present(rc)) rc = ESMF_SUCCESS
end subroutine ESMF_ArraySMMStoreInd4R8