ESMF_ArraySMMStoreInd8I8TP Subroutine

private subroutine ESMF_ArraySMMStoreInd8I8TP(srcArray, dstArray, routehandle, transposeRoutehandle, factorList, factorIndexList, keywordEnforcer, ignoreUnmatchedIndices, srcTermProcessing, pipelineDepth, rc)

Arguments

Type IntentOptional Attributes Name
type(ESMF_Array), intent(inout) :: srcArray
type(ESMF_Array), intent(inout) :: dstArray
type(ESMF_RouteHandle), intent(inout) :: routehandle
type(ESMF_RouteHandle), intent(inout) :: transposeRoutehandle
integer(kind=ESMF_KIND_I8), intent(in), target :: factorList(:)
integer(kind=ESMF_KIND_I8), intent(in) :: factorIndexList(:,:)
type(ESMF_KeywordEnforcer), optional :: keywordEnforcer
logical, intent(in), optional :: ignoreUnmatchedIndices
integer, intent(inout), optional :: srcTermProcessing
integer, intent(inout), optional :: pipelineDepth
integer, intent(out), optional :: rc

Source Code

  subroutine ESMF_ArraySMMStoreInd8I8TP(srcArray, dstArray, routehandle, &
    transposeRoutehandle, factorList, factorIndexList, keywordEnforcer, &
    ignoreUnmatchedIndices, srcTermProcessing, pipelineDepth, rc)
!
! !ARGUMENTS:
    type(ESMF_Array),              intent(inout)           :: srcArray
    type(ESMF_Array),              intent(inout)           :: dstArray
    type(ESMF_RouteHandle),        intent(inout)           :: routehandle
    type(ESMF_RouteHandle),        intent(inout)           :: transposeRoutehandle
    integer(ESMF_KIND_I8), target, intent(in)              :: factorList(:)
    integer(ESMF_KIND_I8),         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
    integer(ESMF_KIND_I8), pointer  :: opt_factorList(:)  ! helper variable
    integer                         :: len_factorList     ! helper variable
    type(ESMF_InterArray)           :: factorIndexListArg ! helper variable
    integer                         :: tupleSize, i       ! helper variable
    integer(ESMF_KIND_I8), allocatable :: transposeFIL(:,:)  ! 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(farray2DI8=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_ArraySMMStoreInd8(srcArray, dstArray, routehandle, &
      ESMF_TYPEKIND_I8, 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
    
    ! Compute the transposeRoutehandle
    ! Construct the transpose of the factorIndexList
    tupleSize = size(factorIndexList,1)
    allocate(transposeFIL(tupleSize, len_factorList))
    if (tupleSize==2) then
      do i=1, len_factorList
        transposeFIL(1,i)=factorIndexList(2,i)
        transposeFIL(2,i)=factorIndexList(1,i)
      enddo
    else if (tupleSize==4) then
      do i=1, len_factorList
        transposeFIL(1,i)=factorIndexList(3,i)
        transposeFIL(2,i)=factorIndexList(4,i)
        transposeFIL(3,i)=factorIndexList(1,i)
        transposeFIL(4,i)=factorIndexList(2,i)
      enddo
    endif
    ! wrap transposeFIL
    factorIndexListArg = &
      ESMF_InterArrayCreate(farray2DI8=transposeFIL, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
      ESMF_CONTEXT, rcToReturn=rc)) return
    ! Call into the C++ interface, which will sort out optional arguments
    call c_ESMC_ArraySMMStoreInd8(dstArray, srcArray, transposeRoutehandle, &
      ESMF_TYPEKIND_I8, 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
    deallocate(transposeFIL)
    ! Mark transposeRoutehandle object as being created
    call ESMF_RouteHandleSetInitCreated(transposeRoutehandle, 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_ArraySMMStoreInd8I8TP