ESMF_AttributeGetAttPackMeshI4 Subroutine

private subroutine ESMF_AttributeGetAttPackMeshI4(target, name, attpack, value, defaultvalue, attnestflag, isPresent, rc)

Arguments

Type IntentOptional Attributes Name
type(ESMF_Mesh), intent(in) :: target
character(len=*), intent(in) :: name
type(ESMF_AttPack), intent(inout) :: attpack
integer(kind=ESMF_KIND_I4), intent(out) :: value
integer(kind=ESMF_KIND_I4), intent(in), optional :: defaultvalue
type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
logical, intent(out), optional :: isPresent
integer, intent(out), optional :: rc

Source Code

subroutine ESMF_AttributeGetAttPackMeshI4(target, name, attpack, value, defaultvalue, attnestflag, isPresent, rc)
  ! 39.11.7
  type(ESMF_Mesh), intent(in) :: target
  character(len=*), intent(in) :: name
  type(ESMF_AttPack), intent(inout) :: attpack
  integer(ESMF_KIND_I4), intent(out) :: value
  integer(ESMF_KIND_I4), intent(in), optional :: defaultvalue
  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
  logical, intent(out), optional :: isPresent
  integer, intent(out), optional :: rc

  integer :: localrc
  character(:), allocatable :: key
  type(ESMF_AttNest_Flag) :: local_attnestflag
  logical :: debug = .false.
  integer :: size
  logical :: is_array, local_isPresent

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_RC_NOT_IMPL
  ! Check object initialization
  ESMF_INIT_CHECK_DEEP(ESMF_MeshGetInit, target, rc)

  if (present(attnestflag)) then
    local_attnestflag = attnestflag
  else
    local_attnestflag = ESMF_ATTR_DEFAULT_ATTNEST
  end if

  key = attpack%formatKey(name=name, rc=localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (debug) then
    call ESMF_LogWrite(ESMF_METHOD//": key="//TRIM(key), rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return
  end if

  local_isPresent = ESMF_InfoIsPresent(attpack%getPayload(), key, attnestflag=local_attnestflag, &
    isPointer=.true., rc=localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (present(isPresent)) then
    isPresent = local_isPresent
  endif

  ! For Attribute, we support scalar to array stuff for single element arrays.
  ! Check if the target is an array with size 1. Operate on it as if it were a
  ! scalar.
  if (local_isPresent) then
    call ESMF_InfoGetArrayMeta(attpack%getPayload(), key, is_array, size, &
      attnestflag=local_attnestflag, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return
  else
    ! Supply some default values for array and size checks to allow logical test
    ! for scalar-array implicit conversion.
    is_array = .false.
    size = 0
  endif

  if (local_isPresent .and. is_array .and. size==1) then
    call ESMF_InfoGet(attpack%getPayload(), key, value, idx=1, default=defaultvalue, &
      attnestflag=local_attnestflag, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return
  else
    call ESMF_InfoGet(attpack%getPayload(), key, value, default=defaultvalue, &
      attnestflag=local_attnestflag, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return
  end if

  if (present(rc)) rc = ESMF_SUCCESS
end subroutine ESMF_AttributeGetAttPackMeshI4