ESMF_AttributeGetObjStateR8 Subroutine

private subroutine ESMF_AttributeGetObjStateR8(target, name, value, defaultvalue, convention, purpose, attnestflag, isPresent, rc)

Arguments

Type IntentOptional Attributes Name
type(ESMF_State), intent(in) :: target
character(len=*), intent(in) :: name
real(kind=ESMF_KIND_R8), intent(out) :: value
real(kind=ESMF_KIND_R8), intent(in), optional :: defaultvalue
character(len=*), intent(in), optional :: convention
character(len=*), intent(in), optional :: purpose
type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
logical, intent(out), optional :: isPresent
integer, intent(out), optional :: rc

Source Code

subroutine ESMF_AttributeGetObjStateR8(target, name, value, defaultvalue, convention, purpose, attnestflag, isPresent, rc)
  type(ESMF_State), intent(in) :: target
  character(len=*), intent(in) :: name
  real(ESMF_KIND_R8), intent(out) :: value
  real(ESMF_KIND_R8), intent(in), optional :: defaultvalue
  character(len=*), optional, intent(in) :: convention
  character(len=*), optional, intent(in) :: purpose
  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_InfoDescribe) :: eidesc
  type(ESMF_Info) :: info
  type(ESMF_AttNest_Flag) :: local_attnestflag
  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_StateGetInit, target, rc)

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

  call ESMF_InfoFormatKey(key, name, localrc, convention=convention, purpose=purpose)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  info = eidesc%GetInfo(target, rc=localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  local_isPresent = ESMF_InfoIsPresent(info, 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_InfoGet(info, key=key, size=size, isArray=is_array, &
      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(info, 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(info, key, value, default=defaultvalue, &
      attnestflag=local_attnestflag, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return
  end if

  deallocate(key)

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