ESMF_DistGridRegDecompSetCubic Subroutine

public subroutine ESMF_DistGridRegDecompSetCubic(regDecomp, deCount, keywordEnforcer, rc)

Arguments

Type IntentOptional Attributes Name
integer, intent(out), target :: regDecomp(:)
integer, intent(in), optional :: deCount
type(ESMF_KeywordEnforcer), optional :: keywordEnforcer
integer, intent(out), optional :: rc

Source Code

  subroutine ESMF_DistGridRegDecompSetCubic(regDecomp, deCount, keywordEnforcer, rc)
!
! !ARGUMENTS:
    integer,        target, intent(out)           :: regDecomp(:)
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
    integer,                intent(in),  optional :: deCount
    integer,                intent(out), optional :: rc
!
! !DESCRIPTION:
!   Construct a regular decomposition argument for DistGrid that is most cubic
!   in {\tt dimCount} dimensions, and multiplies out to {\tt deCount} DEs. The
!   factorization is stable monotonic decreasing, ensuring that earlier entries
!   in {\tt regDecomp} are larger or equal to the later entires.
!
!   The arguments are:
!   \begin{description}
!   \item[regDecomp]
!     The regular decomposition description being constructed.
!   \item[{[deCount]}]
!     The number of DEs. Defaults to {\tt petCount}.
!   \item[{[rc]}]
!     Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!   \end{description}
!
!EOP
!------------------------------------------------------------------------------
    integer               :: localrc      ! local return code
    type(ESMF_InterArray) :: regDecompArg ! helper variable
    integer               :: optDeCount   ! helper variable
    type(ESMF_VM)         :: vm           ! helper variable

    ! initialize return code; assume routine not implemented
    localrc = ESMF_RC_NOT_IMPL
    if (present(rc)) rc = ESMF_RC_NOT_IMPL

    ! create InterArrays
    regDecompArg = ESMF_InterArrayCreate(regDecomp, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
      ESMF_CONTEXT, rcToReturn=rc)) return

    ! Deal with optional arguments
    if (present(deCount)) then
      optDeCount = deCount
    else
      call ESMF_VMGetCurrent(vm, rc=localrc)
      if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
        ESMF_CONTEXT, rcToReturn=rc)) return
      call ESMF_VMGet(vm, petCount=optDeCount, rc=localrc)
      if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
        ESMF_CONTEXT, rcToReturn=rc)) return
    endif
    
    ! call into the C++ interface, which will sort out optional arguments
    call c_ESMC_DistGridRDSetCubic(regDecompArg, optDeCount, localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
      ESMF_CONTEXT, rcToReturn=rc)) return
      
    ! return successfully
    if (present(rc)) rc = ESMF_SUCCESS
 
  end subroutine ESMF_DistGridRegDecompSetCubic