ESMF_Info.F90 Source File


Source Code

! $Id$
!
! Earth System Modeling Framework
! Copyright (c) 2002-2023, University Corporation for Atmospheric Research,
! Massachusetts Institute of Technology, Geophysical Fluid Dynamics
! Laboratory, University of Michigan, National Centers for Environmental
! Prediction, Los Alamos National Laboratory, Argonne National Laboratory,
! NASA Goddard Space Flight Center.
! Licensed under the University of Illinois-NCSA License.
!
!==============================================================================
#define ESMF_FILENAME "ESMF_Info.F90"
!==============================================================================

!==============================================================================
! DO NOT EDIT THIS FILE DIRECTLY. IT IS GENERATED FROM A JINJA2 TEMPLATE FILE.
!   - Template files located in scripts/jinja2_templating/templates
!   - The template file is: ESMF_Info.jinja2
!   - All code edits must be done in that template file then re-generated
!   - See scripts/jinja2_templating/README.md for guidance
!==============================================================================

module ESMF_InfoMod

!==============================================================================
!
! This file contains the Fortran wrapper code for the C++ implementation of
! the Info class.
!
!------------------------------------------------------------------------------

! INCLUDES
#include "ESMF.h"

! =============================================================================
!BOPI
! !MODULE: ESMF_InfoMod
!

!   Fortran API wrapper of C++ implemenation of Info
!
!------------------------------------------------------------------------------

! !USES:
use ESMF_UtilTypesMod     ! ESMF utility types
use ESMF_InitMacrosMod    ! ESMF initializer macros
use ESMF_BaseMod          ! ESMF base class
use ESMF_LogErrMod        ! ESMF error handling
use ESMF_VMMod
use ESMF_HConfigMod
use iso_c_binding

implicit none

include "ESMF_InfoCDef.F90"
include "ESMF_InfoCDefGeneric.F90"

! =============================================================================

type ESMF_Info
  type(C_PTR) :: ptr = C_NULL_PTR
  logical :: is_view = .false.
end type ESMF_Info

! =============================================================================

!BOPI
! !INTERFACE:
interface ESMF_InfoCreate
! !PRIVATE MEMBER FUNCTIONS:
  module procedure ESMF_InfoCreateEmpty
  module procedure ESMF_InfoCreateByKey
  module procedure ESMF_InfoCreateByParse
  module procedure ESMF_InfoCreateFromInfo
! !DESCRIPTION:
! This interface provides a single entry point for \texttt{ESMF\_Info} create
! functions.
!EOPI
end interface ESMF_InfoCreate

interface ESMF_InfoGet
  module procedure ESMF_InfoInquire
  module procedure ESMF_InfoGetI4
  module procedure ESMF_InfoGetI8
  module procedure ESMF_InfoGetR4
  module procedure ESMF_InfoGetR8
  module procedure ESMF_InfoGetLG
  module procedure ESMF_InfoGetCH
  module procedure ESMF_InfoGetArrayI4
  module procedure ESMF_InfoGetArrayI8
  module procedure ESMF_InfoGetArrayR4
  module procedure ESMF_InfoGetArrayR8
  module procedure ESMF_InfoGetArrayLG
  module procedure ESMF_InfoGetArrayCH
end interface ESMF_InfoGet

interface ESMF_InfoGetAlloc
  module procedure ESMF_InfoGetArrayI4Alloc
  module procedure ESMF_InfoGetArrayI8Alloc
  module procedure ESMF_InfoGetArrayR4Alloc
  module procedure ESMF_InfoGetArrayR8Alloc
  module procedure ESMF_InfoGetArrayLGAlloc
  module procedure ESMF_InfoGetArrayCHAlloc
end interface ESMF_InfoGetAlloc

interface ESMF_InfoSet
  module procedure ESMF_InfoSetI4
  module procedure ESMF_InfoSetI8
  module procedure ESMF_InfoSetR4
  module procedure ESMF_InfoSetR8
  module procedure ESMF_InfoSetLG
  module procedure ESMF_InfoSetCH
  module procedure ESMF_InfoSetINFO
  module procedure ESMF_InfoSetHConfig
  module procedure ESMF_InfoSetArrayI4
  module procedure ESMF_InfoSetArrayI8
  module procedure ESMF_InfoSetArrayR4
  module procedure ESMF_InfoSetArrayR8
  module procedure ESMF_InfoSetArrayLG
  module procedure ESMF_InfoSetArrayCH
end interface ESMF_InfoSet

!BOP
! !IROUTINE: ESMF_InfoAssignment(=) - Info assignment
!
! !INTERFACE:
!   interface assignment(=)
!   info1 = info2
!
! !ARGUMENTS:
!   type(ESMF_Info) :: info1
!   type(ESMF_Info) :: info2
!
! !STATUS:
! \begin{itemize}
! \item\apiStatusCompatibleVersion{5.2.0r}
! \end{itemize}
!
! !DESCRIPTION:
!   Assign info1 as an alias to the same ESMF Info object in memory
!   as info2. If info2 is invalid, then info1 will be equally invalid after
!   the assignment.
!
!   The arguments are:
!   \begin{description}
!   \item[info1]
!     The {\tt ESMF\_Info} object on the left hand side of the assignment.
!   \item[info2]
!     The {\tt ESMF\_Info} object on the right hand side of the assignment.
!   \end{description}
!
!EOP

!BOP
! !IROUTINE: ESMF_InfoOperator(==) - Info equality operator
!
! !INTERFACE:
interface operator(==)
! !RETURN VALUE:
!   logical :: result
!
! !ARGUMENTS:
!    type(ESMF_Info), intent(in) :: info1
!    type(ESMF_Info), intent(in) :: info2
!
! !DESCRIPTION:
!   Test if the contents of two \texttt{ESMF\_Info} objects are equal.
!
!   The arguments are:
!   \begin{description}
!     \item [info1]
!       The \texttt{ESMF\_Info} object on the left hand side of the operation.
!     \item [info1]
!       The \texttt{ESMF\_Info} object on the right hand side of the operation.
!   \end{description}
!EOP
  procedure ESMF_InfoEqual
end interface operator(==)

!BOP
! !IROUTINE: ESMF_InfoOperator(/=) - Info not equal operator
!
! !INTERFACE:
interface operator(/=)
! !RETURN VALUE:
!   logical :: result
!
! !ARGUMENTS:
!    type(ESMF_Info), intent(in) :: info1
!    type(ESMF_Info), intent(in) :: info2
!
! !DESCRIPTION:
!   Test if the contents of two \texttt{ESMF\_Info} objects are not equal.
!
!   The arguments are:
!   \begin{description}
!     \item [info1]
!       The \texttt{ESMF\_Info} object on the left hand side of the operation.
!     \item [info1]
!       The \texttt{ESMF\_Info} object on the right hand side of the operation.
!   \end{description}
!EOP
  procedure ESMF_InfoNotEqual
end interface operator(/=)

! =============================================================================

private

public operator(==)
public operator(/=)
public ESMF_Info
public ESMF_InfoCreate
public ESMF_InfoDestroy
public ESMF_InfoGet
public ESMF_InfoGetCharAlloc
public ESMF_InfoGetAlloc
public ESMF_InfoSet
public ESMF_InfoRemove
public ESMF_InfoGetFromBase
public ESMF_InfoGetFromPointer
public ESMF_InfoSetNULL
public ESMF_InfoSetDirty
public ESMF_InfoIsSet
public ESMF_InfoIsPresent
public ESMF_InfoPrint
public ESMF_InfoDump
public ESMF_InfoUpdate
public ESMF_InfoReadJSON
public ESMF_InfoWriteJSON
public ESMF_InfoBroadcast
public ESMF_InfoGetTK
public ESMF_InfoGetArrayMeta

public c_info_base_sync
public c_info_copyforattribute
public c_info_copyforattribute_reference

!------------------------------------------------------------------------------
!------------------------------------------------------------------------------

!------------------------------------------------------------------------------
! The following line turns the CVS identifier string into a printable variable.
character(*), parameter, private :: version = '$Id$'
!------------------------------------------------------------------------------

contains  !====================================================================


!------------------------------------------------------------------------------
! Equality Operators ----------------------------------------------------------
!------------------------------------------------------------------------------

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoEqual()"
impure elemental function ESMF_InfoEqual(lhs, rhs) result(is_equal)
  type(ESMF_Info), intent(in) :: lhs
  type(ESMF_Info), intent(in) :: rhs
  logical :: is_equal

  integer :: localrc
  integer(C_INT) :: local_is_equal

  local_is_equal = 0  !false
  call c_info_is_equal(lhs%ptr, rhs%ptr, local_is_equal, localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT)) return

  is_equal = .false.
  if (local_is_equal == 1) is_equal = .true.
end function ESMF_InfoEqual

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoNotEqual()"
impure elemental function ESMF_InfoNotEqual(lhs, rhs) result(is_equal)
  type(ESMF_Info), intent(in) :: lhs
  type(ESMF_Info), intent(in) :: rhs
  logical :: is_equal
  is_equal = .not. ESMF_InfoEqual(lhs, rhs)
end function ESMF_InfoNotEqual

! -----------------------------------------------------------------------------

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoBroadcast()"
!BOP
! !IROUTINE: ESMF_InfoBroadcast - Broadcast Info contents
! \label{esmf_infobroadcast}
!
! !INTERFACE:
subroutine ESMF_InfoBroadcast(info, rootPet, keywordEnforcer, rc)
! !ARGUMENTS:
  type(ESMF_Info), intent(inout) :: info
  integer, intent(in) :: rootPet
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  integer, intent(out), optional :: rc
!
! !DESCRIPTION:
!     Broadcast an \texttt{ESMF\_Info} object collectively across the current VM.
!
!     Users wishing to synchronize via broadcast an attribute hierarchy associated
!     with an ESMF object should consult the \texttt{ESMF\_InfoSync} documentation
!     \ref{esmf_infosync}
!
!     The arguments are:
!     \begin{description}
!     \item [info]
!       The \texttt{ESMF\_Info} object that is the source (on \textit{rootPet}) or the
!       destination object to populate (on all other PETs). On destination PETs,
!       the structure of \textit{info} is overwritten with data from \textit{rootPet}.
!     \item [rootPet]
!       The root PET identifier.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

  integer :: localrc

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE

  call c_info_broadcast(info%ptr, rootPet, localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

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

!------------------------------------------------------------------------------
! Create ----------------------------------------------------------------------
!------------------------------------------------------------------------------

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoCreateEmpty()"
!BOP
! !IROUTINE: ESMF_InfoCreate - Create a new Info object
!
! !INTERFACE:
  ! Private name; call using ESMF_InfoCreate()
function ESMF_InfoCreateEmpty(rc)
! !ARGUMENTS:
  integer, intent(out), optional :: rc
! !RETURN VALUE:
  type(ESMF_Info) :: ESMF_InfoCreateEmpty
!
! !DESCRIPTION:
!     Create an \texttt{ESMF\_Info} object. This object must be destroyed using
!     \texttt{ESMF\_InfoDestroy} to free its memory allocation
!
!     The arguments are:
!     \begin{description}
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

  integer :: localrc

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE

  ESMF_InfoCreateEmpty%ptr = c_info_create(localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (present(rc)) rc = ESMF_SUCCESS
end function ESMF_InfoCreateEmpty

! -----------------------------------------------------------------------------

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoCreateByKey()"
!BOP
! !IROUTINE: ESMF_InfoCreate - Create a new Info object using a key
!
! !INTERFACE:
  ! Private name; call using ESMF_InfoCreate()
function ESMF_InfoCreateByKey(info, key, keywordEnforcer, rc)
! !ARGUMENTS:
  type(ESMF_Info), intent(in) :: info
  character(len=*), intent(in) :: key
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  integer, intent(out), optional :: rc
! !RETURN VALUE:
  type(ESMF_Info) :: ESMF_InfoCreateByKey
!
! !DESCRIPTION:
!     Create an \texttt{ESMF\_Info} object from a location in \textit{info}
!     defined by \textit{key}. Returned object is a deep copy. The value associated
!     with \texttt{key} must be a nested object (i.e. a collection of key/value
!     pairs).
!
!     The arguments are:
!     \begin{description}
!     \item [info]
!       The \texttt{ESMF\_Info} object providing source data.
!     \item [key]
!       String key to access in \texttt{ESMF\_Info} storage. See section \ref{info_key_format}
!       for an overview of the key format.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

  integer :: localrc

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE

  ESMF_InfoCreateByKey%ptr = c_info_create_by_key(info%ptr, trim(key)//C_NULL_CHAR, localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (present(rc)) rc = ESMF_SUCCESS
end function ESMF_InfoCreateByKey

! -----------------------------------------------------------------------------

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoCreateFromInfo()"
!BOP
! !IROUTINE: ESMF_InfoCreate - Create an Info object from another Info object
!
! !INTERFACE:
  ! Private name; call using ESMF_InfoCreate()
function ESMF_InfoCreateFromInfo(info, keywordEnforcer, rc)
! !ARGUMENTS:
  type(ESMF_Info), intent(in) :: info
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  integer, intent(out), optional :: rc
! !RETURN VALUE:
  type(ESMF_Info) :: ESMF_InfoCreateFromInfo
!
! !DESCRIPTION:
!     Create an \texttt{ESMF\_Info} object from another \texttt{ESMF\_Info} object.
!     The returned object is a deep copy of the source object.
!
!     The arguments are:
!     \begin{description}
!     \item [info]
!       The \texttt{ESMF\_Info} object acting as the source data.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

  integer :: localrc

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE

  ESMF_InfoCreateFromInfo%ptr = c_info_copy(info%ptr, localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (present(rc)) rc = ESMF_SUCCESS
end function ESMF_InfoCreateFromInfo

! -----------------------------------------------------------------------------

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoCreateByParse()"
!BOP
! !IROUTINE: ESMF_InfoCreate - Create a new Info object by string parsing
!
! !INTERFACE:
  ! Private name; call using ESMF_InfoCreate()
function ESMF_InfoCreateByParse(jsonString, keywordEnforcer, rc)
! !ARGUMENTS:
  character(len=*), intent(in) :: jsonString
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  integer, intent(out), optional :: rc
! !RETURN VALUE:
  type(ESMF_Info) :: ESMF_InfoCreateByParse
!
! !DESCRIPTION:
!     Create an \texttt{ESMF\_Info} object by parsing a JSON-formatted string.
!
!     The arguments are:
!     \begin{description}
!     \item [jsonString]
!       The string to parse.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

  integer :: localrc

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE

  ESMF_InfoCreateByParse%ptr = c_info_create_by_parse(trim(jsonString)//C_NULL_CHAR, localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (present(rc)) rc = ESMF_SUCCESS
end function ESMF_InfoCreateByParse

!------------------------------------------------------------------------------
! Destroy ---------------------------------------------------------------------
!------------------------------------------------------------------------------

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoDestroy()"
!BOP
! !IROUTINE: ESMF_InfoDestroy - Destroy an Info object
!
! !INTERFACE:
subroutine ESMF_InfoDestroy(info, keywordEnforcer, rc)
! !ARGUMENTS:
  type(ESMF_Info), intent(inout) :: info
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  integer, intent(out), optional :: rc
!
! !DESCRIPTION:
!     Destroy an \texttt{ESMF\_Info} object. Destroying an \texttt{ESMF\_Info}
!     object created internally by an ESMF object results in an error
!
!     The arguments are:
!     \begin{description}
!     \item [info]
!       Target \texttt{ESMF\_Info} object.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

  integer :: localrc

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE

  if (info%is_view) then
    if (ESMF_LogFoundError(localrc, msg="Object is a view and may not be destroyed. Destroy its host object.", &
      ESMF_CONTEXT, rcToReturn=rc)) return
  endif

  call c_info_destroy(info%ptr, localrc)
  info%ptr = C_NULL_PTR
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

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

!------------------------------------------------------------------------------

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoDump()"
!BOP
! !IROUTINE: ESMF_InfoDump - Dump Info contents to string
!
! !INTERFACE:
function ESMF_InfoDump(info, keywordEnforcer, key, indent, rc) result(output)
! !ARGUMENTS:
  type(ESMF_Info), intent(in) :: info
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  character(*), intent(in), optional :: key
  integer, intent(in), optional :: indent
  integer, intent(out), optional :: rc
! RESULT:
  character(:), allocatable :: output
!
! !DESCRIPTION:
!     Dump the contents of an \texttt{ESMF\_Info} object as a JSON string.
!
!     The arguments are:
!     \begin{description}
!     \item [info]
!       Target \texttt{ESMF\_Info} object.
!     \item [{[key]}]
!       String key to access in \texttt{ESMF\_Info} storage. See section \ref{info_key_format}
!       for an overview of the key format.
!     \item [{[indent]}]
!       Default is 0. Specifying an indentation greater than 0 will result in a
!       "pretty print" for JSON output string (string includes new line breaks).
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

  character(:), allocatable :: l_key
  integer :: dump_length, localrc, local_indent

  localrc = ESMF_RC_NOT_IMPL
  if (present(rc)) rc = ESMF_RC_NOT_IMPL

  if (present(key)) then
    l_key = key
  else
    l_key = ""
  endif
  if (present(indent)) then
    local_indent = indent
  else
    local_indent = 0
  endif

  call c_info_dump_len(info%ptr, dump_length, localrc, local_indent)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT)) return

  allocate(character(dump_length)::output)

  call c_info_dump(info%ptr, output, localrc, local_indent)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT)) return

  if (present(rc)) rc = ESMF_SUCCESS
end function ESMF_InfoDump

!------------------------------------------------------------------------------
! Get (Scalar) ----------------------------------------------------------------
!------------------------------------------------------------------------------

!BOP
! !IROUTINE: ESMF_InfoGet - Get a numeric, logical, or fixed-size character value
!
! !INTERFACE:
!subroutine ESMF_InfoGet(info, key, value, keywordEnforcer, default, idx, attnestflag, rc)
!
! !ARGUMENTS:
!  type(ESMF_Info), intent(in) :: info
!  character(len=*), intent(in) :: key
!  <value>, see below for supported value
!type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
!  <default, optional> see below for supported default value
!  integer, intent(in), optional :: idx
!  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
!  integer, intent(out), optional :: rc
!
! !DESCRIPTION:
!     Get a value from an \texttt{ESMF\_Info} object using a key. If the key is
!     not found, \textit{rc} will not equal \texttt{ESMF\_SUCCESS}. The returned
!     value is always a copy including gets with a \textit{default}.
!
!     Overloaded \textit{value} for the following types:
!     \begin{itemize}
!       \item \texttt{integer(kind=ESMF\_KIND\_I4)}
!       \item \texttt{integer(kind=ESMF\_KIND\_I8)}
!       \item \texttt{real(kind=ESMF\_KIND\_R4)}
!       \item \texttt{real(kind=ESMF\_KIND\_R8)}
!       \item \texttt{logical}
!       \item \texttt{character(:)}
!     \end{itemize}
!
!     The arguments are:
!     \begin{description}
!     \item [info]
!       Target \texttt{ESMF\_Info} object.
!     \item [key]
!       String key to access in \texttt{ESMF\_Info} storage. See section \ref{info_key_format}
!       for an overview of the key format.
!     \item [value]
!       The output value associated with the key.
!     \item [{[default]}]
!       A default value to use if the key is not present in the target \texttt{ESMF\_Info}
!       object. Must be the same typekind and size as \textit{value}.
!     \item [{[idx]}]
!       An integer index to get if the target key's value is a list.
!     \item [{[attnestflag]}]
!       Setting to \texttt{ESMF\_ATTNEST\_ON} triggers a recursive search. The
!       first instance of the key (searching by depth) will be found in the hierarchy.
!       Default is \texttt{ESMF\_ATTNEST\_OFF}.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

!BOP
! !IROUTINE: ESMF_InfoGetCharAlloc - Get an allocatable character value
!
! !INTERFACE:
!subroutine ESMF_InfoGetCharAlloc(info, key, value, keywordEnforcer, default, idx, attnestflag, rc)
!
! !ARGUMENTS:
!  type(ESMF_Info), intent(in) :: info
!  character(len=*), intent(in) :: key
!  character(:), allocatable, intent(out) :: value
!type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
!  character(len=*), intent(in), optional :: default
!  integer, intent(in), optional :: idx
!  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
!  integer, intent(out), optional :: rc
!
! !DESCRIPTION:
!     Get a value from an \texttt{ESMF\_Info} object using a key. If the key is
!     not found, \textit{rc} will not equal \texttt{ESMF\_SUCCESS}. The returned
!     value is always a copy including gets with a \textit{default}.
!
!     The arguments are:
!     \begin{description}
!     \item [info]
!       Target \texttt{ESMF\_Info} object.
!     \item [key]
!       String key to access in \texttt{ESMF\_Info} storage. See section \ref{info_key_format}
!       for an overview of the key format.
!     \item [value]
!       The output value associated with the key.
!     \item [{[default]}]
!       A default value to use if the key is not present in the target \texttt{ESMF\_Info}
!       object. Must be the same typekind and size as \textit{value}.
!     \item [{[idx]}]
!       An integer index to get if the target key's value is a list.
!     \item [{[attnestflag]}]
!       Setting to \texttt{ESMF\_ATTNEST\_ON} triggers a recursive search. The
!       first instance of the key (searching by depth) will be found in the hierarchy.
!       Default is \texttt{ESMF\_ATTNEST\_OFF}.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! NOTE: Documentation stub located above for generic interface compliance.
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoGetR4()"
subroutine ESMF_InfoGetR4(info, key, value, keywordEnforcer, default, idx, attnestflag, rc)
  type(ESMF_Info), intent(in) :: info
  character(len=*), intent(in) :: key
  real(ESMF_KIND_R4), intent(out) :: value
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  real(ESMF_KIND_R4), intent(in), optional :: default
  integer, intent(in), optional :: idx
  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
  integer, intent(out), optional :: rc

  integer :: localrc
  real(C_FLOAT), target :: local_default
  integer(C_INT), target :: local_idx
  type(C_PTR) :: local_default_ptr, local_idx_ptr
  integer(C_INT) :: recursive, strlen_only

  ! Set up local return code
  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE
  recursive = 0 !false
  strlen_only = 0 !false

  ! Handle optional arguments for C ###########################################

  if (present(default)) then
    local_default = default
    local_default_ptr = C_LOC(local_default)
  else
    local_default_ptr = C_NULL_PTR
  end if
  if (present(idx)) then
    local_idx = idx - 1  ! Shift to C (zero-based) indexing
    local_idx_ptr = C_LOC(local_idx)
  else
    local_idx_ptr = C_NULL_PTR
  end if
  if (present(attnestflag)) then
    if (attnestflag%value==ESMF_ATTNEST_ON%value) recursive = 1 !true
  end if

  ! Call C ####################################################################

  call c_info_get_R4(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    value, &
    localrc, &
    local_default_ptr, &
    local_idx_ptr, &
    recursive)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

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

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoGetR8()"
subroutine ESMF_InfoGetR8(info, key, value, keywordEnforcer, default, idx, attnestflag, rc)
  type(ESMF_Info), intent(in) :: info
  character(len=*), intent(in) :: key
  real(ESMF_KIND_R8), intent(out) :: value
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  real(ESMF_KIND_R8), intent(in), optional :: default
  integer, intent(in), optional :: idx
  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
  integer, intent(out), optional :: rc

  integer :: localrc
  real(C_DOUBLE), target :: local_default
  integer(C_INT), target :: local_idx
  type(C_PTR) :: local_default_ptr, local_idx_ptr
  integer(C_INT) :: recursive, strlen_only

  ! Set up local return code
  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE
  recursive = 0 !false
  strlen_only = 0 !false

  ! Handle optional arguments for C ###########################################

  if (present(default)) then
    local_default = default
    local_default_ptr = C_LOC(local_default)
  else
    local_default_ptr = C_NULL_PTR
  end if
  if (present(idx)) then
    local_idx = idx - 1  ! Shift to C (zero-based) indexing
    local_idx_ptr = C_LOC(local_idx)
  else
    local_idx_ptr = C_NULL_PTR
  end if
  if (present(attnestflag)) then
    if (attnestflag%value==ESMF_ATTNEST_ON%value) recursive = 1 !true
  end if

  ! Call C ####################################################################

  call c_info_get_R8(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    value, &
    localrc, &
    local_default_ptr, &
    local_idx_ptr, &
    recursive)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

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

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoGetI4()"
subroutine ESMF_InfoGetI4(info, key, value, keywordEnforcer, default, idx, attnestflag, rc)
  type(ESMF_Info), intent(in) :: info
  character(len=*), intent(in) :: key
  integer(ESMF_KIND_I4), intent(out) :: value
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  integer(ESMF_KIND_I4), intent(in), optional :: default
  integer, intent(in), optional :: idx
  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
  integer, intent(out), optional :: rc

  integer :: localrc
  integer(C_INT), target :: local_default
  integer(C_INT), target :: local_idx
  type(C_PTR) :: local_default_ptr, local_idx_ptr
  integer(C_INT) :: recursive, strlen_only

  ! Set up local return code
  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE
  recursive = 0 !false
  strlen_only = 0 !false

  ! Handle optional arguments for C ###########################################

  if (present(default)) then
    local_default = default
    local_default_ptr = C_LOC(local_default)
  else
    local_default_ptr = C_NULL_PTR
  end if
  if (present(idx)) then
    local_idx = idx - 1  ! Shift to C (zero-based) indexing
    local_idx_ptr = C_LOC(local_idx)
  else
    local_idx_ptr = C_NULL_PTR
  end if
  if (present(attnestflag)) then
    if (attnestflag%value==ESMF_ATTNEST_ON%value) recursive = 1 !true
  end if

  ! Call C ####################################################################

  call c_info_get_I4(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    value, &
    localrc, &
    local_default_ptr, &
    local_idx_ptr, &
    recursive)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

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

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoGetI8()"
subroutine ESMF_InfoGetI8(info, key, value, keywordEnforcer, default, idx, attnestflag, rc)
  type(ESMF_Info), intent(in) :: info
  character(len=*), intent(in) :: key
  integer(ESMF_KIND_I8), intent(out) :: value
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  integer(ESMF_KIND_I8), intent(in), optional :: default
  integer, intent(in), optional :: idx
  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
  integer, intent(out), optional :: rc

  integer :: localrc
  integer(C_LONG), target :: local_default
  integer(C_INT), target :: local_idx
  type(C_PTR) :: local_default_ptr, local_idx_ptr
  integer(C_INT) :: recursive, strlen_only

  ! Set up local return code
  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE
  recursive = 0 !false
  strlen_only = 0 !false

  ! Handle optional arguments for C ###########################################

  if (present(default)) then
    local_default = default
    local_default_ptr = C_LOC(local_default)
  else
    local_default_ptr = C_NULL_PTR
  end if
  if (present(idx)) then
    local_idx = idx - 1  ! Shift to C (zero-based) indexing
    local_idx_ptr = C_LOC(local_idx)
  else
    local_idx_ptr = C_NULL_PTR
  end if
  if (present(attnestflag)) then
    if (attnestflag%value==ESMF_ATTNEST_ON%value) recursive = 1 !true
  end if

  ! Call C ####################################################################

  call c_info_get_I8(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    value, &
    localrc, &
    local_default_ptr, &
    local_idx_ptr, &
    recursive)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

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

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoGetCH()"
subroutine ESMF_InfoGetCH(info, key, value, keywordEnforcer, default, idx, attnestflag, rc)
  type(ESMF_Info), intent(in) :: info
  character(len=*), intent(in) :: key
  character(len=*), intent(out) :: value
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  character(len=*), intent(in), optional :: default
  integer, intent(in), optional :: idx
  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
  integer, intent(out), optional :: rc

  integer :: localrc, vlen
  character(:), allocatable, target :: local_default
  integer(C_INT), target :: local_idx
  type(C_PTR) :: local_default_ptr, local_idx_ptr
  integer(C_INT) :: recursive, strlen_only

  ! Set up local return code
  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE
  recursive = 0 !false
  strlen_only = 0 !false

  ! Handle optional arguments for C ###########################################

  if (present(default)) then
    local_default = trim(default)//C_NULL_CHAR
    local_default_ptr = C_LOC(local_default)
  else
    local_default_ptr = C_NULL_PTR
  end if
  if (present(idx)) then
    local_idx = idx - 1  ! Shift to C (zero-based) indexing
    local_idx_ptr = C_LOC(local_idx)
  else
    local_idx_ptr = C_NULL_PTR
  end if
  if (present(attnestflag)) then
    if (attnestflag%value==ESMF_ATTNEST_ON%value) recursive = 1 !true
  end if

  ! Call C ####################################################################

    vlen = LEN(value)
  call c_info_get_CH(info%ptr, trim(key)//C_NULL_CHAR, value, vlen, &
    localrc, local_default_ptr, local_idx_ptr, recursive, strlen_only)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

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

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoGetCharAlloc()"
subroutine ESMF_InfoGetCharAlloc(info, key, value, keywordEnforcer, default, idx, attnestflag, rc)
  type(ESMF_Info), intent(in) :: info
  character(len=*), intent(in) :: key
  character(:), allocatable, intent(out) :: value
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  character(len=*), intent(in), optional :: default
  integer, intent(in), optional :: idx
  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
  integer, intent(out), optional :: rc

  integer :: localrc, vlen
  character(:), allocatable, target :: local_default
  integer(C_INT), target :: local_idx
  type(C_PTR) :: local_default_ptr, local_idx_ptr
  integer(C_INT) :: recursive, strlen_only

  ! Set up local return code
  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE
  recursive = 0 !false
  strlen_only = 0 !false

  ! Handle optional arguments for C ###########################################

  if (present(default)) then
    local_default = trim(default)//C_NULL_CHAR
    local_default_ptr = C_LOC(local_default)
  else
    local_default_ptr = C_NULL_PTR
  end if
  if (present(idx)) then
    local_idx = idx - 1  ! Shift to C (zero-based) indexing
    local_idx_ptr = C_LOC(local_idx)
  else
    local_idx_ptr = C_NULL_PTR
  end if
  if (present(attnestflag)) then
    if (attnestflag%value==ESMF_ATTNEST_ON%value) recursive = 1 !true
  end if

  ! Call C ####################################################################

    strlen_only = 1 !true
    call c_info_get_CH(info%ptr, trim(key)//C_NULL_CHAR, value, vlen, &
      localrc, local_default_ptr, local_idx_ptr, recursive, strlen_only)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return
    strlen_only = 0 !false
    allocate(character(len=vlen) :: value)
  call c_info_get_CH(info%ptr, trim(key)//C_NULL_CHAR, value, vlen, &
    localrc, local_default_ptr, local_idx_ptr, recursive, strlen_only)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

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

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoGetLG()"
subroutine ESMF_InfoGetLG(info, key, value, keywordEnforcer, default, idx, attnestflag, rc)
  type(ESMF_Info), intent(in) :: info
  character(len=*), intent(in) :: key
  logical, intent(inout) :: value
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  logical, intent(in), optional :: default
  integer, intent(in), optional :: idx
  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
  integer, intent(out), optional :: rc

  integer :: localrc
  logical(C_BOOL), target :: local_default
  logical(C_BOOL) :: local_value
  integer(C_INT), target :: local_idx
  type(C_PTR) :: local_default_ptr, local_idx_ptr
  integer(C_INT) :: recursive, strlen_only

  ! Set up local return code
  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE
  recursive = 0 !false
  strlen_only = 0 !false

  ! Handle optional arguments for C ###########################################

  if (present(default)) then
    local_default = default
    local_default_ptr = C_LOC(local_default)
  else
    local_default_ptr = C_NULL_PTR
  end if
  if (present(idx)) then
    local_idx = idx - 1  ! Shift to C (zero-based) indexing
    local_idx_ptr = C_LOC(local_idx)
  else
    local_idx_ptr = C_NULL_PTR
  end if
  if (present(attnestflag)) then
    if (attnestflag%value==ESMF_ATTNEST_ON%value) recursive = 1 !true
  end if

  ! Call C ####################################################################

  call c_info_get_LG(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    local_value, &
    localrc, &
    local_default_ptr, &
    local_idx_ptr, &
    recursive)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  value = local_value

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

!------------------------------------------------------------------------------
! GetArray --------------------------------------------------------------------
!------------------------------------------------------------------------------

!BOP
! !IROUTINE: ESMF_InfoGet - Get a list
!
! !INTERFACE:
!subroutine ESMF_InfoGet(info, key, values, keywordEnforcer, itemCount, attnestflag, scalarToArray, rc)
! !ARGUMENTS:
!  type(ESMF_Info), intent(in) :: info
!  character(len=*), intent(in) :: key
!  <values>, see below for supported values
!type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
!  integer, intent(out), optional :: itemCount
!  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
!  logical, intent(in), optional :: scalarToArray
!  integer, intent(out), optional :: rc
!
! !DESCRIPTION:
!     Get a value list from an \texttt{ESMF\_Info} object using a key. If the key
!     is not found, \textit{rc} will not equal \texttt{ESMF\_SUCCESS}. The returned
!     value is always a copy.
!
!     The length of \textit{values} must match its length in storage.
!
!     Overloaded \textit{values} for the following types:
!     \begin{itemize}
!       \item \texttt{integer(kind=ESMF\_KIND\_I4), dimension(:)}
!       \item \texttt{integer(kind=ESMF\_KIND\_I8), dimension(:)}
!       \item \texttt{real(kind=ESMF\_KIND\_R4), dimension(:)}
!       \item \texttt{real(kind=ESMF\_KIND\_R8), dimension(:)}
!       \item \texttt{logical, dimension(:)}
!       \item \texttt{character(:), dimension(:)}
!     \end{itemize}
!
!     The arguments are:
!     \begin{description}
!     \item [info]
!       Target \texttt{ESMF\_Info} object.
!     \item [key]
!       String key to access in \texttt{ESMF\_Info} storage. See section \ref{info_key_format}
!       for an overview of the key format.
!     \item [values]
!       The output value list associated with the key.
!     \item [{[itemCount]}]
!       The number of items in \textit{values}.
!     \item [{[attnestflag]}]
!       Default is \texttt{ESMF\_ATTNEST\_OFF}. Setting to \texttt{ESMF\_ATTNEST\_ON}
!       triggers a recursive search. The first instance of the key will be found
!       in the hierarchy.
!     \item [{[scalarToArray]}]
!       Default is false. If true, allow conversion of scalar values in storage
!       to single-valued lists.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

!BOP
! !IROUTINE: ESMF_InfoGetAlloc - Get an allocatable list
!
! !INTERFACE:
!subroutine ESMF_InfoGetAlloc(info, key, values, keywordEnforcer, itemCount, attnestflag, scalarToArray, rc)
! !ARGUMENTS:
!  type(ESMF_Info), intent(in) :: info
!  character(len=*), intent(in) :: key
!  <values>, see below for supported values
!type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
!  integer, intent(out), optional :: itemCount
!  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
!  logical, intent(in), optional :: scalarToArray
!  integer, intent(out), optional :: rc
!
! !DESCRIPTION:
!     Get a value list from an \texttt{ESMF\_Info} object using a key. If the key
!     is not found, \textit{rc} will not equal \texttt{ESMF\_SUCCESS}. The returned
!     value is always a copy.
!
!     Overloaded \textit{values} for the following types:
!     \begin{itemize}
!       \item \texttt{integer(kind=ESMF\_KIND\_I4), dimension(:), allocatable}
!       \item \texttt{integer(kind=ESMF\_KIND\_I8), dimension(:), allocatable}
!       \item \texttt{real(kind=ESMF\_KIND\_R4), dimension(:), allocatable}
!       \item \texttt{real(kind=ESMF\_KIND\_R8), dimension(:), allocatable}
!       \item \texttt{logical, dimension(:), allocatable}
!       \item \texttt{character(:), dimension(:), allocatable}
!     \end{itemize}
!
!     The arguments are:
!     \begin{description}
!     \item [info]
!       Target \texttt{ESMF\_Info} object.
!     \item [key]
!       String key to access in \texttt{ESMF\_Info} storage. See section \ref{info_key_format}
!       for an overview of the key format.
!     \item [values]
!       The output value list associated with the key.
!     \item [{[itemCount]}]
!       The number of items in \textit{values}.
!     \item [{[attnestflag]}]
!       Default is \texttt{ESMF\_ATTNEST\_OFF}. Setting to \texttt{ESMF\_ATTNEST\_ON}
!       triggers a recursive search. The first instance of the key will be found
!       in the hierarchy.
!     \item [{[scalarToArray]}]
!       Default is false. If true, allow conversion of scalar values in storage
!       to single-valued lists.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! NOTE: Documentation stubs located above for generic interface compliance.
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoGetArrayR4()"
subroutine ESMF_InfoGetArrayR4(info, key, values, keywordEnforcer, itemCount, attnestflag, scalarToArray, rc)
  type(ESMF_Info), intent(in) :: info
  character(len=*), intent(in) :: key
  real(ESMF_KIND_R4), dimension(:), intent(out) :: values
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  integer, intent(out), optional :: itemCount
  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
  logical, intent(in), optional :: scalarToArray
  integer, intent(out), optional :: rc

  integer :: localrc
  integer(C_INT) :: recursive, local_itemCount, expected_size, local_scalarToArray

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE
  recursive = 0 !false

  if (present(attnestflag)) then
    if (attnestflag%value==ESMF_ATTNEST_ON%value) recursive = 1 !true
  end if

  local_scalarToArray = 0 !false
  if (present(scalarToArray)) then
    if (scalarToArray) local_scalarToArray = 1 !true
  end if

  ! This is used to validate the size of the outgoing array given the size in
  ! storage.
  expected_size = size(values)

  call c_info_get_array_R4(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    values, &
    local_itemCount, &
    localrc, &
    recursive, &
    local_scalarToArray, &
    expected_size)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (present(itemCount)) itemCount = local_itemCount
  if (present(rc)) rc = ESMF_SUCCESS
end subroutine ESMF_InfoGetArrayR4

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoGetArrayR4Alloc()"
subroutine ESMF_InfoGetArrayR4Alloc(info, key, values, keywordEnforcer, itemCount, attnestflag, scalarToArray, rc)
  type(ESMF_Info), intent(in) :: info
  character(len=*), intent(in) :: key
  real(ESMF_KIND_R4), dimension(:),  allocatable, intent(out) :: values
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  integer, intent(out), optional :: itemCount
  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
  logical, intent(in), optional :: scalarToArray
  integer, intent(out), optional :: rc

  integer :: localrc
  integer(C_INT) :: recursive, local_itemCount, expected_size, local_scalarToArray
  logical :: is_array

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE
  recursive = 0 !false

  if (present(attnestflag)) then
    if (attnestflag%value==ESMF_ATTNEST_ON%value) recursive = 1 !true
  end if

  local_scalarToArray = 0 !false
  if (present(scalarToArray)) then
    if (scalarToArray) local_scalarToArray = 1 !true
  end if

  ! Set the value to a negative one to indicate we are doing an allocatable call
  ! into the storage layer, and the size should not be checked.
  expected_size = -1

  ! Get the array size from the info store
  call ESMF_InfoGetArrayMeta(info, key, is_array, local_itemcount, attnestflag=attnestflag, &
    rc=localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  ! Allocate the outgoing storage array and call into C to fill the array
  allocate(values(local_itemCount))

  call c_info_get_array_R4(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    values, &
    local_itemCount, &
    localrc, &
    recursive, &
    local_scalarToArray, &
    expected_size)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (present(itemCount)) itemCount = local_itemCount
  if (present(rc)) rc = ESMF_SUCCESS
end subroutine ESMF_InfoGetArrayR4Alloc

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoGetArrayR8()"
subroutine ESMF_InfoGetArrayR8(info, key, values, keywordEnforcer, itemCount, attnestflag, scalarToArray, rc)
  type(ESMF_Info), intent(in) :: info
  character(len=*), intent(in) :: key
  real(ESMF_KIND_R8), dimension(:), intent(out) :: values
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  integer, intent(out), optional :: itemCount
  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
  logical, intent(in), optional :: scalarToArray
  integer, intent(out), optional :: rc

  integer :: localrc
  integer(C_INT) :: recursive, local_itemCount, expected_size, local_scalarToArray

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE
  recursive = 0 !false

  if (present(attnestflag)) then
    if (attnestflag%value==ESMF_ATTNEST_ON%value) recursive = 1 !true
  end if

  local_scalarToArray = 0 !false
  if (present(scalarToArray)) then
    if (scalarToArray) local_scalarToArray = 1 !true
  end if

  ! This is used to validate the size of the outgoing array given the size in
  ! storage.
  expected_size = size(values)

  call c_info_get_array_R8(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    values, &
    local_itemCount, &
    localrc, &
    recursive, &
    local_scalarToArray, &
    expected_size)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (present(itemCount)) itemCount = local_itemCount
  if (present(rc)) rc = ESMF_SUCCESS
end subroutine ESMF_InfoGetArrayR8

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoGetArrayR8Alloc()"
subroutine ESMF_InfoGetArrayR8Alloc(info, key, values, keywordEnforcer, itemCount, attnestflag, scalarToArray, rc)
  type(ESMF_Info), intent(in) :: info
  character(len=*), intent(in) :: key
  real(ESMF_KIND_R8), dimension(:),  allocatable, intent(out) :: values
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  integer, intent(out), optional :: itemCount
  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
  logical, intent(in), optional :: scalarToArray
  integer, intent(out), optional :: rc

  integer :: localrc
  integer(C_INT) :: recursive, local_itemCount, expected_size, local_scalarToArray
  logical :: is_array

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE
  recursive = 0 !false

  if (present(attnestflag)) then
    if (attnestflag%value==ESMF_ATTNEST_ON%value) recursive = 1 !true
  end if

  local_scalarToArray = 0 !false
  if (present(scalarToArray)) then
    if (scalarToArray) local_scalarToArray = 1 !true
  end if

  ! Set the value to a negative one to indicate we are doing an allocatable call
  ! into the storage layer, and the size should not be checked.
  expected_size = -1

  ! Get the array size from the info store
  call ESMF_InfoGetArrayMeta(info, key, is_array, local_itemcount, attnestflag=attnestflag, &
    rc=localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  ! Allocate the outgoing storage array and call into C to fill the array
  allocate(values(local_itemCount))

  call c_info_get_array_R8(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    values, &
    local_itemCount, &
    localrc, &
    recursive, &
    local_scalarToArray, &
    expected_size)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (present(itemCount)) itemCount = local_itemCount
  if (present(rc)) rc = ESMF_SUCCESS
end subroutine ESMF_InfoGetArrayR8Alloc

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoGetArrayI4()"
subroutine ESMF_InfoGetArrayI4(info, key, values, keywordEnforcer, itemCount, attnestflag, scalarToArray, rc)
  type(ESMF_Info), intent(in) :: info
  character(len=*), intent(in) :: key
  integer(ESMF_KIND_I4), dimension(:), intent(out) :: values
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  integer, intent(out), optional :: itemCount
  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
  logical, intent(in), optional :: scalarToArray
  integer, intent(out), optional :: rc

  integer :: localrc
  integer(C_INT) :: recursive, local_itemCount, expected_size, local_scalarToArray

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE
  recursive = 0 !false

  if (present(attnestflag)) then
    if (attnestflag%value==ESMF_ATTNEST_ON%value) recursive = 1 !true
  end if

  local_scalarToArray = 0 !false
  if (present(scalarToArray)) then
    if (scalarToArray) local_scalarToArray = 1 !true
  end if

  ! This is used to validate the size of the outgoing array given the size in
  ! storage.
  expected_size = size(values)

  call c_info_get_array_I4(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    values, &
    local_itemCount, &
    localrc, &
    recursive, &
    local_scalarToArray, &
    expected_size)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (present(itemCount)) itemCount = local_itemCount
  if (present(rc)) rc = ESMF_SUCCESS
end subroutine ESMF_InfoGetArrayI4

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoGetArrayI4Alloc()"
subroutine ESMF_InfoGetArrayI4Alloc(info, key, values, keywordEnforcer, itemCount, attnestflag, scalarToArray, rc)
  type(ESMF_Info), intent(in) :: info
  character(len=*), intent(in) :: key
  integer(ESMF_KIND_I4), dimension(:),  allocatable, intent(out) :: values
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  integer, intent(out), optional :: itemCount
  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
  logical, intent(in), optional :: scalarToArray
  integer, intent(out), optional :: rc

  integer :: localrc
  integer(C_INT) :: recursive, local_itemCount, expected_size, local_scalarToArray
  logical :: is_array

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE
  recursive = 0 !false

  if (present(attnestflag)) then
    if (attnestflag%value==ESMF_ATTNEST_ON%value) recursive = 1 !true
  end if

  local_scalarToArray = 0 !false
  if (present(scalarToArray)) then
    if (scalarToArray) local_scalarToArray = 1 !true
  end if

  ! Set the value to a negative one to indicate we are doing an allocatable call
  ! into the storage layer, and the size should not be checked.
  expected_size = -1

  ! Get the array size from the info store
  call ESMF_InfoGetArrayMeta(info, key, is_array, local_itemcount, attnestflag=attnestflag, &
    rc=localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  ! Allocate the outgoing storage array and call into C to fill the array
  allocate(values(local_itemCount))

  call c_info_get_array_I4(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    values, &
    local_itemCount, &
    localrc, &
    recursive, &
    local_scalarToArray, &
    expected_size)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (present(itemCount)) itemCount = local_itemCount
  if (present(rc)) rc = ESMF_SUCCESS
end subroutine ESMF_InfoGetArrayI4Alloc

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoGetArrayI8()"
subroutine ESMF_InfoGetArrayI8(info, key, values, keywordEnforcer, itemCount, attnestflag, scalarToArray, rc)
  type(ESMF_Info), intent(in) :: info
  character(len=*), intent(in) :: key
  integer(ESMF_KIND_I8), dimension(:), intent(out) :: values
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  integer, intent(out), optional :: itemCount
  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
  logical, intent(in), optional :: scalarToArray
  integer, intent(out), optional :: rc

  integer :: localrc
  integer(C_INT) :: recursive, local_itemCount, expected_size, local_scalarToArray

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE
  recursive = 0 !false

  if (present(attnestflag)) then
    if (attnestflag%value==ESMF_ATTNEST_ON%value) recursive = 1 !true
  end if

  local_scalarToArray = 0 !false
  if (present(scalarToArray)) then
    if (scalarToArray) local_scalarToArray = 1 !true
  end if

  ! This is used to validate the size of the outgoing array given the size in
  ! storage.
  expected_size = size(values)

  call c_info_get_array_I8(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    values, &
    local_itemCount, &
    localrc, &
    recursive, &
    local_scalarToArray, &
    expected_size)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (present(itemCount)) itemCount = local_itemCount
  if (present(rc)) rc = ESMF_SUCCESS
end subroutine ESMF_InfoGetArrayI8

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoGetArrayI8Alloc()"
subroutine ESMF_InfoGetArrayI8Alloc(info, key, values, keywordEnforcer, itemCount, attnestflag, scalarToArray, rc)
  type(ESMF_Info), intent(in) :: info
  character(len=*), intent(in) :: key
  integer(ESMF_KIND_I8), dimension(:),  allocatable, intent(out) :: values
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  integer, intent(out), optional :: itemCount
  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
  logical, intent(in), optional :: scalarToArray
  integer, intent(out), optional :: rc

  integer :: localrc
  integer(C_INT) :: recursive, local_itemCount, expected_size, local_scalarToArray
  logical :: is_array

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE
  recursive = 0 !false

  if (present(attnestflag)) then
    if (attnestflag%value==ESMF_ATTNEST_ON%value) recursive = 1 !true
  end if

  local_scalarToArray = 0 !false
  if (present(scalarToArray)) then
    if (scalarToArray) local_scalarToArray = 1 !true
  end if

  ! Set the value to a negative one to indicate we are doing an allocatable call
  ! into the storage layer, and the size should not be checked.
  expected_size = -1

  ! Get the array size from the info store
  call ESMF_InfoGetArrayMeta(info, key, is_array, local_itemcount, attnestflag=attnestflag, &
    rc=localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  ! Allocate the outgoing storage array and call into C to fill the array
  allocate(values(local_itemCount))

  call c_info_get_array_I8(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    values, &
    local_itemCount, &
    localrc, &
    recursive, &
    local_scalarToArray, &
    expected_size)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (present(itemCount)) itemCount = local_itemCount
  if (present(rc)) rc = ESMF_SUCCESS
end subroutine ESMF_InfoGetArrayI8Alloc

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoGetArrayCH()"
subroutine ESMF_InfoGetArrayCH(info, key, values, keywordEnforcer, itemCount, attnestflag, scalarToArray, rc)
  type(ESMF_Info), intent(in) :: info
  character(len=*), intent(in) :: key
  character(len=*), dimension(:), intent(out) :: values
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  integer, intent(out), optional :: itemCount
  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
  logical, intent(in), optional :: scalarToArray
  integer, intent(out), optional :: rc

  integer :: localrc
  integer(C_INT) :: recursive, local_itemCount, expected_size, local_scalarToArray
  integer :: ii
  logical :: is_array

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE
  recursive = 0 !false

  if (present(attnestflag)) then
    if (attnestflag%value==ESMF_ATTNEST_ON%value) recursive = 1 !true
  end if

  local_scalarToArray = 0 !false
  if (present(scalarToArray)) then
    if (scalarToArray) local_scalarToArray = 1 !true
  end if

  ! Set the value to a negative one to indicate we are doing an allocatable call
  ! into the storage layer, and the size should not be checked.
  expected_size = -1

  ! Get the array size from the info store
  call ESMF_InfoGetArrayMeta(info, key, is_array, local_itemcount, attnestflag=attnestflag, &
    rc=localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (.not. is_array .and. local_scalarToArray == 0) then
    if (ESMF_LogFoundError(ESMF_RC_ATTR_WRONGTYPE, &
      msg="Array requested but type in JSON storage is not an array. Key is: "//TRIM(key), &
      ESMF_CONTEXT, rcToReturn=rc)) return
  end if

  if (local_itemCount /= SIZE(values)) then
    if (ESMF_LogFoundError(ESMF_RC_ATTR_ITEMSOFF, msg="values allocation size does not match size in Info storage", ESMF_CONTEXT, rcToReturn=rc)) return
  end if
  do ii=1,local_itemCount
    if (.not. is_array .and. local_scalarToArray == 1) then
      call ESMF_InfoGetCH(info, key, values(ii), attnestflag=attnestflag, rc=localrc)
    else
      call ESMF_InfoGetCH(info, key, values(ii), idx=ii, attnestflag=attnestflag, rc=localrc)
    end if
  enddo
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (present(itemCount)) itemCount = local_itemCount
  if (present(rc)) rc = ESMF_SUCCESS
end subroutine ESMF_InfoGetArrayCH

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoGetArrayCHAlloc()"
subroutine ESMF_InfoGetArrayCHAlloc(info, key, values, keywordEnforcer, itemCount, attnestflag, scalarToArray, rc)
  type(ESMF_Info), intent(in) :: info
  character(len=*), intent(in) :: key
  character(len=*), dimension(:),  allocatable, intent(out) :: values
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  integer, intent(out), optional :: itemCount
  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
  logical, intent(in), optional :: scalarToArray
  integer, intent(out), optional :: rc

  integer :: localrc
  integer(C_INT) :: recursive, local_itemCount, expected_size, local_scalarToArray
  integer :: ii
  logical :: is_array

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE
  recursive = 0 !false

  if (present(attnestflag)) then
    if (attnestflag%value==ESMF_ATTNEST_ON%value) recursive = 1 !true
  end if

  local_scalarToArray = 0 !false
  if (present(scalarToArray)) then
    if (scalarToArray) local_scalarToArray = 1 !true
  end if

  ! Set the value to a negative one to indicate we are doing an allocatable call
  ! into the storage layer, and the size should not be checked.
  expected_size = -1

  ! Get the array size from the info store
  call ESMF_InfoGetArrayMeta(info, key, is_array, local_itemcount, attnestflag=attnestflag, &
    rc=localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (.not. is_array .and. local_scalarToArray == 0) then
    if (ESMF_LogFoundError(ESMF_RC_ATTR_WRONGTYPE, &
      msg="Array requested but type in JSON storage is not an array. Key is: "//TRIM(key), &
      ESMF_CONTEXT, rcToReturn=rc)) return
  end if

  ! Allocate the outgoing storage array and call into C to fill the array
  allocate(values(local_itemCount))
  do ii=1,local_itemCount
    if (.not. is_array .and. local_scalarToArray == 1) then
      call ESMF_InfoGetCH(info, key, values(ii), attnestflag=attnestflag, rc=localrc)
    else
      call ESMF_InfoGetCH(info, key, values(ii), idx=ii, attnestflag=attnestflag, rc=localrc)
    end if
  enddo
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (present(itemCount)) itemCount = local_itemCount
  if (present(rc)) rc = ESMF_SUCCESS
end subroutine ESMF_InfoGetArrayCHAlloc

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoGetArrayLG()"
subroutine ESMF_InfoGetArrayLG(info, key, values, keywordEnforcer, itemCount, attnestflag, scalarToArray, rc)
  type(ESMF_Info), intent(in) :: info
  character(len=*), intent(in) :: key
  logical, dimension(:), intent(out) :: values
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  integer, intent(out), optional :: itemCount
  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
  logical, intent(in), optional :: scalarToArray
  integer, intent(out), optional :: rc

  integer :: localrc
  integer(C_INT) :: recursive, local_itemCount, expected_size, local_scalarToArray
  integer :: ii
  logical(C_BOOL), dimension(:), allocatable :: local_values

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE
  recursive = 0 !false

  if (present(attnestflag)) then
    if (attnestflag%value==ESMF_ATTNEST_ON%value) recursive = 1 !true
  end if

  local_scalarToArray = 0 !false
  if (present(scalarToArray)) then
    if (scalarToArray) local_scalarToArray = 1 !true
  end if

  ! This is used to validate the size of the outgoing array given the size in
  ! storage.
  expected_size = size(values)

  ! Store boolean integers for transfer from C
  allocate(local_values(SIZE(values)))

  call c_info_get_array_LG(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    local_values, &
    local_itemCount, &
    localrc, &
    recursive, &
    local_scalarToArray, &
    expected_size)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  ! Transfer to logical storage from boolean integers
  do ii=1,SIZE(values)
    values(ii) = local_values(ii)
  enddo
  deallocate(local_values)

  if (present(itemCount)) itemCount = local_itemCount
  if (present(rc)) rc = ESMF_SUCCESS
end subroutine ESMF_InfoGetArrayLG

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoGetArrayLGAlloc()"
subroutine ESMF_InfoGetArrayLGAlloc(info, key, values, keywordEnforcer, itemCount, attnestflag, scalarToArray, rc)
  type(ESMF_Info), intent(in) :: info
  character(len=*), intent(in) :: key
  logical, dimension(:),  allocatable, intent(out) :: values
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  integer, intent(out), optional :: itemCount
  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
  logical, intent(in), optional :: scalarToArray
  integer, intent(out), optional :: rc

  integer :: localrc
  integer(C_INT) :: recursive, local_itemCount, expected_size, local_scalarToArray
  integer :: ii
  logical :: is_array
  logical(C_BOOL), dimension(:), allocatable :: local_values

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE
  recursive = 0 !false

  if (present(attnestflag)) then
    if (attnestflag%value==ESMF_ATTNEST_ON%value) recursive = 1 !true
  end if

  local_scalarToArray = 0 !false
  if (present(scalarToArray)) then
    if (scalarToArray) local_scalarToArray = 1 !true
  end if

  ! Set the value to a negative one to indicate we are doing an allocatable call
  ! into the storage layer, and the size should not be checked.
  expected_size = -1

  ! Get the array size from the info store
  call ESMF_InfoGetArrayMeta(info, key, is_array, local_itemcount, attnestflag=attnestflag, &
    rc=localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  ! Allocate the outgoing storage array and call into C to fill the array
  allocate(values(local_itemCount))

  ! Store boolean integers for transfer from C
  allocate(local_values(SIZE(values)))

  call c_info_get_array_LG(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    local_values, &
    local_itemCount, &
    localrc, &
    recursive, &
    local_scalarToArray, &
    expected_size)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  ! Transfer to logical storage from boolean integers
  do ii=1,SIZE(values)
    values(ii) = local_values(ii)
  enddo
  deallocate(local_values)

  if (present(itemCount)) itemCount = local_itemCount
  if (present(rc)) rc = ESMF_SUCCESS
end subroutine ESMF_InfoGetArrayLGAlloc

!------------------------------------------------------------------------------

!BOP
! !IROUTINE: ESMF_InfoGet - Inquire an Info object for metadata
!
! !INTERFACE:
  ! Private name; call using ESMF_InfoGet()
!subroutine ESMF_InfoInquire(info, keywordEnforcer, size, key, jsonType, isArray, &
!  isDirty, idx, typekind, ikey, isPresent, isStructured, isNull, rc)
! !ARGUMENTS:
!  type(ESMF_Info), intent(in) :: info
!type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
!  integer, intent(out), optional :: size
!  character(len=*), intent(in), optional :: key
!  character(len=*), intent(out), optional :: jsonType
!  logical, intent(out), optional :: isArray
!  logical, intent(out), optional :: isDirty
!  integer, intent(in), optional :: idx
!  type(ESMF_TypeKind_Flag), intent(out), optional :: typekind
!  character(len=*), intent(out), optional :: ikey
!  logical, intent(out), optional :: isPresent
!  logical, intent(out), optional :: isStructured
!  logical, intent(out), optional :: isNull
!  integer, intent(out), optional :: rc
!
! !DESCRIPTION:
!     Inquire an \texttt{ESMF\_Info} object for metadata.
!
!     The arguments are:
!     \begin{description}
!     \item [info]
!       Target \texttt{ESMF\_Info} object.
!     \item [{[size]}]
!       Returns the size of the target. The following rules apply:
!       \begin{itemize}
!         \item If the target is an object, return the number of key-value pairs.
!         \item If the target is a scalar, return \texttt{1}.
!         \item If the target is an array, return its length.
!       \end{itemize}
!     \item [{[key]}]
!       If provided, use this location as the origin instead of root. See section
!       \ref{info_key_format} for an overview of the key format.
!     \item [{[jsonType]}]
!       Returns the JSON object type name \cite{json_for_modern_cpp_typename}.
!     \item [{[isArray]}]
!       Returns true if the target is an array.
!     \item [{[isDirty]}]
!       Returns true if the \texttt{ESMF\_Info} object should be synchronized during
!       an \texttt{ESMF\_InfoSync} operation.
!     \item [{[idx]}]
!       An integer index to use. This will index into an object type providing
!       the primary mechanism for iteration.
!     \item [{[typekind]}]
!       Get the ESMF typekind for the target. The minimum typekind required to
!       hold the value is returned.
!       See section \ref{const:typekind} for valid values.
!     \item [{[ikey]}]
!       If present, this will be set to the key's name for the current inquire.
!       Useful when iterating using an index. This does \textit{not} return the full key
!       path if nested.
!     \item [{[isPresent]}]
!       Returns true if the \textit{key} exists in storage. If no \textit{key}
!       is provided, this will return true.
!     \item [{[isStructured]}]
!       Returns true if the target is structured \cite{json_for_modern_cpp_is_structured}.
!       This means it is either an object (a map) or an array.
!     \item [{[isNull]}]
!       Returns true if the target is null \cite{json_for_modern_cpp_null}.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP
! Private parameter documentation:
!     \item [{[attrCount]}]
!       Returns the count of of \texttt{ESMF\_Attribute} objects excluding \texttt{ESMF\_AttPack}
!       conventions and purposes. Does not descend the object hierarchy when counting.
!     \item [{[attrCountTotal]}]
!       Returns the count of of \texttt{ESMF\_Attribute} objects excluding \texttt{ESMF\_AttPack}
!       conventions and purposes. Descends the object hierarchy when counting.
!     \item [{[attPackCount]}]
!       Returns the count of \texttt{ESMF\_AttPack} conventions and purposes.
!       Does not descend the object hierarchy when counting.
!     \item [{[attPackCountTotal]}]
!       Returns the count of \texttt{ESMF\_AttPack} conventions and purposes.
!       Does descend the object hierarchy when counting.
!     \item [{[attrCompliance]}]
!       Default is false. If true, treat the JSON scheme as an \textit{ESMF\_Attribute}
!       schema.

!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! NOTE: Documentation stub above to make some parameters private.
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoInquire()"
subroutine ESMF_InfoInquire(info, keywordEnforcer, size, key, attrCount, attrCountTotal, jsonType, &
  isArray, isDirty, attPackCount, attPackCountTotal, attnestflag, idx, typekind, &
  ikey, isPresent, isStructured, isNull, attrCompliance, rc)

  type(ESMF_Info), intent(in) :: info
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  integer(C_INT), intent(out), optional :: size
  character(len=*), intent(in), optional :: key
  integer(C_INT), intent(out), optional :: attrCount
  integer(C_INT), intent(out), optional :: attrCountTotal
  character(len=*), intent(out), optional :: jsonType
  logical, intent(out), optional :: isArray
  logical, intent(out), optional :: isDirty
  integer(C_INT), intent(out), optional :: attPackCount
  integer(C_INT), intent(out), optional :: attPackCountTotal
  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
  integer(C_INT), intent(in), optional :: idx
  type(ESMF_TypeKind_Flag), intent(out), optional :: typekind
  character(len=*), intent(out), optional :: ikey
  logical, intent(out), optional :: isPresent
  logical, intent(out), optional :: isStructured
  logical, intent(out), optional :: isNull
  logical, intent(in), optional :: attrCompliance
  integer, intent(out), optional :: rc

  integer :: localrc, esmc_typekind, local_size
  type(ESMF_Info) :: inq, inq2
  character(:), allocatable :: local_key
  integer :: recursive, local_attrCompliance
  integer(C_INT), target :: local_idx
  type(C_PTR) :: local_idx_ptr

  if (present(rc)) rc = ESMF_FAILURE
  localrc = ESMF_FAILURE
  recursive = 0 !false
  local_attrCompliance = 0 !false

  if (present(key)) then
    local_key = trim(key)
  else
    local_key = ""
  end if
  if (present(attnestflag)) then
    if (attnestflag%value==ESMF_ATTNEST_ON%value) recursive = 1 !true
  end if
  if (present(idx)) then
    local_idx = idx - 1  ! Shift to C (zero-based) indexing
    local_idx_ptr = C_LOC(local_idx)
  else
    local_idx_ptr = C_NULL_PTR
  end if
  if (present(attrCompliance)) then
    if (attrCompliance) local_attrCompliance = 1 !true
  end if

  if (present(isPresent)) then
    if (LEN(key) > 0) then
      isPresent = ESMF_InfoIsPresent(info, local_key, isPointer=.true., &
        attnestflag=attnestflag, rc=localrc)
      if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return
    else
      isPresent = .true.
    end if
  end if

  inq = ESMF_InfoCreate(rc=localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  call c_info_inquire(info%ptr, inq%ptr, local_key//C_NULL_CHAR, recursive, &
   local_idx_ptr, local_attrCompliance, localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (present(size)) then
      call ESMF_InfoGet(inq, "size", size, rc=localrc)
      if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return
  end if
  if (present(attrCount)) then
    call ESMF_InfoGet(inq, "attrCount", attrCount, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return
  end if
  if (present(attrCountTotal)) then
      call ESMF_InfoGet(inq, "attrCountTotal", attrCountTotal, rc=localrc)
      if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return
  end if
  if (present(jsonType)) then
    call ESMF_InfoGet(inq, "jsonType", jsonType, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return
  end if
  if (present(isArray)) then
    call ESMF_InfoGet(inq, "isArray", isArray, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return
  end if
  if (present(isDirty)) then
    call ESMF_InfoGet(inq, "isDirty", isDirty, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return
  end if
  if (present(attPackCount)) then
    call ESMF_InfoGet(inq, "attPackCount", attPackCount, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return
  end if
  if (present(attPackCountTotal)) then
      call ESMF_InfoGet(inq, "attPackCountTotal", attPackCountTotal, rc=localrc)
      if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return
  end if
  if (present(typekind)) then
    call ESMF_InfoGet(inq, "ESMC_TypeKind_Flag", esmc_typekind, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return
    typekind = ESMF_TypeKind_Flag(esmc_typekind)
  end if
  if (present(ikey)) then
    call ESMF_InfoGet(inq, "key", ikey, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return
  end if
  if (present(idx)) then
    if (present(isPresent)) then
      if (isPresent) then
          inq2 = ESMF_InfoCreate(rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

          call c_info_inquire(info%ptr, inq2%ptr, local_key//C_NULL_CHAR, recursive, &
            C_NULL_PTR, local_attrCompliance, localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

          call ESMF_InfoGet(inq2, "size", local_size, rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

          if (idx > local_size) then
            isPresent = .false.
          end if

          call ESMF_InfoDestroy(inq2, rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return
       end if
    end if
  end if
  if (present(isStructured)) then
    call ESMF_InfoGet(inq, "isStructured", isStructured, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return
  end if
  if (present(isNull)) then
    call ESMF_InfoGet(inq, "isNull", isNull, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return
  end if

  call ESMF_InfoDestroy(inq, rc=rc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

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

!------------------------------------------------------------------------------
! GetHandle Documentation -----------------------------------------------------
!------------------------------------------------------------------------------

!BOP
! !IROUTINE: ESMF_InfoGetFromHost - Get an Info handle from an ESMF object
!
! !INTERFACE:
!subroutine ESMF_InfoGetFromHost(host, info, keywordEnforcer, rc)
! !ARGUMENTS:
!  type(ESMF_*), intent(inout) :: host
!  type(ESMF_Info), intent(out) :: info
!type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
!  integer, intent(out), optional :: rc
!
! !DESCRIPTION:
!     Get an \texttt{ESMF\_Info} object handle from a host ESMF object. The returned
!     handle should not be destroyed.
!
!     The arguments are:
!     \begin{description}
!     \item [host]
!       Target ESMF object. Overloaded for:
!       \begin{itemize}
!         \item \texttt{ESMF\_Array}
!         \item \texttt{ESMF\_ArrayBundle}
!         \item \texttt{ESMF\_CplComp}
!         \item \texttt{ESMF\_GridComp}
!         \item \texttt{ESMF\_SciComp}
!         \item \texttt{ESMF\_DistGrid}
!         \item \texttt{ESMF\_Field}
!         \item \texttt{ESMF\_FieldBundle}
!         \item \texttt{ESMF\_Grid}
!         \item \texttt{ESMF\_State}
!         \item \texttt{ESMF\_LocStream}
!         \item \texttt{ESMF\_Mesh}
!       \end{itemize}
!     \item [info]
!       Outgoing \texttt{ESMF\_Info} object.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

! =============================================================================

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoGetFromBase()"
!BOPI
! !IROUTINE: ESMF_InfoGetFromBase
!
! !INTERFACE:
subroutine ESMF_InfoGetFromBase(base, info, keywordEnforcer, rc)
! !ARGUMENTS:
  type(ESMF_Base), intent(in) :: base
  type(ESMF_Info), intent(out) :: info
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  integer, intent(out), optional :: rc
!
! !DESCRIPTION:
!     Get an \texttt{ESMF\_Info} handle from an \texttt{ESMF\_Base} object.
!
!     The arguments are:
!     \begin{description}
!     \item [base]
!       Target \texttt{ESMF\_Base} object.
!     \item [info]
!       Outgoing \texttt{ESMF\_Info} object.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOPI

  if (present(rc)) rc = ESMF_FAILURE
  info%ptr = c_info_base_get(base%this%ptr)
  if (present(rc)) rc = ESMF_SUCCESS
end subroutine ESMF_InfoGetFromBase

! -----------------------------------------------------------------------------

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoGetFromPointer()"
!BOPI
! !IROUTINE: ESMF_InfoGetFromPointer
!
! !INTERFACE:
subroutine ESMF_InfoGetFromPointer(ptr, info, keywordEnforcer, rc)
! !ARGUMENTS:
  type(ESMF_Pointer), intent(in) :: ptr
  type(ESMF_Info), intent(out) :: info
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  integer, intent(out), optional :: rc
!
! !DESCRIPTION:
!     Get an \texttt{ESMF\_Info} handle from an \texttt{ESMF\Pointer} object.
!
!     The arguments are:
!     \begin{description}
!     \item [base]
!       Target \texttt{ESMF\_Pointer} object.
!     \item [info]
!       Outgoing \texttt{ESMF\_Info} object.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOPI

  if (present(rc)) rc = ESMF_FAILURE
  info%ptr = c_info_base_get(ptr%ptr)
  if (present(rc)) rc = ESMF_SUCCESS
end subroutine ESMF_InfoGetFromPointer

!------------------------------------------------------------------------------

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoGetTK()"
!BOP
! !IROUTINE: ESMF_InfoGetTK - Retrieve the ESMF TypeKind for a key
!
! !INTERFACE:
function ESMF_InfoGetTK(info, key, keywordEnforcer, attnestflag, rc) result(typekind)
! !ARGUMENTS:
  type(ESMF_Info), intent(in) :: info
  character(len=*), intent(in) :: key
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
  integer, intent(out), optional :: rc
! !RETURN VALUE:
  type(ESMF_TypeKind_Flag) :: typekind
!
! !DESCRIPTION:
!     Return the ESMF TypeKind of the value associated with \textit{key}.
!     See section \ref{const:typekind} for valid return values.
!
!     The arguments are:
!     \begin{description}
!     \item [info]
!       Target \texttt{ESMF\_Info} object.
!     \item [key]
!       String key to access in \texttt{ESMF\_Info} storage. See section \ref{info_key_format}
!       for an overview of the key format.
!     \item [{[attnestflag]}]
!       Setting to \texttt{ESMF\_ATTNEST\_ON} triggers a recursive search for
!       \textit{keyParent}. The first instance of the key will be found in the
!       hierarchy. Default is \texttt{ESMF\_ATTNEST\_OFF}.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

  integer :: localrc
  integer(C_INT) :: local_typekind
  integer(C_INT) :: recursive

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE
  recursive = 0 !false
  local_typekind = 0

  if (present(attnestflag)) then
    if (attnestflag%value==ESMF_ATTNEST_ON%value) recursive = 1 !true
  end if

  call c_info_get_tk(info%ptr, trim(key)//C_NULL_CHAR, local_typekind, &
    localrc, recursive)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  typekind = ESMF_TypeKind_Flag(local_typekind)

  if (present(rc)) rc = ESMF_SUCCESS
end function ESMF_InfoGetTK

!------------------------------------------------------------------------------

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoGetArrayMeta()"
!BOP
! !IROUTINE: ESMF_InfoGetArrayMeta - Retrieve array metadata information
!
! !INTERFACE:
subroutine ESMF_InfoGetArrayMeta(info, key, isArray, size, keywordEnforcer, attnestflag, rc)
! !ARGUMENTS:
  type(ESMF_Info), intent(in) :: info
  character(len=*), intent(in) :: key
  logical, intent(out) :: isArray
  integer(C_INT), intent(out) :: size
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
  integer, intent(out), optional :: rc
!
! !DESCRIPTION:
!     Return a value's array status and size using a \textit{key}.
!
!     The arguments are:
!     \begin{description}
!     \item [info]
!       Target \texttt{ESMF\_Info} object.
!     \item [key]
!       String key to access in \texttt{ESMF\_Info} storage. See section \ref{info_key_format}
!       for an overview of the key format.
!     \item [isArray]
!       Set to \texttt{true} if the target is an array in storage.
!     \item [size]
!       Set to the size of the target object in storage (i.e. length of the array).
!     \item [{[attnestflag]}]
!       Setting to \texttt{ESMF\_ATTNEST\_ON} triggers a recursive search for
!       \textit{keyParent}. The first instance of the key will be found in the
!       hierarchy. Default is \texttt{ESMF\_ATTNEST\_OFF}.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

  integer :: localrc
  integer(C_INT) :: local_is_array
  integer(C_INT) :: recursive

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE

  recursive = 0 !false
  if (present(attnestflag)) then
    if (attnestflag%value==ESMF_ATTNEST_ON%value) recursive = 1 !true
  end if

  call c_info_get_array_meta(info%ptr, trim(key)//C_NULL_CHAR, local_is_array, size, recursive, &
    localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (local_is_array == 1) then
    isArray = .true.
  else
    isArray = .false.
  end if

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

!------------------------------------------------------------------------------

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoIsPresent()"
!BOP
! !IROUTINE: ESMF_InfoIsPresent - Check for key presence
!
! !INTERFACE:
function ESMF_InfoIsPresent(info, key, keywordEnforcer, attnestflag, isPointer, rc) result(is_present)
! !ARGUMENTS:
  type(ESMF_Info), intent(in) :: info
  character(len=*), intent(in) :: key
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
  logical, intent(in), optional :: isPointer
  integer, intent(out), optional :: rc
! !RETURN VALUE:
  logical :: is_present
!
! !DESCRIPTION:
!     Return true if \textit{key} exists in \texttt{ESMF\_Info}'s storage.
!
!     The arguments are:
!     \begin{description}
!     \item [info]
!       Target \texttt{ESMF\_Info} object.
!     \item [key]
!       String key to access in \texttt{ESMF\_Info} storage. See section \ref{info_key_format}
!       for an overview of the key format.
!     \item [{[attnestflag]}]
!       Setting to \texttt{ESMF\_ATTNEST\_ON} triggers a recursive search for
!       \textit{keyParent}. The first instance of the key will be found in the
!       hierarchy. Default is \texttt{ESMF\_ATTNEST\_OFF}.
!     \item [{[isPointer]}]
!       Default is true. If true, expect the \textit{key} is using JSON Pointer
!       syntax (see section \ref{info_key_format}). Setting to false will trigger
!       a slightly faster search.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

  logical :: local_isPointer
  integer :: localrc
  integer(C_INT) :: isPointer_forC
  integer(C_INT) :: local_is_present
  integer(C_INT) :: recursive

  is_present = .false.
  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE
  recursive = 0 !false
  local_is_present = 0 !false

  if (present(attnestflag)) then
    if (attnestflag%value==ESMF_ATTNEST_ON%value) recursive = 1 !true
  end if
  if (present(isPointer)) then
    local_isPointer = isPointer
  else
    local_isPointer = .true.
  end if

  if (local_isPointer) then
    isPointer_forC = 1 !true
  else
    isPointer_forC = 0 !false
  end if

  call c_info_is_present(info%ptr, trim(key)//C_NULL_CHAR, local_is_present, &
    localrc, recursive, isPointer_forC)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (local_is_present == 1) then
    is_present = .true.
  else
    is_present = .false.
  end if

  if (present(rc)) rc = ESMF_SUCCESS
end function ESMF_InfoIsPresent

!------------------------------------------------------------------------------

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoIsSet()"
!BOP
! !IROUTINE: ESMF_InfoIsSet - Check if a value is null
!
! !INTERFACE:
function ESMF_InfoIsSet(info, key, keywordEnforcer, rc) result(is_set)
! !ARGUMENTS:
  type(ESMF_Info), intent(in) :: info
  character(len=*), intent(in) :: key
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  integer, intent(out), optional :: rc
! !RETURN VALUE:
  logical :: is_set
!
! !DESCRIPTION:
!     Returns true if the target value is not null \cite{json_for_modern_cpp_null}.
!
!     The arguments are:
!     \begin{description}
!     \item [info]
!       Target \texttt{ESMF\_Info} object.
!     \item [key]
!       String key to access in \texttt{ESMF\_Info} storage. See section \ref{info_key_format}
!       for an overview of the key format.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

  integer :: localrc
  integer(C_INT) :: is_set_c

  is_set = .false.
  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE

  call c_info_is_set(info%ptr, trim(key)//C_NULL_CHAR, is_set_c, localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (is_set_c == 1) then
    is_set = .true.
  end if

  if (present(rc)) rc = ESMF_SUCCESS
end function ESMF_InfoIsSet

!------------------------------------------------------------------------------

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoPrint()"
!BOP
! !IROUTINE: ESMF_InfoPrint - Print contents of an Info object
!
! !INTERFACE:
subroutine ESMF_InfoPrint(info, keywordEnforcer, indent, preString, unit, rc)
! !ARGUMENTS:
  type(ESMF_Info), intent(in) :: info
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  character(*), intent(in), optional :: preString
  character(*), intent(out), optional :: unit
  integer, intent(in), optional :: indent
  integer, intent(out), optional :: rc
!
! !DESCRIPTION:
!     Print \texttt{ESMF\_Info} contents in JSON format.
!
!     The arguments are:
!     \begin{description}
!     \item [info]
!       Target \texttt{ESMF\_Info} object.
!     \item [{[indent]}]
!       Default is 0. Specify a "pretty print" indentation for the JSON output string.
!     \item [{[preString]}]
!       Optionally prepended string. Default to empty string.
!     \item [{[unit]}]
!       Internal unit, i.e. a string. Default to printing to stdout.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

  integer :: localrc, local_indent
  character(:), allocatable :: output, local_preString

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE

  if (present(indent)) then
    local_indent = indent
  else
    local_indent = 4
  end if
  if (present(preString)) then
    local_preString = preString
  else
    local_preString = ""
  endif

  output = ESMF_InfoDump(info, indent=local_indent, rc=localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (present(unit)) then
    unit = local_preString//output
  else
    print *, local_preString//output
  endif

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

!------------------------------------------------------------------------------

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoReadJSON()"
!BOP
! !IROUTINE: ESMF_InfoReadJSON - Read JSON data from file
!
! !INTERFACE:
function ESMF_InfoReadJSON(filename, keywordEnforcer, rc) result(info_r)
! !ARGUMENTS:
  character(len=*), intent(in) :: filename
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  integer, intent(out), optional :: rc
! !RETURN VALUE:
  type(ESMF_Info) :: info_r
!
! !DESCRIPTION:
!     Read JSON data from a file and return a new \texttt{ESMF\_Info} object.
!
!     The arguments are:
!     \begin{description}
!     \item [filename]
!       Path to the input file.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

  integer :: localrc

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE

  info_r = ESMF_InfoCreate(rc=localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  call c_info_read_json(info_r%ptr, trim(filename)//C_NULL_CHAR, localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (present(rc)) rc = ESMF_SUCCESS
end function ESMF_InfoReadJSON

!------------------------------------------------------------------------------

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoRemove()"
!BOP
! !IROUTINE: ESMF_InfoRemove - Remove a key-value pair from an Info object
!
! !INTERFACE:
subroutine ESMF_InfoRemove(info, keyParent, keywordEnforcer, keyChild, attnestflag, rc)
! !ARGUMENTS:
  type(ESMF_Info), intent(inout) :: info
  character(len=*), intent(in) :: keyParent
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  character(len=*), intent(in), optional :: keyChild
  type(ESMF_AttNest_Flag), intent(in), optional :: attnestflag
  integer, intent(out), optional :: rc
!
! !DESCRIPTION:
!     Remove a key-value pair from an \texttt{ESMF\_Info} object.
!
!     The arguments are:
!     \begin{description}
!     \item [info]
!       Target \texttt{ESMF\_Info} object.
!     \item [keyParent]
!       String key to identify the parent location for the removal. If no \textit{keyChild}
!       is specified, then the root location is assumed. See section \ref{info_key_format}
!       for an overview of the key format.
!     \item [{[keyChild]}]
!       String key to identify the value for the removal. This \textit{may not}
!       be a path.
!     \item [{[attnestflag]}]
!       Setting to \texttt{ESMF\_ATTNEST\_ON} triggers a recursive search for
!       \textit{keyParent}. The first instance of the key will be found in the
!       hierarchy. Default is \texttt{ESMF\_ATTNEST\_OFF}.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

  integer :: localrc
  character(:), allocatable :: localkeyChild
  logical(C_BOOL) :: recursive

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE
  recursive = .false.

  if (present(keyChild)) then
    localkeyChild = keyChild
  else
    localkeyChild = ""
  end if
  if (present(attnestflag)) then
    if (attnestflag%value==ESMF_ATTNEST_ON%value) recursive = .true.
  end if

  call c_info_erase(info%ptr, trim(keyParent)//C_NULL_CHAR, &
                     trim(localkeyChild)//C_NULL_CHAR, recursive, localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

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

!------------------------------------------------------------------------------
! Set (Scalar) ----------------------------------------------------------------
!------------------------------------------------------------------------------

!BOP
! !IROUTINE: ESMF_InfoSet - Set a value
!
! !INTERFACE:
!subroutine ESMF_InfoSet(info, key, value, keywordEnforcer, force, idx, pkey, rc)
!
! !ARGUMENTS:
!  type(ESMF_Info), intent(inout) :: info
!  character(len=*), intent(in) :: key
!  <value>, see below for supported value
!type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
!  logical, intent(in), optional :: force
!  integer, intent(in), optional :: idx
!  character(len=*), intent(in), optional :: pkey
!  integer, intent(out), optional :: rc
!
! !DESCRIPTION:
!     Set a \textit{value} in an \texttt{ESMF\_Info} object using a key.
!
!     Overloaded \textit{value} for the following types:
!     \begin{itemize}
!       \item \texttt{integer(kind=ESMF\_KIND\_I4)}
!       \item \texttt{integer(kind=ESMF\_KIND\_I8)}
!       \item \texttt{real(kind=ESMF\_KIND\_R4)}
!       \item \texttt{real(kind=ESMF\_KIND\_R8)}
!       \item \texttt{logical}
!       \item \texttt{character(:)}
!     \end{itemize}
!
!     The arguments are:
!     \begin{description}
!     \item [info]
!       Target \texttt{ESMF\_Info} object.
!     \item [key]
!       String key to access in \texttt{ESMF\_Info} storage. See section \ref{info_key_format}
!       for an overview of the key format.
!     \item [value]
!       The input value associated with the key.
!     \item [{[force]}]
!       Default is true. When true, insert the key even if it already exists in
!       storage. If false, \textit{rc} will not return {\tt ESMF\_SUCCESS} if the
!       key already exists.
!     \item [{[idx]}]
!       An integer index to set if the target key's value is a list.
!     \item [{[pkey]}]
!       Use this key's location as the origin for the set call.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! NOTE: Documentation stub located above for generic interface compliance.
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoSetR4()"
subroutine ESMF_InfoSetR4(info, key, value, keywordEnforcer, force, idx, pkey, rc)
  type(ESMF_Info), intent(inout) :: info
  character(len=*), intent(in) :: key
  real(ESMF_KIND_R4), intent(in) :: value
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  logical, intent(in), optional :: force
  integer, intent(in), optional :: idx
  character(len=*), intent(in), optional :: pkey
  integer, intent(out), optional :: rc

  integer :: localrc
  logical(C_BOOL) :: local_force
  integer(C_INT), target :: local_idx
  type(C_PTR) :: local_idx_ptr
  character(:), allocatable :: local_pkey

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE

  if (present(force)) then
    local_force = force
  else
    local_force = .true.
  end if
  if (present(idx)) then
    local_idx = idx - 1  ! Shift to C (zero-based) indexing
    local_idx_ptr = C_LOC(local_idx)
  else
    local_idx_ptr = C_NULL_PTR
  end if
  if (present(pkey)) then
    local_pkey = TRIM(pkey)//C_NULL_CHAR
  else
    local_pkey = ""//C_NULL_CHAR
  end if

  call c_info_set_R4(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    value, &
    local_force, &
    localrc, &
    local_idx_ptr, &
    local_pkey)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (present(rc)) rc = ESMF_SUCCESS
end subroutine ESMF_InfoSetR4
#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoSetR8()"
subroutine ESMF_InfoSetR8(info, key, value, keywordEnforcer, force, idx, pkey, rc)
  type(ESMF_Info), intent(inout) :: info
  character(len=*), intent(in) :: key
  real(ESMF_KIND_R8), intent(in) :: value
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  logical, intent(in), optional :: force
  integer, intent(in), optional :: idx
  character(len=*), intent(in), optional :: pkey
  integer, intent(out), optional :: rc

  integer :: localrc
  logical(C_BOOL) :: local_force
  integer(C_INT), target :: local_idx
  type(C_PTR) :: local_idx_ptr
  character(:), allocatable :: local_pkey

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE

  if (present(force)) then
    local_force = force
  else
    local_force = .true.
  end if
  if (present(idx)) then
    local_idx = idx - 1  ! Shift to C (zero-based) indexing
    local_idx_ptr = C_LOC(local_idx)
  else
    local_idx_ptr = C_NULL_PTR
  end if
  if (present(pkey)) then
    local_pkey = TRIM(pkey)//C_NULL_CHAR
  else
    local_pkey = ""//C_NULL_CHAR
  end if

  call c_info_set_R8(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    value, &
    local_force, &
    localrc, &
    local_idx_ptr, &
    local_pkey)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (present(rc)) rc = ESMF_SUCCESS
end subroutine ESMF_InfoSetR8
#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoSetI4()"
subroutine ESMF_InfoSetI4(info, key, value, keywordEnforcer, force, idx, pkey, rc)
  type(ESMF_Info), intent(inout) :: info
  character(len=*), intent(in) :: key
  integer(ESMF_KIND_I4), intent(in) :: value
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  logical, intent(in), optional :: force
  integer, intent(in), optional :: idx
  character(len=*), intent(in), optional :: pkey
  integer, intent(out), optional :: rc

  integer :: localrc
  logical(C_BOOL) :: local_force
  integer(C_INT), target :: local_idx
  type(C_PTR) :: local_idx_ptr
  character(:), allocatable :: local_pkey

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE

  if (present(force)) then
    local_force = force
  else
    local_force = .true.
  end if
  if (present(idx)) then
    local_idx = idx - 1  ! Shift to C (zero-based) indexing
    local_idx_ptr = C_LOC(local_idx)
  else
    local_idx_ptr = C_NULL_PTR
  end if
  if (present(pkey)) then
    local_pkey = TRIM(pkey)//C_NULL_CHAR
  else
    local_pkey = ""//C_NULL_CHAR
  end if

  call c_info_set_I4(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    value, &
    local_force, &
    localrc, &
    local_idx_ptr, &
    local_pkey)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (present(rc)) rc = ESMF_SUCCESS
end subroutine ESMF_InfoSetI4
#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoSetI8()"
subroutine ESMF_InfoSetI8(info, key, value, keywordEnforcer, force, idx, pkey, rc)
  type(ESMF_Info), intent(inout) :: info
  character(len=*), intent(in) :: key
  integer(ESMF_KIND_I8), intent(in) :: value
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  logical, intent(in), optional :: force
  integer, intent(in), optional :: idx
  character(len=*), intent(in), optional :: pkey
  integer, intent(out), optional :: rc

  integer :: localrc
  logical(C_BOOL) :: local_force
  integer(C_INT), target :: local_idx
  type(C_PTR) :: local_idx_ptr
  character(:), allocatable :: local_pkey

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE

  if (present(force)) then
    local_force = force
  else
    local_force = .true.
  end if
  if (present(idx)) then
    local_idx = idx - 1  ! Shift to C (zero-based) indexing
    local_idx_ptr = C_LOC(local_idx)
  else
    local_idx_ptr = C_NULL_PTR
  end if
  if (present(pkey)) then
    local_pkey = TRIM(pkey)//C_NULL_CHAR
  else
    local_pkey = ""//C_NULL_CHAR
  end if

  call c_info_set_I8(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    value, &
    local_force, &
    localrc, &
    local_idx_ptr, &
    local_pkey)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (present(rc)) rc = ESMF_SUCCESS
end subroutine ESMF_InfoSetI8
#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoSetCH()"
subroutine ESMF_InfoSetCH(info, key, value, keywordEnforcer, force, idx, pkey, rc)
  type(ESMF_Info), intent(inout) :: info
  character(len=*), intent(in) :: key
  character(len=*), intent(in) :: value
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  logical, intent(in), optional :: force
  integer, intent(in), optional :: idx
  character(len=*), intent(in), optional :: pkey
  integer, intent(out), optional :: rc

  integer :: localrc
  logical(C_BOOL) :: local_force
  integer(C_INT), target :: local_idx
  type(C_PTR) :: local_idx_ptr
  character(:), allocatable :: local_pkey

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE

  if (present(force)) then
    local_force = force
  else
    local_force = .true.
  end if
  if (present(idx)) then
    local_idx = idx - 1  ! Shift to C (zero-based) indexing
    local_idx_ptr = C_LOC(local_idx)
  else
    local_idx_ptr = C_NULL_PTR
  end if
  if (present(pkey)) then
    local_pkey = TRIM(pkey)//C_NULL_CHAR
  else
    local_pkey = ""//C_NULL_CHAR
  end if

  call c_info_set_CH(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    trim(value)//C_NULL_CHAR, &
    local_force, &
    localrc, &
    local_idx_ptr, &
    local_pkey)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  if (present(rc)) rc = ESMF_SUCCESS
end subroutine ESMF_InfoSetCH
#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoSetLG()"
subroutine ESMF_InfoSetLG(info, key, value, keywordEnforcer, force, idx, pkey, rc)
  type(ESMF_Info), intent(inout) :: info
  character(len=*), intent(in) :: key
  logical, intent(in) :: value
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  logical, intent(in), optional :: force
  integer, intent(in), optional :: idx
  character(len=*), intent(in), optional :: pkey
  integer, intent(out), optional :: rc

  integer :: localrc
  logical(C_BOOL) :: local_force, local_value
  integer(C_INT), target :: local_idx
  type(C_PTR) :: local_idx_ptr
  character(:), allocatable :: local_pkey

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE

  if (present(force)) then
    local_force = force
  else
    local_force = .true.
  end if
  if (present(idx)) then
    local_idx = idx - 1  ! Shift to C (zero-based) indexing
    local_idx_ptr = C_LOC(local_idx)
  else
    local_idx_ptr = C_NULL_PTR
  end if
  if (present(pkey)) then
    local_pkey = TRIM(pkey)//C_NULL_CHAR
  else
    local_pkey = ""//C_NULL_CHAR
  end if

  local_value = value
  call c_info_set_LG(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    local_value, &
    local_force, &
    localrc, &
    local_idx_ptr, &
    local_pkey)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

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

! -----------------------------------------------------------------------------

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoSetINFO()"
!BOP
! !IROUTINE: ESMF_InfoSet - Set a key to the contents of an Info object
!
! !INTERFACE:
  ! Private name; call using ESMF_InfoSet
subroutine ESMF_InfoSetINFO(info, key, value, keywordEnforcer, force, rc)
! !ARGUMENTS:
  type(ESMF_Info), intent(inout) :: info
  character(len=*), intent(in) :: key
  type(ESMF_Info), intent(in) :: value
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  logical, intent(in), optional :: force
  integer, intent(out), optional :: rc
!
! !DESCRIPTION:
!     Set a value to the contents of an \texttt{ESMF\_Info} object. A copy of
!     the source contents is made.
!
!     The arguments are:
!     \begin{description}
!     \item [info]
!       Target \texttt{ESMF\_Info} object.
!     \item [key]
!       String key to access in \texttt{ESMF\_Info} storage. See section \ref{info_key_format}
!       for an overview of the key format.
!     \item [value]
!       The \texttt{ESMF\_Info} object to use as source data.
!     \item [{[force]}]
!       Default is true. When true, insert the key even if it already exists in
!       storage. If false, \textit{rc} will not return {\tt ESMF\_SUCCESS} if the
!       key already exists.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

  integer :: localrc
  logical(C_BOOL) :: local_force

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE

  if (present(force)) then
    local_force = force
  else
    local_force = .true.
  end if

  call c_info_set_INFO(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    value%ptr, &
    local_force, &
    localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

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

! -----------------------------------------------------------------------------

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoSetHConfig()"
!BOP
! !IROUTINE: ESMF_InfoSet - Set contents from a HConfig object
!
! !INTERFACE:
  ! Private name; call using ESMF_InfoSet
recursive subroutine ESMF_InfoSetHConfig(info, value, keywordEnforcer, keyPrefix, force, rc)
! !ARGUMENTS:
  type(ESMF_Info), intent(inout) :: info
  type(ESMF_HConfig), intent(in) :: value
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  character(len=*), intent(in), optional :: keyPrefix
  logical, intent(in), optional :: force
  integer, intent(out), optional :: rc
!
! !DESCRIPTION:
!     The provided \texttt{ESMF\_HConfig} object is expected to be a {\em map}.
!     An error is returned if this condition is not met. Each key-value pair
!     held by the \texttt{ESMF\_HConfig} object is added to the
!     \texttt{ESMF\_Info} object. A copy of the source contents is made.
!
!     Transfer of {\em scalar}, {\em sequence}, and {\em map} values
!     from \texttt{ESMF\_HConfig} to \texttt{ESMF\_Info} are supported.
!     Maps are treated recursively. Sequences are restricted to scalar elements
!     of the same typekind.
!
!     The keys of any map provided by the \texttt{ESMF\_HConfig} object must
!     be of scalar type. Keys are interpreted as strings when transferred to the
!     \texttt{ESMF\_Info} object. YAML merge keys "<<" are supported.
!
!     When existing keys in {\tt info} are overridden by this operation, the
!     typekind of the associated value element is allowed to change.
!
!     The arguments are:
!     \begin{description}
!     \item [info]
!       Target \texttt{ESMF\_Info} object.
!     \item [value]
!       The \texttt{ESMF\_HConfig} object to use as source data.
!     \item [{[keyPrefix]}]
!       If provided, prepend {\tt keyPrefix} to each of the keys found in the
!       {\tt value} map.
!     \item [{[force]}]
!       Default is true. When true, insert the key even if it already exists in
!       storage. If false, \textit{rc} will not return {\tt ESMF\_SUCCESS} if the
!       key already exists.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

  integer                             :: localrc
  logical                             :: isFlag
  type(ESMF_HConfigIter)              :: hconfigIterBegin, hconfigIterEnd
  type(ESMF_HConfigIter)              :: hconfigIter
  character(len=:), allocatable       :: key, tag, fullKey, msgString
  character(len=:), allocatable       :: valueStr
  integer(ESMF_KIND_I4)               :: valueInt
  real(ESMF_KIND_R4)                  :: valueFloat
  logical                             :: valueBool
  type(ESMF_HConfig)                  :: valueHConfig
  type(ESMF_Info)                     :: valueInfo
  character(len=:), allocatable       :: valueStrSeq(:)
  integer(ESMF_KIND_I4), allocatable  :: valueIntSeq(:)
  real(ESMF_KIND_R4), allocatable     :: valueFloatSeq(:)
  logical, allocatable                :: valueBoolSeq(:)

  if (present(rc)) rc = ESMF_SUCCESS

  isFlag = ESMF_HConfigIsNull(value, rc=localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
    ESMF_CONTEXT, rcToReturn=rc)) return

  if (isFlag) return  ! noop for NULL value

  isFlag = ESMF_HConfigIsDefined(value, rc=localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
    ESMF_CONTEXT, rcToReturn=rc)) return

  if (.not.isFlag) return  ! noop for not defined value

  isFlag = ESMF_HConfigIsMap(value, rc=localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
    ESMF_CONTEXT, rcToReturn=rc)) return

  if (.not.isFlag) then
    call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, &
      msg="Value must be HConfig map", &
      ESMF_CONTEXT, rcToReturn=rc)
    return
  endif

  hconfigIterBegin = ESMF_HConfigIterBegin(value, rc=localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
    ESMF_CONTEXT, rcToReturn=rc)) return

  hconfigIterEnd = ESMF_HConfigIterEnd(value, rc=localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
    ESMF_CONTEXT, rcToReturn=rc)) return

  hconfigIter = hconfigIterBegin
  do while (ESMF_HConfigIterLoop(hconfigIter, hconfigIterBegin, hconfigIterEnd, &
    rc=localrc))
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
      ESMF_CONTEXT, rcToReturn=rc)) return

    key = ESMF_HConfigAsStringMapKey(hconfigIter, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
      ESMF_CONTEXT, rcToReturn=rc)) return

    tag = ESMF_HConfigGetTagMapVal(hconfigIter, rc=localrc)
    if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
      ESMF_CONTEXT, rcToReturn=rc)) return

    if (key=="<<" .and. tag=="tag:yaml.org,2002:map") then
      ! dealing with YAML merge key -> recursivey handle it

      valueHConfig = ESMF_HConfigCreateAtMapVal(hconfigIter, rc=localrc)
      if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
        ESMF_CONTEXT, rcToReturn=rc)) return
      call ESMF_InfoSet(info, valueHConfig, keyPrefix=keyPrefix, force=force, &
        rc=localrc)
      if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
        ESMF_CONTEXT, rcToReturn=rc)) return
      ! clean-up
      call ESMF_HConfigDestroy(valueHConfig, rc=localrc)
      if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
        ESMF_CONTEXT, rcToReturn=rc)) return

    else
      ! regular key

      ! determine full key to be used for adding into info
      if (present(keyPrefix)) then
        fullKey=trim(keyPrefix)//"/"//key
      else
        fullKey=key
      endif

      ! set entry at full key to null to prevent conflict if typekind changes
      call ESMF_InfoSetNull(info, key=fullKey, rc=localrc)
      if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
        ESMF_CONTEXT, rcToReturn=rc)) return

      if (tag=="tag:yaml.org,2002:str") then
        valueStr = ESMF_HConfigAsStringMapVal(hconfigIter, rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
        call ESMF_InfoSet(info, key=fullKey, value=valueStr, force=force, &
          rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
      else if (tag=="tag:yaml.org,2002:bool") then
        valueBool = ESMF_HConfigAsLogicalMapVal(hconfigIter, rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
        call ESMF_InfoSet(info, key=fullKey, value=valueBool, force=force, &
          rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
      else if (tag=="tag:yaml.org,2002:int") then
        valueInt = ESMF_HConfigAsI4MapVal(hconfigIter, rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
        call ESMF_InfoSet(info, key=fullKey, value=valueInt, force=force, &
          rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
      else if (tag=="tag:yaml.org,2002:float") then
        valueFloat = ESMF_HConfigAsR4MapVal(hconfigIter, rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
        call ESMF_InfoSet(info, key=fullKey, value=valueFloat, force=force, &
          rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
      else if (tag=="tag:yaml.org,2002:map") then
        ! ESMF_Info supports maps recursively... go for it...
        valueHConfig = ESMF_HConfigCreateAtMapVal(hconfigIter, rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
        valueInfo = ESMF_InfoCreate(rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
        ! recursive call to set up the info object from hconfig
        call ESMF_InfoSet(valueInfo, valueHConfig, force=force, rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
        ! insert info under the respective key
        call ESMF_InfoSet(info, key=fullKey, value=valueInfo, force=force, &
          rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
        ! clean-up
        call ESMF_HConfigDestroy(valueHConfig, rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
        call ESMF_InfoDestroy(valueInfo, rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
      else if (tag=="tag:yaml.org,2002:seq") then
        ! ESMF_Info supports sequences only supported as 1d vectors same typekind
        ! ...detect the typekind by looking at the first element
        tag = ESMF_HConfigGetTagMapVal(hconfigIter, index=1, rc=localrc)
        if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
          ESMF_CONTEXT, rcToReturn=rc)) return
        if (tag=="tag:yaml.org,2002:str") then
          valueStrSeq = ESMF_HConfigAsStringSeqMapVal(hconfigIter, &
            stringLen=ESMF_MAXSTR, rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
          call ESMF_InfoSet(info, key=fullKey, values=valueStrSeq, force=force, &
            rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
        else if (tag=="tag:yaml.org,2002:bool") then
          valueBoolSeq = ESMF_HConfigAsLogicalSeqMapVal(hconfigIter, rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
          call ESMF_InfoSet(info, key=fullKey, values=valueBoolSeq, force=force, &
            rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
        else if (tag=="tag:yaml.org,2002:int") then
          valueIntSeq = ESMF_HConfigAsI4SeqMapVal(hconfigIter, rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
          call ESMF_InfoSet(info, key=fullKey, values=valueIntSeq, force=force, &
            rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
        else if (tag=="tag:yaml.org,2002:float") then
          valueFloatSeq = ESMF_HConfigAsR4SeqMapVal(hconfigIter, rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
          call ESMF_InfoSet(info, key=fullKey, values=valueFloatSeq, force=force, &
            rc=localrc)
          if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
        else
          msgString = "Unsupported typekind for sequence conversion, tag="//tag
          call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, &
            msg="Unsupported typekind for sequence conversion", &
            ESMF_CONTEXT, rcToReturn=rc)
          return
        endif
      else
        msgString = "Unsupported typekind, tag="//tag
        call ESMF_LogSetError(rcToCheck=ESMF_RC_ARG_BAD, &
          msg=msgString, &
          ESMF_CONTEXT, rcToReturn=rc)
        return
      endif

    endif

  enddo

end subroutine ESMF_InfoSetHConfig

!------------------------------------------------------------------------------
! SetArray --------------------------------------------------------------------
!------------------------------------------------------------------------------

!BOP
! !IROUTINE: ESMF_InfoSet - Set a value list
!
! !INTERFACE:
!subroutine ESMF_InfoSet(info, key, values, keywordEnforcer, force, pkey, rc)
!
! !ARGUMENTS:
!  type(ESMF_Info), intent(inout) :: info
!  character(len=*), intent(in) :: key
!  <values>, see below for supported values
!type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
!  logical, intent(in), optional :: force
!  character(len=*), intent(in), optional :: pkey
!  integer, intent(out), optional :: rc
!
! !DESCRIPTION:
!     Set a value list in an \texttt{ESMF\_Info} object using a key. List values
!     are initialized to null.
!
!     Overloaded \textit{values} for the following types:
!     \begin{itemize}
!       \item \texttt{integer(kind=ESMF\_KIND\_I4), dimension(:)}
!       \item \texttt{integer(kind=ESMF\_KIND\_I8), dimension(:)}
!       \item \texttt{real(kind=ESMF\_KIND\_R4), dimension(:)}
!       \item \texttt{real(kind=ESMF\_KIND\_R8), dimension(:)}
!       \item \texttt{logical, dimension(:)}
!       \item \texttt{character(:), dimension(:)}
!     \end{itemize}
!
!     The arguments are:
!     \begin{description}
!     \item [info]
!       Target \texttt{ESMF\_Info} object.
!     \item [key]
!       String key to access in \texttt{ESMF\_Info} storage. See section \ref{info_key_format}
!       for an overview of the key format.
!     \item [values]
!       The input value list associated with the key.
!     \item [{[force]}]
!       Default is true. When true, insert the key even if it already exists in
!       storage. If false, \textit{rc} will not return {\tt ESMF\_SUCCESS} if the
!       key already exists.
!     \item [{[pkey]}]
!       Use this key's location as the origin for the set call. Used primarily
!       for recursive requirements related to \texttt{ESMF\_Attribute}.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! NOTE: Documentation stub located above for generic interface compliance.
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoSetArrayR4()"
subroutine ESMF_InfoSetArrayR4(info, key, values, keywordEnforcer, force, pkey, rc)
  type(ESMF_Info), intent(inout) :: info
  character(len=*), intent(in) :: key
  real(ESMF_KIND_R4), dimension(:), intent(in) :: values
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  logical, intent(in), optional :: force
  character(len=*), intent(in), optional :: pkey
  integer, intent(out), optional :: rc

  integer :: localrc
  logical(C_BOOL) :: local_force
  character(:), allocatable :: local_pkey
  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE

  if (present(force)) then
    local_force = force
  else
    local_force = .true.
  end if
  if (present(pkey)) then
    local_pkey = TRIM(pkey)//C_NULL_CHAR
  else
    local_pkey = ""//C_NULL_CHAR
  end if

  call c_info_set_array_R4(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    values, &
    SIZE(values), &
    local_force, &
    localrc, &
    local_pkey)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

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

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoSetArrayR8()"
subroutine ESMF_InfoSetArrayR8(info, key, values, keywordEnforcer, force, pkey, rc)
  type(ESMF_Info), intent(inout) :: info
  character(len=*), intent(in) :: key
  real(ESMF_KIND_R8), dimension(:), intent(in) :: values
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  logical, intent(in), optional :: force
  character(len=*), intent(in), optional :: pkey
  integer, intent(out), optional :: rc

  integer :: localrc
  logical(C_BOOL) :: local_force
  character(:), allocatable :: local_pkey
  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE

  if (present(force)) then
    local_force = force
  else
    local_force = .true.
  end if
  if (present(pkey)) then
    local_pkey = TRIM(pkey)//C_NULL_CHAR
  else
    local_pkey = ""//C_NULL_CHAR
  end if

  call c_info_set_array_R8(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    values, &
    SIZE(values), &
    local_force, &
    localrc, &
    local_pkey)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

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

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoSetArrayI4()"
subroutine ESMF_InfoSetArrayI4(info, key, values, keywordEnforcer, force, pkey, rc)
  type(ESMF_Info), intent(inout) :: info
  character(len=*), intent(in) :: key
  integer(ESMF_KIND_I4), dimension(:), intent(in) :: values
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  logical, intent(in), optional :: force
  character(len=*), intent(in), optional :: pkey
  integer, intent(out), optional :: rc

  integer :: localrc
  logical(C_BOOL) :: local_force
  character(:), allocatable :: local_pkey
  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE

  if (present(force)) then
    local_force = force
  else
    local_force = .true.
  end if
  if (present(pkey)) then
    local_pkey = TRIM(pkey)//C_NULL_CHAR
  else
    local_pkey = ""//C_NULL_CHAR
  end if

  call c_info_set_array_I4(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    values, &
    SIZE(values), &
    local_force, &
    localrc, &
    local_pkey)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

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

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoSetArrayI8()"
subroutine ESMF_InfoSetArrayI8(info, key, values, keywordEnforcer, force, pkey, rc)
  type(ESMF_Info), intent(inout) :: info
  character(len=*), intent(in) :: key
  integer(ESMF_KIND_I8), dimension(:), intent(in) :: values
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  logical, intent(in), optional :: force
  character(len=*), intent(in), optional :: pkey
  integer, intent(out), optional :: rc

  integer :: localrc
  logical(C_BOOL) :: local_force
  character(:), allocatable :: local_pkey
  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE

  if (present(force)) then
    local_force = force
  else
    local_force = .true.
  end if
  if (present(pkey)) then
    local_pkey = TRIM(pkey)//C_NULL_CHAR
  else
    local_pkey = ""//C_NULL_CHAR
  end if

  call c_info_set_array_I8(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    values, &
    SIZE(values), &
    local_force, &
    localrc, &
    local_pkey)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

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

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoSetArrayCH()"
subroutine ESMF_InfoSetArrayCH(info, key, values, keywordEnforcer, force, pkey, rc)
  type(ESMF_Info), intent(inout) :: info
  character(len=*), intent(in) :: key
  character(len=*), dimension(:), intent(in) :: values
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  logical, intent(in), optional :: force
  character(len=*), intent(in), optional :: pkey
  integer, intent(out), optional :: rc

  integer :: localrc
  logical(C_BOOL) :: local_force
  integer :: ii
  integer(C_INT) :: idx
  character(:), allocatable :: local_pkey
  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE

  if (present(force)) then
    local_force = force
  else
    local_force = .true.
  end if
  if (present(pkey)) then
    local_pkey = TRIM(pkey)//C_NULL_CHAR
  else
    local_pkey = ""//C_NULL_CHAR
  end if

  ! Allocate storage in C
  call c_info_set_array_CH(info%ptr, trim(key)//C_NULL_CHAR, &
    SIZE(values), local_force, localrc, local_pkey)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  ! Set each character element in the underlying store
  do ii=1,SIZE(values)
    call ESMF_InfoSetCH(info, key, values(ii), idx=ii, pkey=local_pkey, rc=localrc)
  enddo
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

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

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoSetArrayLG()"
subroutine ESMF_InfoSetArrayLG(info, key, values, keywordEnforcer, force, pkey, rc)
  type(ESMF_Info), intent(inout) :: info
  character(len=*), intent(in) :: key
  logical, dimension(:), intent(in) :: values
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  logical, intent(in), optional :: force
  character(len=*), intent(in), optional :: pkey
  integer, intent(out), optional :: rc

  integer :: localrc
  logical(C_BOOL) :: local_force
  integer :: ii
  logical(C_BOOL), dimension(:), allocatable :: local_values
  character(:), allocatable :: local_pkey
  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE

  if (present(force)) then
    local_force = force
  else
    local_force = .true.
  end if
  if (present(pkey)) then
    local_pkey = TRIM(pkey)//C_NULL_CHAR
  else
    local_pkey = ""//C_NULL_CHAR
  end if

  allocate(local_values(SIZE(values)))
  do ii=1,SIZE(values)
    local_values(ii) = values(ii)
  enddo

  call c_info_set_array_LG(&
    info%ptr, &
    trim(key)//C_NULL_CHAR, &
    local_values, &
    SIZE(values), &
    local_force, &
    localrc, &
    local_pkey)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

  deallocate(local_values)

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


!------------------------------------------------------------------------------

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoSetNULL()"
!BOP
! !IROUTINE: ESMF_InfoSetNULL - Set a value to null
!
! !INTERFACE:
subroutine ESMF_InfoSetNULL(info, key, keywordEnforcer, force, rc)
! !ARGUMENTS:
  type(ESMF_Info), intent(inout) :: info
  character(len=*), intent(in) :: key
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  logical, intent(in), optional :: force
  integer, intent(out), optional :: rc
!
! !DESCRIPTION:
!     Set a value to null \cite{json_for_modern_cpp_null}.
!
!     The arguments are:
!     \begin{description}
!     \item [info]
!       Target \texttt{ESMF\_Info} object.
!     \item [key]
!       String key to access in \texttt{ESMF\_Info} storage. See section \ref{info_key_format}
!       for an overview of the key format.
!     \item [{[force]}]
!       Default is true. When true, insert the key even if it already exists in
!       storage. If false, \textit{rc} will not return {\tt ESMF\_SUCCESS} if the
!       key already exists.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

  integer :: localrc
  logical(C_BOOL) :: local_force

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE

  if (present(force)) then
    local_force = force
  else
    local_force = .true.
  end if

  call c_info_set_NULL(info%ptr, trim(key)//C_NULL_CHAR, local_force, localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

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

!------------------------------------------------------------------------------

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoSetDirty()"
!BOPI
! !IROUTINE: ESMF_InfoSetDirty - Change dirty state for synchronize operations
!
! !INTERFACE:
subroutine ESMF_InfoSetDirty(info, flag, keywordEnforcer, rc)
! !ARGUMENTS:
  type(ESMF_Info), intent(inout) :: info
  logical, intent(in) :: flag
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  integer, intent(out), optional :: rc
! !DESCRIPTION:
!     Change the dirty state of an \texttt{ESMF\_Info} object.
!
!     The arguments are:
!     \begin{description}
!     \item [info]
!       Target \texttt{ESMF\_Info} object.
!     \item [flag]
!       Logical value to set.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOPI

  integer :: localrc, local_flag

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE

  if (flag) then
    local_flag = 1 !true
  else
    local_flag = 0 !false
  end if

  call c_info_set_dirty(info%ptr, local_flag, localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

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

!------------------------------------------------------------------------------
! Sync Documentation ----------------------------------------------------------
!------------------------------------------------------------------------------

!BOP
! !IROUTINE: ESMF_InfoSync - Synchronize Info contents across a VM
! \label{esmf_infosync}
!
! !INTERFACE:
!subroutine ESMF_InfoSync(host, rootPet, vm, keywordEnforcer, markClean, &
!   rc)
! !ARGUMENTS:
!  type(ESMF_*), intent(inout) :: host
!  integer, intent(in) :: rootPet
!  type(ESMF_VM), intent(in) :: vm
!type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
!  logical, intent(in), optional :: markClean
!  integer, intent(out), optional :: rc
!
! !DESCRIPTION:
!     Synchronize \texttt{ESMF\_Info} contents collectively across the current VM.
!     Contents on the \textit{rootPet} are set as the contents on matching objects
!     sharing the VM. An attempt is made to optimize by only communicating updated
!     contents (i.e. something set or modified). This subroutine will traverse
!     the ESMF object hierarchy associated with \texttt{host} (i.e. Arrays in
!     an ArrayBundle, Fields in a FieldBundle, etc.).
!
!     Users interested in broadcasting only the \texttt{ESMF\_Info} object should
!     consult the \texttt{ESMF\_InfoBroadcast} documentation \ref{esmf_infobroadcast}.
!
!     The arguments are:
!     \begin{description}
!     \item [host]
!       Target ESMF object. Overloaded for:
!       \begin{itemize}
!         \item \texttt{ESMF\_State}
!         \item \texttt{ESMF\_CplComp}
!         \item \texttt{ESMF\_GridComp}
!         \item \texttt{ESMF\_SciComp}
!         \item \texttt{ESMF\_Field}
!         \item \texttt{ESMF\_FieldBundle}
!       \end{itemize}
!     \item [rootPet]
!       The root PET to use for the synchronization.
!     \item [vm]
!       The VM to synchronize across.
!     \item [{[markClean]}]
!       Default is false. If true, mark changed \texttt{ESMF\_Info} contents as
!       clean once synchronized. These contents will no longer be broadcast in
!       consecutive calls to \texttt{ESMF\_InfoSync}.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

!------------------------------------------------------------------------------

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoUpdate()"
!BOP
! !IROUTINE: ESMF_InfoUpdate - Update the contents of an Info object
!
! !INTERFACE:
subroutine ESMF_InfoUpdate(lhs, rhs, keywordEnforcer, recursive, overwrite, rc)
! !ARGUMENTS:
  type(ESMF_Info), intent(inout) :: lhs
  type(ESMF_Info), intent(in) :: rhs
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  logical, intent(in), optional :: recursive
  logical, intent(in), optional :: overwrite
  integer, intent(out), optional :: rc
!
! !DESCRIPTION:
!     Update the contents of \textit{lhs} using the contents of \textit{rhs}. The
!     operation inserts or overwrites any key in \textit{lhs} if it exists in \textit{rhs}.
!     Otherwise, the contents of \textit{lhs} is left unaltered. See the JSON
!     documentation for implementation details \cite{json_for_modern_cpp_update}.
!     If \textit{recursive} is \textit{.true.} (default is \texttt{.false.}),
!     nested objects will be updated by their component key/values. Otherwise,
!     the first instance or top-level key is replaced without the child contents
!     being updated element-by-element.
!
!     The arguments are:
!     \begin{description}
!     \item [lhs]
!       The \texttt{ESMF\_Info} object to update.
!     \item [rhs]
!       The \texttt{ESMF\_Info} object whose contents are used to update \textit{lhs}.
!     \item [{[recursive]}]
!       Default is \texttt{.false.}. If \texttt{.true.}, descend into nested objects
!       and recursively update the contents.
!     \item [{[overwrite]}]
!       Default is \texttt{.false.}. If \texttt{.true.}, key-values that exist
!       in \textit{lhs} will be overwritten by key-values in \textit{rhs}. Flag
!       is only applicable when \textit{recursive} is \texttt{.true.}.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

  integer :: localrc
  integer :: recursive_int, overwrite_int

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE
  recursive_int = 0  ! .false.
  if (present(recursive)) then
    if (recursive) then
      recursive_int = 1  ! .true.
    end if
  end if
  overwrite_int = 0  ! .false.
  if (present(overwrite)) then
    if (overwrite) then
      overwrite_int = 1  ! .true.
    end if
  end if

  call c_info_update(lhs%ptr, rhs%ptr, recursive_int, overwrite_int, localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

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

!------------------------------------------------------------------------------

#undef  ESMF_METHOD
#define ESMF_METHOD "ESMF_InfoWriteJSON()"
!BOP
! !IROUTINE: ESMF_InfoWriteJSON - Write Info contents to file
!
! !INTERFACE:
subroutine ESMF_InfoWriteJSON(info, filename, keywordEnforcer, rc)
! !ARGUMENTS:
  type(ESMF_Info), intent(in) :: info
  character(len=*), intent(in) :: filename
type(ESMF_KeywordEnforcer), optional:: keywordEnforcer ! must use keywords below
  integer, intent(out), optional :: rc
!
! !DESCRIPTION:
!     Write \texttt{ESMF\_Info} contents to file using the JSON format.
!
!     The arguments are:
!     \begin{description}
!     \item [info]
!       Target \texttt{ESMF\_Info} object.
!     \item [filename]
!       Path to the output file.
!     \item [{[rc]}]
!       Return code; equals {\tt ESMF\_SUCCESS} if there are no errors.
!     \end{description}
!EOP

  integer :: localrc

  localrc = ESMF_FAILURE
  if (present(rc)) rc = ESMF_FAILURE

  call c_info_write_json(info%ptr, trim(filename)//C_NULL_CHAR, localrc)
  if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, ESMF_CONTEXT, rcToReturn=rc)) return

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

end module ESMF_InfoMod  !=====================================================