ESMF_ArrayRedistStoreNF Subroutine

private subroutine ESMF_ArrayRedistStoreNF(srcArray, dstArray, routehandle, keywordEnforcer, srcToDstTransposeMap, ignoreUnmatchedIndices, pipelineDepth, rc)

Arguments

Type IntentOptional Attributes Name
type(ESMF_Array), intent(in) :: srcArray
type(ESMF_Array), intent(inout) :: dstArray
type(ESMF_RouteHandle), intent(inout) :: routehandle
type(ESMF_KeywordEnforcer), optional :: keywordEnforcer
integer, intent(in), optional :: srcToDstTransposeMap(:)
logical, intent(in), optional :: ignoreUnmatchedIndices
integer, intent(inout), optional :: pipelineDepth
integer, intent(out), optional :: rc

Source Code

  subroutine ESMF_ArrayRedistStoreNF(srcArray, dstArray, routehandle, &
    keywordEnforcer, srcToDstTransposeMap, ignoreUnmatchedIndices, &
    pipelineDepth, rc)
!
! !ARGUMENTS:
    type(ESMF_Array),       intent(in)              :: srcArray
    type(ESMF_Array),       intent(inout)           :: dstArray
    type(ESMF_RouteHandle), intent(inout)           :: routehandle
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
    integer,                intent(in),    optional :: srcToDstTransposeMap(:)
    logical,                intent(in),    optional :: ignoreUnmatchedIndices
    integer,                intent(inout), optional :: pipelineDepth
    integer,                intent(out),   optional :: rc
!
! !STATUS:
! \begin{itemize}
! \item\apiStatusCompatibleVersion{5.2.0r}
! \item\apiStatusModifiedSinceVersion{5.2.0r}
! \begin{description}
! \item[6.1.0] Added argument {\tt pipelineDepth}.
!              The new argument provide access to the tuning parameter
!              affecting the sparse matrix execution.
! \item[7.0.0] Added argument {\tt transposeRoutehandle} to allow a handle to
!              the transposed redist operation to be returned.\newline
!              Added argument {\tt ignoreUnmatchedIndices} to support situations 
!              where not all elements between source and destination Arrays 
!              match.
! \item[7.1.0r] Removed argument {\tt transposeRoutehandle} and provide it
!              via interface overloading instead. This allows argument 
!              {\tt srcArray} to stay strictly intent(in) for this entry point.
! \end{description}
! \end{itemize}
!
! !DESCRIPTION:
! \label{ArrayRedistStoreNF}
! {\tt ESMF\_ArrayRedistStore()} is a collective method across all PETs of the
! current Component. The interface of the method is overloaded, allowing 
! -- in principle -- each PET to call into {\tt ESMF\_ArrayRedistStore()}
! through a different entry point. Restrictions apply as to which combinations
! are sensible. All other combinations result in ESMF run time errors. The
! complete semantics of the {\tt ESMF\_ArrayRedistStore()} method, as provided
! through the separate entry points shown in \ref{ArrayRedistStoreTK} and
! \ref{ArrayRedistStoreNF}, is described in the following paragraphs as a whole.
!
! Store an Array redistribution operation from {\tt srcArray} to {\tt dstArray}.
! Interface \ref{ArrayRedistStoreTK} allows PETs to specify a {\tt factor}
! argument. PETs not specifying a {\tt factor} argument call into interface
! \ref{ArrayRedistStoreNF}. If multiple PETs specify the {\tt factor} argument,
! its type and kind, as well as its value must match across all PETs. If none
! of the PETs specify a {\tt factor} argument the default will be a factor of
! 1. The resulting factor is applied to all of the source data during
! redistribution, allowing scaling of the data, e.g. for unit transformation.
!
! Both {\tt srcArray} and {\tt dstArray} are interpreted as sequentialized 
! vectors. The sequence is defined by the order of DistGrid dimensions and the
! order of tiles within the DistGrid or by user-supplied arbitrary sequence
! indices. See section \ref{Array:SparseMatMul} for details on the definition
! of {\em sequence indices}.
!
! Source Array, destination Array, and the factor may be of different
! <type><kind>. Further, source and destination Arrays may differ in shape,
! however, the number of elements must match. 
!
! The default redistribution operation, when {\tt srcToDstTransposeMap} is not
! specified, corresponds to the identity mapping: each element of the
! sequentialized source Array is copied to the sequentialized
! destination Array element in order.
!
! \begin{sloppypar}
! If the {\tt srcToDstTransposeMap} argument is provided it must be identical
! across all PETs. The {\tt srcToDstTransposeMap} allows source and destination
! Array dimensions to be transposed during the redistribution. To support this
! option, the number of source and destination Array dimensions must be equal
! and the size of the associated dimensions must match.
! See section \ref{Array:Redist:TransposeMode} for more details about the
! use of the {\tt srcToDstTransposeMap} argument.
! \end{sloppypar}
!
! It is erroneous to specify the identical Array object for {\tt srcArray} and
! {\tt dstArray} arguments. 
!
!   The routine returns an {\tt ESMF\_RouteHandle} that can be used to call 
!   {\tt ESMF\_ArrayRedist()} on any pair of Arrays that matches 
!   {\tt srcArray} and {\tt dstArray} in {\em type}, {\em kind}, and 
!   memory layout of the {\em distributed} dimensions. However, the size,
!   number, and index order of {\em undistributed} dimensions may be different.
!   See section \ref{RH:Reusability} for a more detailed discussion of
!   RouteHandle reusability.
!
! This call is {\em collective} across the current VM.  
!
!   \begin{description}
!
!   \item [srcArray]
!     {\tt ESMF\_Array} with source data.
!
!   \item [dstArray]
!     {\tt ESMF\_Array} with destination data. The data in this Array may be
!     destroyed by this call.
!
!   \item [routehandle]
!     Handle to the precomputed Route.
!
!   \item [{[srcToDstTransposeMap]}]
!     A list with as many entries as there are dimensions in {\tt srcArray}, or
!     {\tt tileCount} times this many entries.
!     Each entry maps the corresponding {\tt srcArray} dimension against the
!     specified {\tt dstArray} dimension. Mixing distributed and
!     undistributed dimensions is supported.
!     Negative entries reverse the order of elements along the specified
!     dimension when going from source to destination.
!     When providing $rank \times tileCount$ elements in
!     {\tt srcToDstTransposeMap},  each block of size {\tt rank} is associated
!     with the corresponding tile (in order), and interpreted as the
!     tile-specific transpose map.
!
!   \item [{[ignoreUnmatchedIndices]}]
!     A logical flag that affects the behavior for when not all elements match
!     between the {\tt srcArray} and {\tt dstArray} side. The default setting
!     is {\tt .false.}, indicating that it is an error when such a situation is 
!     encountered. Setting {\tt ignoreUnmatchedIndices} to {\tt .true.} ignores
!     unmatched indices.
!
!   \item [{[pipelineDepth]}]
!     The {\tt pipelineDepth} parameter controls how many messages a PET
!     may have outstanding during a redist exchange. Larger values
!     of {\tt pipelineDepth} typically lead to better performance. However,
!     on some systems too large a value may lead to performance degradation,
!     or runtime errors.
!
!     Note that the pipeline depth has no effect on the bit-for-bit
!     reproducibility of the results. However, it may affect the performance
!     reproducibility of the exchange.
!
!     The {\tt ESMF\_ArraySMMStore()} method implements an auto-tuning scheme
!     for the {\tt pipelineDepth} parameter. The intent on the 
!     {\tt pipelineDepth} argument is "{\tt inout}" in order to 
!     support both overriding and accessing the auto-tuning parameter.
!     If an argument $>= 0$ is specified, it is used for the 
!     {\tt pipelineDepth} parameter, and the auto-tuning phase is skipped.
!     In this case the {\tt pipelineDepth} argument is not modified on
!     return. If the provided argument is $< 0$, the {\tt pipelineDepth}
!     parameter is determined internally using the auto-tuning scheme. In this
!     case the {\tt pipelineDepth} argument is re-set to the internally
!     determined value on return. Auto-tuning is also used if the optional 
!     {\tt pipelineDepth} argument is omitted.
!
!   \item [{[rc]}]
!     Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!   \end{description}
!
!EOP
!------------------------------------------------------------------------------
    integer               :: localrc                  ! local return code
    type(ESMF_InterArray) :: srcToDstTransposeMapArg  ! index helper
    type(ESMF_Logical)    :: opt_ignoreUnmatched      ! helper variable
    integer, allocatable  :: dstToSrcTransposeMap(:)  ! helper variable
    integer               :: i

    ! 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)

    ! Deal with (optional) array arguments
    srcToDstTransposeMapArg = ESMF_InterArrayCreate(srcToDstTransposeMap, &
      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_ArrayRedistStoreNF(srcArray, dstArray, routehandle, &
      srcToDstTransposeMapArg, opt_ignoreUnmatched, pipelineDepth, 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
    
    ! garbage collection
    call ESMF_InterArrayDestroy(srcToDstTransposeMapArg, 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_ArrayRedistStoreNF