Subsections

3 API


3.1 Generic Component: NUOPC_Driver


MODULE:

  module NUOPC_Driver


DESCRIPTION:
Component that drives and coordinates initialization of its child components: Model, Mediator, and Connector components. For every Driver time step the same run sequence, i.e. sequence of Model, Mediator, and Connector Run methods is called. The run sequence is fully customizable. The default run sequence implements explicit time stepping.


SUPER:

  ESMF_GridComp


USE DEPENDENCIES:

  use ESMF


SETSERVICES:

  subroutine SetServices(driver, rc)
    type(ESMF_GridComp)   :: driver
    integer, intent(out)  :: rc


SEMANTIC SPECIALIZATION LABELS:

3.1.1 NUOPC_DriverAddComp - Add a GridComp child to a Driver


INTERFACE:

   ! Private name; call using NUOPC_DriverAddComp()
   recursive subroutine NUOPC_DriverAddGridComp(driver, compLabel, &
     compSetServicesRoutine, compSetVMRoutine, petList, info, config, comp, rc)
ARGUMENTS:
     type(ESMF_GridComp)                               :: driver
     character(len=*),    intent(in)                   :: compLabel
 #if defined (__NVCOMPILER) || defined (__PGI) || defined (ESMF_COMPILER_AOCC)
     interface
       recursive subroutine compSetServicesRoutine(gridcomp, rc)
         use ESMF
         implicit none
         type(ESMF_GridComp)        :: gridcomp ! must not be optional
         integer, intent(out)       :: rc       ! must not be optional
       end subroutine
     end interface
     interface
       recursive subroutine compSetVMRoutine(gridcomp, rc)
         use ESMF
         implicit none
         type(ESMF_GridComp)        :: gridcomp ! must not be optional
         integer, intent(out)       :: rc       ! must not be optional
       end subroutine
     end interface
     optional                                          :: compSetVMRoutine
 #else
     abstract interface
       recursive subroutine SetServicesRoutine(gridcomp, rc)
         use ESMF
         implicit none
         type(ESMF_GridComp)        :: gridcomp ! must not be optional
         integer, intent(out)       :: rc       ! must not be optional
       end subroutine
       recursive subroutine SetVMRoutine(gridcomp, rc)
         use ESMF
         implicit none
         type(ESMF_GridComp)        :: gridcomp ! must not be optional
         integer, intent(out)       :: rc       ! must not be optional
       end subroutine
     end interface
     procedure(SetServicesRoutine)                     :: compSetServicesRoutine
     procedure(SetVMRoutine),                 optional :: compSetVMRoutine
 #endif
     integer,             intent(in),         optional :: petList(:)
     type(ESMF_Info),     intent(in),         optional :: info
     type(ESMF_Config),   intent(in),         optional :: config
     type(ESMF_GridComp), intent(out),        optional :: comp
     integer,             intent(out),        optional :: rc
DESCRIPTION:

Create and add a GridComp (i.e. Model, Mediator, or Driver) as a child component to a Driver. The component is created on the provided petList, or by default across all of the Driver PETs.

The specified compSetServicesRoutine() is called back immediately after the new child component has been created internally. Very little around the component is set up at that time (e.g. NUOPC component attributes will not be available). The routine should therefore be very light weight, with the sole purpose of setting the entry points of the component - typically by deriving from a generic component followed by the appropriate specilizations.

If provided, the compSetVMRoutine() is called back before the compSetServicesRoutine(). This allows the child component to set aspects of its own VM, such as threading or the PE distribution among PETs.

The info argument can be used to pass custom attributes to the child component. These attributes are available on the component when compSetVMRoutine() and compSetServicesRoutine() are called. The attributes provided in info are copied onto the child component. This allows the same info object to be used for multiple child components without conflict.

The compLabel must uniquely identify the child component within the context of the Driver component.

If the comp argument is specified, it will reference the newly created component on return.

3.1.2 NUOPC_DriverAddComp - Add a GridComp child from shared object to a Driver


INTERFACE:

   ! Private name; call using NUOPC_DriverAddComp()
   recursive subroutine NUOPC_DriverAddGridCompSO(driver, compLabel, &
     sharedObj, petList, info, config, comp, rc)
ARGUMENTS:
     type(ESMF_GridComp)                        :: driver
     character(len=*),    intent(in)            :: compLabel
     character(len=*),    intent(in),  optional :: sharedObj
     integer,             intent(in),  optional :: petList(:)
     type(ESMF_Info),     intent(in),  optional :: info
     type(ESMF_Config),   intent(in),  optional :: config
     type(ESMF_GridComp), intent(out), optional :: comp
     integer,             intent(out), optional :: rc
DESCRIPTION:

Create and add a GridComp (i.e. Model, Mediator, or Driver) as a child component to a Driver. The component is created on the provided petList, or by default across all of the Driver PETs.

The SetServices() routine in the sharedObj is called back immediately after the new child component has been created internally. Very little around the component is set up at that time (e.g. NUOPC component attributes will not be available). The routine should therefore be very light weight, with the sole purpose of setting the entry points of the component - typically by deriving from a generic component followed by the appropriate specilizations.

The info argument can be used to pass custom attributes to the child component. These attributes are available on the component when compSetVMRoutine() and compSetServicesRoutine() are called. The attributes provided in info are copied onto the child component. This allows the same info object to be used for multiple child components without conflict.

The compLabel must uniquely identify the child component within the context of the Driver component.

If the comp argument is specified, it will reference the newly created component on return.

3.1.3 NUOPC_DriverAddComp - Add a CplComp child to a Driver


INTERFACE:

   ! Private name; call using NUOPC_DriverAddComp()
   recursive subroutine NUOPC_DriverAddCplComp(driver, srcCompLabel, &
     dstCompLabel, compSetServicesRoutine, compSetVMRoutine, petList, info, &
     config, comp, rc)
ARGUMENTS:
     type(ESMF_GridComp)                               :: driver
     character(len=*),    intent(in)                   :: srcCompLabel
     character(len=*),    intent(in)                   :: dstCompLabel
 #if defined (__NVCOMPILER) || defined (__PGI) || defined (ESMF_COMPILER_AOCC)
     interface
       recursive subroutine compSetServicesRoutine(cplcomp, rc)
         use ESMF
         implicit none
         type(ESMF_CplComp)         :: cplcomp  ! must not be optional
         integer, intent(out)       :: rc       ! must not be optional
       end subroutine
     end interface
     interface
       recursive subroutine compSetVMRoutine(cplcomp, rc)
         use ESMF
         implicit none
         type(ESMF_CplComp)         :: cplcomp  ! must not be optional
         integer, intent(out)       :: rc       ! must not be optional
       end subroutine
     end interface
     optional                                          :: compSetVMRoutine
 #else
     abstract interface
       recursive subroutine SetServicesRoutine(cplcomp, rc)
         use ESMF
         implicit none
         type(ESMF_CplComp)         :: cplcomp  ! must not be optional
         integer, intent(out)       :: rc       ! must not be optional
       end subroutine
       recursive subroutine SetVMRoutine(cplcomp, rc)
         use ESMF
         implicit none
         type(ESMF_CplComp)         :: cplcomp  ! must not be optional
         integer, intent(out)       :: rc       ! must not be optional
       end subroutine
     end interface
     procedure(SetServicesRoutine)                     :: compSetServicesRoutine
     procedure(SetVMRoutine),                 optional :: compSetVMRoutine
 #endif
     integer, target,     intent(in),         optional :: petList(:)
     type(ESMF_Info),     intent(in),         optional :: info
     type(ESMF_Config),   intent(in),         optional :: config
     type(ESMF_CplComp),  intent(out),        optional :: comp
     integer,             intent(out),        optional :: rc
DESCRIPTION:

Create and add a CplComp (i.e. Connector) as a child component to a Driver. The component is created on the provided petList, or by default across the union of PETs of the components indicated by srcCompLabel and dstCompLabel.

The specified SetServices() routine is called back immediately after the new child component has been created internally. Very little around the component is set up at that time (e.g. NUOPC component attributes will not be available). The routine should therefore be very light weight, with the sole purpose of setting the entry points of the component - typically by deriving from a generic component followed by the appropriate specilizations.

The info argument can be used to pass custom attributes to the child component. These attributes are available on the component when compSetVMRoutine() and compSetServicesRoutine() are called. The attributes provided in info are copied onto the child component. This allows the same info object to be used for multiple child components without conflict.

The compLabel must uniquely identify the child component within the context of the Driver component.

If the comp argument is specified, it will reference the newly created component on return.

3.1.4 NUOPC_DriverAddRunElement - Add RunElement for Model, Mediator, or Driver


INTERFACE:

   ! Private name; call using NUOPC_DriverAddRunElement()
   recursive subroutine NUOPC_DriverAddRunElementMPL(driver, slot, compLabel, &
     phaseLabel, relaxedflag, rc)
ARGUMENTS:
     type(ESMF_GridComp)                        :: driver
     integer,             intent(in)            :: slot
     character(len=*),    intent(in)            :: compLabel
 -- The following arguments require argument keyword syntax (e.g. rc=rc). --
     character(len=*),    intent(in),  optional :: phaseLabel
     logical,             intent(in),  optional :: relaxedflag
     integer,             intent(out), optional :: rc
DESCRIPTION:

Add an element associated with a Model, Mediator, or Driver component to the run sequence of the Driver. The component must have been added to the Driver, and associated with compLabel prior to this call.

If phaseLabel was not specified, the first entry in the RunPhaseMap attribute of the referenced component will be used to determine the run phase of the added element.

By default an error is returned if no component is associated with the specified compLabel. This error can be suppressed by setting relaxedflag=.true., and no entry will be added to the run sequence.

The slot number identifies the run sequence time slot in case multiple sequences are available. Slots start counting from 1.

3.1.5 NUOPC_DriverAddRunElement - Add RunElement for Connector


INTERFACE:

   ! Private name; call using NUOPC_DriverAddRunElement()
   recursive subroutine NUOPC_DriverAddRunElementCPL(driver, slot, srcCompLabel,&
     dstCompLabel, phaseLabel, relaxedflag, rc)
ARGUMENTS:
     type(ESMF_GridComp)                        :: driver
     integer,             intent(in)            :: slot
     character(len=*),    intent(in)            :: srcCompLabel
     character(len=*),    intent(in)            :: dstCompLabel
 -- The following arguments require argument keyword syntax (e.g. rc=rc). --
     character(len=*),    intent(in),  optional :: phaseLabel
     logical,             intent(in),  optional :: relaxedflag
     integer,             intent(out), optional :: rc
DESCRIPTION:

Add an element associated with a Connector component to the run sequence of the Driver. The component must have been added to the Driver, and associated with srcCompLabel and dstCompLabel prior to this call.

If phaseLabel was not specified, the first entry in the RunPhaseMap attribute of the referenced component will be used to determine the run phase of the added element.

By default an error is returned if no component is associated with the specified compLabel. This error can be suppressed by setting relaxedflag=.true., and no entry will be added to the run sequence.

The slot number identifies the run sequence time slot in case multiple sequences are available. Slots start counting from 1.

3.1.6 NUOPC_DriverAddRunElement - Add RunElement that links to another slot


INTERFACE:

   ! Private name; call using NUOPC_DriverAddRunElement()
   recursive subroutine NUOPC_DriverAddRunElementL(driver, slot, linkSlot, rc)
ARGUMENTS:
     type(ESMF_GridComp)                        :: driver
     integer,             intent(in)            :: slot
     integer,             intent(in)            :: linkSlot
     integer,             intent(out), optional :: rc
DESCRIPTION:

Add an element to the run sequence of the Driver that links to the time slot indicated by linkSlot.

3.1.7 NUOPC_DriverEgestRunSequence - Egest the run sequence as FreeFormat


INTERFACE:

   recursive subroutine NUOPC_DriverEgestRunSequence(driver, freeFormat, rc)
ARGUMENTS:
     type(ESMF_GridComp)                           :: driver
     type(NUOPC_FreeFormat), intent(out)           :: freeFormat
     integer,                intent(out), optional :: rc
DESCRIPTION:

Egest the run sequence stored in the driver as a FreeFormat object. It is the caller's responsibility to destroy the created freeFormat object.

3.1.8 NUOPC_DriverGet - Get info from a Driver


INTERFACE:

   ! Private name; call using NUOPC_DriverGet()
   recursive subroutine NUOPC_DriverGet(driver, slotCount, parentClock, &
     importState, exportState, rc)
ARGUMENTS:
     type(ESMF_GridComp)                        :: driver
     integer,             intent(out), optional :: slotCount
     type(ESMF_Clock),    intent(out), optional :: parentClock
     type(ESMF_State),    intent(out), optional :: importState
     type(ESMF_State),    intent(out), optional :: exportState
     integer,             intent(out), optional :: rc
DESCRIPTION:

Access Driver information.

3.1.9 NUOPC_DriverGetComp - Get a GridComp child from a Driver


INTERFACE:

   ! Private name; call using NUOPC_DriverGetComp()
   recursive subroutine NUOPC_DriverGetGridComp(driver, compLabel, comp, petList, &
     importState, exportState, relaxedflag, rc)
ARGUMENTS:
     type(ESMF_GridComp)                        :: driver
     character(len=*),    intent(in)            :: compLabel
     type(ESMF_GridComp), intent(out), optional :: comp
     integer,             pointer,     optional :: petList(:)
     type(ESMF_State),    intent(out), optional :: importState
     type(ESMF_State),    intent(out), optional :: exportState
     logical,             intent(in),  optional :: relaxedflag
     integer,             intent(out), optional :: rc
DESCRIPTION:

Query the Driver for a GridComp (i.e. Model, Mediator, or Driver) child component that was added under compLabel.

If provided, the petList argument will be associated with the petList that was used to create the referenced component. This is an internal allocation owned by the library. This pointer must not be deallocated by the user!

By default an error is returned if no component is associated with the specified compLabel. This error can be suppressed by setting relaxedflag=.true., and unassociated arguments will be returned.

3.1.10 NUOPC_DriverGetComp - Get a CplComp child from a Driver


INTERFACE:

   ! Private name; call using NUOPC_DriverGetComp()
   recursive subroutine NUOPC_DriverGetCplComp(driver, srcCompLabel, &
     dstCompLabel, comp, petList, relaxedflag, rc)
ARGUMENTS:
     type(ESMF_GridComp)                        :: driver
     character(len=*),    intent(in)            :: srcCompLabel
     character(len=*),    intent(in)            :: dstCompLabel
     type(ESMF_CplComp),  intent(out), optional :: comp
     integer,             pointer    , optional :: petList(:)
     logical,             intent(in),  optional :: relaxedflag
     integer,             intent(out), optional :: rc
DESCRIPTION:

Query the Driver for a CplComp (i.e. Connector) child component that was added under compLabel.

If provided, the petList argument will be associated with the petList that was used to create the referenced component. This is an internal allocation owned by the library. This pointer must not be deallocated by the user!

By default an error is returned if no component is associated with the specified compLabel. This error can be suppressed by setting relaxedflag=.true., and unassociated arguments will be returned.

3.1.11 NUOPC_DriverGetComp - Get all the GridComp child components from a Driver


INTERFACE:

   ! Private name; call using NUOPC_DriverGetComp()
   recursive subroutine NUOPC_DriverGetAllGridComp(driver, compList, petLists, &
     rc)
ARGUMENTS:
     type(ESMF_GridComp)                        :: driver
     type(ESMF_GridComp), pointer, optional     :: compList(:)
     type(ESMF_PtrInt1D), pointer, optional     :: petLists(:)
     integer,             intent(out), optional :: rc
DESCRIPTION:

Get all the GridComp (i.e. Model, Mediator, or Driver) child components from a Driver.

The incoming compList and petLists arguments must enter unassociated. This means that the user code must explicitly call nullify() or use the => null() syntax on the variables passed in as the actual arguments.

On return it becomes the responsibility of the caller to deallocate associated compList and petLists arguments:

     if (associated(compList)) deallocate(compList)
     if (associated(petLists)) deallocate(petLists)

Notice that the petLists(i)%ptr members, if associated, are pointing to an internal allocation owned by the library. These pointers must not be deallocated by the user!

3.1.12 NUOPC_DriverGetComp - Get all the CplComp child components from a Driver


INTERFACE:

   ! Private name; call using NUOPC_DriverGetComp()
   recursive subroutine NUOPC_DriverGetAllCplComp(driver, compList, petLists, rc)
ARGUMENTS:
     type(ESMF_GridComp)                        :: driver
     type(ESMF_CplComp),  pointer               :: compList(:)
     type(ESMF_PtrInt1D), pointer, optional     :: petLists(:)
     integer,             intent(out), optional :: rc
DESCRIPTION:

Get all the CplComp (i.e. Connector) child components from a Driver.

The incoming compList and petLists arguments must enter unassociated. This means that the user code must explicitly call nullify() or use the => null() syntax on the variables passed in as the actual arguments.

On return it becomes the responsibility of the caller to deallocate associated compList and petLists arguments:

     if (associated(compList)) deallocate(compList)
     if (associated(petLists)) deallocate(petLists)

Notice that the petLists(i)%ptr members, if associated, are pointing to an internal allocation owned by the library. These pointers must not be deallocated by the user!

3.1.13 NUOPC_DriverIngestRunSequence - Ingest the run sequence from FreeFormat


INTERFACE:

   ! Private name; call using NUOPC_DriverIngestRunSequence()
   recursive subroutine NUOPC_DriverIngestRunSequenceFF(driver, freeFormat, &
     autoAddConnectors, rc)
ARGUMENTS:
     type(ESMF_GridComp)                           :: driver
     type(NUOPC_FreeFormat), intent(in),  target   :: freeFormat
     logical,                intent(in),  optional :: autoAddConnectors
     integer,                intent(out), optional :: rc
DESCRIPTION:

Ingest the run sequence from a FreeFormat object and replace the run sequence currently held by the driver. Every line in freeFormat corresponds to either a component run sequence element, or is part of a time loop or alarm block defintion. Anything following a '#' character on a line is considered a comment, and ignored for the purpose of ingesting run sequence elements.

Component run sequence elements define the run method of a single component. The lines are interpreted sequentially, however, components will execute concurrently as long as this is not prevented by data-dependencies or overlapping petLists.

Each line specifies the precise run method phase for a single component instance. For model, mediator, and driver components the format is this:

     compLabel [phaseLabel]
Here compLabel is the label by which the component instance is known to the driver. It is optionally followed a phaseLabel identifying a specific run phase. An example of calling the run phase of the ATM instance that contains the "fast" processes, and is labeled fast:

     ATM fast
By default, i.e. without phaseLabel, the first registered run method of the component is used.

The format for connector components is different. It looks like this:

     srcCompLabel -> dstCompLabel [connectionOptions]
A connector instance is uniquely known by the two components it connects, i.e. by srcCompLabel and dstCompLabel. The syntax requires that the token -> be specified between source and destination. Optionally connectionOptions can be supplied using the format discussed under section 2.4.5. The connection options are set as attribute ConnectionOptions on the respective connector component.

An example of executing the connector instance that transfers fields from the ATM component to the OCN component, using redistribution for remapping:

     ATM -> OCN :remapMethod=redist

By default autoAddConnectors is .false., which means that all components referenced in the freeFormat run sequence, including connectors, must already be available as child components of the driver component. An error will be returned if this is not the case. However, when autoAddConnectors is set to .true., connector components encountered in the run sequence that are no already present in the driver will be added automatically. The default NUOPC_Connector implementation is used for all automatically added connector instances.

Lines that contain a time loop definition have the general format:

     @{timeStep|*}[:runDuration]
       ...
       ...
     @
Both timeStep and runDuration are numbers in units of seconds. Time loops can be nested and concatenated.

A wildcard "*" character can be specified in place of an actual timeStep number. In this case the timeStep of the associated run clock object is set to be equal to the timeStep of the time loop one level up in the loop nesting hierarchy. If a wildcard time step is used for a single outer time loop in the run sequence, then the associated run clock is identical to the driver clock and must be set explicitly by the driver code, or its parent component.

The runDuration specification is optional. If omitted, the duration of the associated run clock is set to the timeStep of the time loop one level up in the loop nesting hierarchy. This ensures that for a single nested time loop, the loop returns to the parent loop level at the appropriate time.

A simple example of a single time loop with one hour timestep:

     @3600
       ...
       ...
     @
Each time loop has its own associated clock object. NUOPC manages these clock objects, i.e. their creation and destruction, as well as startTime, endTime, timeStep adjustments during the execution. The outer most time loop of the run sequence is a special case. It uses the driver clock itself. If a single outer most loop is defined in the run sequence provided by freeFormat, this loop becomes the driver loop level directly. Therefore, setting the timeStep or runDuration for the outer most time loop results modifiying the driver clock itself. However, for cases with concatenated loops on the upper level of the run sequence in freeFormat, a single outer loop is added automatically during ingestion, and the driver clock is used for this loop instead.

A more complex run sequence example, that shows component run sequence elements outside of time loops, a nested time loop, time step wildcards, explicit duration specifications, and concatenated time loops:

     @100:800
       ATM -> OCN
       OCN -> ATM
       ATM
       OCN
       @*
         OCN -> EXTOCN
         EXTOCN
       @
     @
     ATM -> EXTATM
     EXTATM
     @100:1000
       ATM -> OCN
       OCN -> ATM
       ATM
       OCN
     @
Here the timeStep of the first time loop is explicitly chosen at $100s$. The runDuration is explicitly set to $800s$. The first time loop steps the current time forward for $800s$, for each iteration executing ATM-OCN coupling, followed by the nested loop that calls the OCN -> EXTOCN and EXTOCN components. The nested loop uses a wildcard timeStep and therefore is identical to the parent loop level timeStep of $100s$. The nested runDuration is not specified and therefore also defaults to the parent time step of $100s$. In other words, the nested loop is executed exactly once for every parent loop iteration.

After $800s$ the first time loop is exited, and followed by explicit calls to ATM -> EXTAMT and EXTATM components. Finally the second time loop is entered for another $1000s$ runDuration. The timeStep is again explicitly set to $100s$. The second time loop only implements ATM-OCN coupling, and no coupling to EXTOCN is implemented. Finally, after $1800s$ the sequence returns to the driver level loop.

Lines that contain an alarm block definition have the general format:

     @@{alarmTime|*}
       ...
       ...
     @@
The alarmTime is a number in units of seconds, and indicates at which interval the alarm will ring. The first ring time of an alarm is the current time of the parent clock.

Specification of the wildcard character * sets the alarmTime equal to the timeStep of the parentClock.

When an alarm rings, the entire alarm block is executed once.

Nesting of time loops and alarm blocks is supported.

3.1.14 NUOPC_DriverIngestRunSequence - Ingest the run sequence from HConfig


INTERFACE:

   ! Private name; call using NUOPC_DriverIngestRunSequence()
   recursive subroutine NUOPC_DriverIngestRunSequenceHC(driver, hconfig, &
     autoAddConnectors, rc)
ARGUMENTS:
     type(ESMF_GridComp)                           :: driver
     type(ESMF_HConfig),     intent(in)            :: hconfig
     logical,                intent(in),  optional :: autoAddConnectors
     integer,                intent(out), optional :: rc
DESCRIPTION:

Ingest the run sequence from a HConfig object and replace the run sequence currently held by the driver. The provided hconfig must be a scalar, or else an error is returned. The scalar is interpreted as a string, broken into lines at the newline character. Each line is subsequently interpreted according to the rules described under the FreeFormat version of the NUOPC_DriverIngestRunSequence() interface.

To preserve newline characters in run sequences expressed in YAML block notation, it is important to use literals indicated by the '|' character in YAML. For example:

   # A simple run sequence example as a YAML block literal
   --- |
    @900:1800          # comments are ignored
      MED
      MED -> ATM       # any line can have a comment
      MED -> OCN
      ATM
      OCN
      ATM -> MED
      OCN -> MED
    @

Notice the leading whitespace character(s) on each line of the block literal string. YAML requires at least one (1) leading whitespace character for strings in block notation.

3.1.15 NUOPC_DriverNewRunSequence - Replace the run sequence in a Driver


INTERFACE:

   recursive subroutine NUOPC_DriverNewRunSequence(driver, slotCount, rc)
ARGUMENTS:
     type(ESMF_GridComp)                        :: driver
     integer,             intent(in)            :: slotCount
     integer,             intent(out), optional :: rc
DESCRIPTION:

Replace the current run sequence of the Driver with a new one that has slotCount slots. Each slot uses its own clock for time keeping.

3.1.16 NUOPC_DriverPrint - Print internal Driver information


INTERFACE:

   recursive subroutine NUOPC_DriverPrint(driver, orderflag, rc)
ARGUMENTS:
     type(ESMF_GridComp)                        :: driver
     logical,             intent(in),  optional :: orderflag
     integer,             intent(out), optional :: rc
DESCRIPTION:

Print internal Driver information. If orderflag is provided and set to .true., the output is ordered from lowest to highest PET. Setting this flag makes the method collective.

3.1.17 NUOPC_DriverSetRunSequence - Set internals of RunSequence slot


INTERFACE:

   ! Private name; call using NUOPC_DriverSetRunSequence()
   recursive subroutine NUOPC_DriverSetRunSequence(driver, slot, clock, alarm, rc)
ARGUMENTS:
     type(ESMF_GridComp)                        :: driver
     integer,             intent(in)            :: slot
     type(ESMF_Clock),    intent(in)            :: clock
     type(ESMF_Alarm),    intent(in),  optional :: alarm
     integer,             intent(out), optional :: rc
DESCRIPTION:

Set the clock in the run sequence under slot of the Driver.


3.2 Generic Component: NUOPC_ModelBase


MODULE:

  module NUOPC_ModelBase


DESCRIPTION:
Partial specialization of a component with a default explicit time dependency. Each time the Run method is called the component steps one timeStep forward on the passed in parent clock. The component flags incompatibility during Run if the current time of the incoming clock does not match the current time of the internal clock.


SUPER:

  ESMF_GridComp


USE DEPENDENCIES:

  use ESMF


SETSERVICES:

  subroutine SetServices(modelBase, rc)
    type(ESMF_GridComp)   :: modelBase
    integer, intent(out)  :: rc


SEMANTIC SPECIALIZATION LABELS:


3.3 Generic Component: NUOPC_Model


MODULE:

  module NUOPC_Model


DESCRIPTION:
Model component with a default explicit time dependency. Each time the Run method is called the model integrates one timeStep forward on the passed in parent clock. The internal clock is advanced at the end of each Run call. The component flags incompatibility during Run if the current time of the incoming clock does not match the current time of the internal clock.


SUPER:

  NUOPC_ModelBase


USE DEPENDENCIES:

  use ESMF


SETSERVICES:

  subroutine SetServices(model, rc)
    type(ESMF_GridComp)   :: model
    integer, intent(out)  :: rc


SEMANTIC SPECIALIZATION LABELS:

3.3.1 NUOPC_ModelGet - Get info from a Model


INTERFACE:

   subroutine NUOPC_ModelGet(model, driverClock, modelClock, &
     importState, exportState, rc)
ARGUMENTS:
     type(ESMF_GridComp)                        :: model
     type(ESMF_Clock),    intent(out), optional :: driverClock
     type(ESMF_Clock),    intent(out), optional :: modelClock
     type(ESMF_State),    intent(out), optional :: importState
     type(ESMF_State),    intent(out), optional :: exportState
     integer,             intent(out), optional :: rc
DESCRIPTION:

Access Model information.


3.4 Generic Component: NUOPC_Mediator


MODULE:

  module NUOPC_Mediator


DESCRIPTION:
Mediator component with a default explicit time dependency. Each time the Run method is called, the time stamp on the imported Fields must match the current time (on both the incoming and internal clock). Before returning, the Mediator time stamps the exported Fields with the same current time, before advancing the internal clock one timeStep forward.


SUPER:

  NUOPC_ModelBase


USE DEPENDENCIES:

  use ESMF


SETSERVICES:

  subroutine SetServices(mediator, rc)
    type(ESMF_GridComp)   :: mediator
    integer, intent(out)  :: rc


SEMANTIC SPECIALIZATION LABELS:

3.4.1 NUOPC_MediatorGet - Get info from a Mediator


INTERFACE:

   subroutine NUOPC_MediatorGet(mediator, driverClock, mediatorClock, &
     importState, exportState, rc)
ARGUMENTS:
     type(ESMF_GridComp)                        :: mediator
     type(ESMF_Clock),    intent(out), optional :: driverClock
     type(ESMF_Clock),    intent(out), optional :: mediatorClock
     type(ESMF_State),    intent(out), optional :: importState
     type(ESMF_State),    intent(out), optional :: exportState
     integer,             intent(out), optional :: rc
DESCRIPTION:

Access Mediator information.


3.5 Generic Component: NUOPC_Connector


MODULE:

  module NUOPC_Connector


DESCRIPTION:
Component that makes a unidirectional connection between model, mediator, and or driver components. During initialization field pairing is performed between the import and export side according to section 2.4.2, and paired fields are connected. By default the bilinear regrid method is used during Run to transfer data from the connected import Fields to the connected export Fields.


SUPER:

  ESMF_CplComp


USE DEPENDENCIES:

  use ESMF


SETSERVICES:

  subroutine SetServices(connector, rc)
    type(ESMF_CplComp)    :: connector
    integer, intent(out)  :: rc


SEMANTIC SPECIALIZATION LABELS:

3.5.1 NUOPC_ConnectorGet - Get parameters from a Connector


INTERFACE:

   subroutine NUOPC_ConnectorGet(connector, srcFields, dstFields, rh, state, &
     CplSet, cplSetList, srcVM, dstVM, driverClock, rc)
ARGUMENTS:
     type(ESMF_CplComp)                            :: connector
     type(ESMF_FieldBundle), intent(out), optional :: srcFields
     type(ESMF_FieldBundle), intent(out), optional :: dstFields
     type(ESMF_RouteHandle), intent(out), optional :: rh
     type(ESMF_State),       intent(out), optional :: state
     character(*),           intent(in),  optional :: CplSet
     character(ESMF_MAXSTR), pointer,     optional :: cplSetList(:)
     type(ESMF_VM),          intent(out), optional :: srcVM
     type(ESMF_VM),          intent(out), optional :: dstVM
     type(ESMF_Clock),       intent(out), optional :: driverClock
     integer,                intent(out), optional :: rc
DESCRIPTION:

Get parameters from the connector internal state.

The Connector keeps information about the connection that it implements in its internal state. When customizing a Connector, it is often necessary to access and sometimes modify these data objects.

The arguments are:

connector
The Connector component.
[srcFields]
The FieldBundle under which the Connector keeps track of all connected source side fields. The order in which the fields are stored in srcFields is significant, as it corresponds to the order of fields in dstFields. Consequently, when accessing and modifying the fields inside of srcFields, it is important to use the itemorderflag=ESMF_ITEMORDER_ADDORDER option to ESMF_FieldBundleGet().
[dstFields]
The FieldBundle under which the Connector keeps track of all connected destination side fields. The order in which the fields are stored in dstFields is significant, as it corresponds to the order of fields in srcFields. Consequently, when accessing and modifying the fields inside of dstFields, it is important to use the itemorderflag=ESMF_ITEMORDER_ADDORDER option to ESMF_FieldBundleGet().
[rh]
The RouteHandle that the Connector uses to move data from srcFields to dstFields.
[state]
A State object that the Connector keeps to make customization of the Connector more convenient. The generic Connector code handles creation and destruction of state, but does not access it directly for information.
[CplSet]
If present, all of the returned information is specific to the specified coupling set.
[cplSetList]
The list of coupling sets currently known to the Connector. This argument must enter the call unassociated or an error is returned. This means that the user code must explicitly call nullify() or use the => null() syntax on the variable passed in as cplSetList argument. On return, the cplSetList argument will be associated, potentially of size zero. The responsibility for deallocation transfers to the caller.
[srcVM]
The VM of the source side component.
[dstVM]
The VM of the destination side component.
[driverClock]
The Clock object used by the current RunSequence level to drive this component.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.5.2 NUOPC_ConnectorSet - Set parameters in a Connector


INTERFACE:

   subroutine NUOPC_ConnectorSet(connector, srcFields, dstFields, rh, state, &
     CplSet, srcVM, dstVM, rc)
ARGUMENTS:
     type(ESMF_CplComp)                            :: connector
     type(ESMF_FieldBundle), intent(in),  optional :: srcFields
     type(ESMF_FieldBundle), intent(in),  optional :: dstFields
     type(ESMF_RouteHandle), intent(in),  optional :: rh
     type(ESMF_State),       intent(in),  optional :: state
     character(*),           intent(in),  optional :: CplSet
     type(ESMF_VM),          intent(in),  optional :: srcVM
     type(ESMF_VM),          intent(in),  optional :: dstVM
     integer,                intent(out), optional :: rc
DESCRIPTION:

Set parameters in the connector internal state.

The Connector keeps information about the connection that it implements in its internal state. When customizing a Connector, it is often necessary to access and sometimes modify these data objects.

The arguments are:

connector
The Connector component.
[srcFields]
The FieldBundle under which the Connector keeps track of all connected source side fields. The order in which the fields are stored in srcFields is significant, as it corresponds to the order of fields in dstFields. Consequently, when setting srcFields, it is important to add them in the same order as for dstFields.
[dstFields]
The FieldBundle under which the Connector keeps track of all connected destination side fields. The order in which the fields are stored in dstFields is significant, as it corresponds to the order of fields in srcFields. Consequently, when setting dstFields, it is important to add them in the same order as for srcFields.
[rh]
The RouteHandle that the Connector uses to move data from srcFields to dstFields.
[state]
A State object that the Connector keeps to make customization of the Connector more convenient. Only in very rare cases would the user want to replace the state that is managed by the generic Connector implementation. If state is set by this call, the user essentially claims ownership of the previous state object, and becomes responsible for its destruction. Ownership of the new state is transferred to the Connector and must not be explicitly destroyed by the user code.
[CplSet]
If present, all of the passed in information is set under the specified coupling set.
[srcVM]
The VM of the source side component.
[dstVM]
The VM of the destination side component.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.6 General Generic Component Methods

3.6.1 NUOPC_CompAreServicesSet - Check if SetServices was called


INTERFACE:

   ! Private name; call using NUOPC_CompAreServicesSet() 
   function NUOPC_GridCompAreServicesSet(comp, rc)
RETURN VALUE:
     logical :: NUOPC_GridCompAreServicesSet
ARGUMENTS:
     type(ESMF_GridComp), intent(in)            :: comp
     integer,             intent(out), optional :: rc
DESCRIPTION:

Return .true. if SetServices has been called for comp. Otherwise return .false..

3.6.2 NUOPC_CompAreServicesSet - Check if SetServices was called


INTERFACE:

   ! Private name; call using NUOPC_CompAreServicesSet() 
   function NUOPC_CplCompAreServicesSet(comp, rc)
RETURN VALUE:
     logical :: NUOPC_CplCompAreServicesSet
ARGUMENTS:
     type(ESMF_CplComp), intent(in)            :: comp
     integer,            intent(out), optional :: rc
DESCRIPTION:

Return .true. if SetServices has been called for comp. Otherwise return .false..

3.6.3 NUOPC_CompAttributeAdd - Add NUOPC GridComp Attributes


INTERFACE:

   ! Private name; call using NUOPC_CompAttributeAdd() 
   subroutine NUOPC_GridCompAttributeAdd(comp, attrList, rc)
ARGUMENTS:
     type(ESMF_GridComp)                       :: comp
     character(len=*),   intent(in)            :: attrList(:)
     integer,            intent(out), optional :: rc
DESCRIPTION:

Add Attributes to the highest level of the standard NUOPC AttPack hierarchy (convention="NUOPC", purpose="Instance").

3.6.4 NUOPC_CompAttributeAdd - Add NUOPC CplComp Attributes


INTERFACE:

   ! Private name; call using NUOPC_CompAttributeAdd() 
   subroutine NUOPC_CplCompAttributeAdd(comp, attrList, rc)
ARGUMENTS:
     type(ESMF_CplComp)                        :: comp
     character(len=*),   intent(in)            :: attrList(:)
     integer,            intent(out), optional :: rc
DESCRIPTION:

Add Attributes to the highest level of the standard NUOPC AttPack hierarchy (convention="NUOPC", purpose="Instance").

3.6.5 NUOPC_CompAttributeEgest - Egest NUOPC GridComp Attributes in FreeFormat


INTERFACE:

   ! Private name; call using NUOPC_CompAttributeEgest() 
   subroutine NUOPC_GridCompAttributeEge(comp, freeFormat, rc)
ARGUMENTS:
     type(ESMF_GridComp),    intent(in)            :: comp
     type(NUOPC_FreeFormat), intent(out)           :: freeFormat
     integer,                intent(out), optional :: rc
DESCRIPTION:

Egest the Attributes of the highest level of the standard NUOPC AttPack hierarchy (convention="NUOPC", purpose="Instance") as a FreeFormat object. It is the caller's responsibility to destroy the created freeFormat object.

3.6.6 NUOPC_CompAttributeEgest - Egest NUOPC CplComp Attributes in FreeFormat


INTERFACE:

   ! Private name; call using NUOPC_CompAttributeEgest() 
   subroutine NUOPC_CplCompAttributeEge(comp, freeFormat, rc)
ARGUMENTS:
     type(ESMF_CplComp),     intent(in)            :: comp
     type(NUOPC_FreeFormat), intent(out)           :: freeFormat
     integer,                intent(out), optional :: rc
DESCRIPTION:

Egest the Attributes of the highest level of the standard NUOPC AttPack hierarchy (convention="NUOPC", purpose="Instance") as a FreeFormat object. It is the caller's responsibility to destroy the created freeFormat object.

3.6.7 NUOPC_CompAttributeGet - Get a NUOPC GridComp Attribute - string


INTERFACE:

   ! Private name; call using NUOPC_CompAttributeGet() 
   subroutine NUOPC_GridCompAttributeGet(comp, name, value, isPresent, isSet, rc)
ARGUMENTS:
     type(ESMF_GridComp), intent(in)            :: comp
     character(*),        intent(in)            :: name
     character(*),        intent(out)           :: value
     logical,             intent(out), optional :: isPresent
     logical,             intent(out), optional :: isSet
     integer,             intent(out), optional :: rc
DESCRIPTION:

Access the attribute name inside of comp using the convention NUOPC and purpose Instance.

This call assumes to find a scalar value. An error is returned otherwise.

This call concverts to a string value, regardless of the actual attribute storage.

Unless isPresent and isSet are provided, return with error if the attribute is not present or not set, respectively.

The arguments are:

comp
The ESMF_GridComp object to be queried.
name
The name of the queried attribute.
value
The value of the queried attribute.
[isPresent]
Set to .true. if the queried attribute is present, .false. otherwise.
[isSet]
Set to .true. if the queried attribute is set, .false. otherwise.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.6.8 NUOPC_CompAttributeGet - Get a NUOPC CplComp Attribute - string


INTERFACE:

   ! Private name; call using NUOPC_CompAttributeGet() 
   subroutine NUOPC_CplCompAttributeGet(comp, name, value, isPresent, isSet, rc)
ARGUMENTS:
     type(ESMF_CplComp),  intent(in)            :: comp
     character(*),        intent(in)            :: name
     character(*),        intent(out)           :: value
     logical,             intent(out), optional :: isPresent
     logical,             intent(out), optional :: isSet
     integer,             intent(out), optional :: rc
DESCRIPTION:

Access the attribute name inside of comp using the convention NUOPC and purpose Instance.

This call assumes to find a scalar value. An error is returned otherwise.

This call concverts to a string value, regardless of the actual attribute storage.

Unless isPresent and isSet are provided, return with error if the attribute is not present or not set, respectively.

The arguments are:

comp
The ESMF_CplComp object to be queried.
name
The name of the queried attribute.
value
The value of the queried attribute.
[isPresent]
Set to .true. if the queried attribute is present, .false. otherwise.
[isSet]
Set to .true. if the queried attribute is set, .false. otherwise.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.6.9 NUOPC_CompAttributeGet - Get a NUOPC GridComp Attribute - integer


INTERFACE:

   ! Private name; call using NUOPC_CompAttributeGet() 
   subroutine NUOPC_GridCompAttributeGetI(comp, name, value, isPresent, isSet, rc)
ARGUMENTS:
     type(ESMF_GridComp), intent(in)            :: comp
     character(*),        intent(in)            :: name
     integer,             intent(out)           :: value
     logical,             intent(out), optional :: isPresent
     logical,             intent(out), optional :: isSet
     integer,             intent(out), optional :: rc
DESCRIPTION:

Access the attribute name inside of comp using the convention NUOPC and purpose Instance.

Unless isPresent and isSet are provided, return with error if the attribute is not present or not set, respectively.

The arguments are:

comp
The ESMF_GridComp object to be queried.
name
The name of the queried attribute.
value
The value of the queried attribute.
[isPresent]
Set to .true. if the queried attribute is present, .false. otherwise.
[isSet]
Set to .true. if the queried attribute is set, .false. otherwise.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.6.10 NUOPC_CompAttributeGet - Get a NUOPC CplComp Attribute - integer


INTERFACE:

   ! Private name; call using NUOPC_CompAttributeGet() 
   subroutine NUOPC_CplCompAttributeGetI(comp, name, value, isPresent, isSet, rc)
ARGUMENTS:
     type(ESMF_CplComp),  intent(in)            :: comp
     character(*),        intent(in)            :: name
     integer,             intent(out)           :: value
     logical,             intent(out), optional :: isPresent
     logical,             intent(out), optional :: isSet
     integer,             intent(out), optional :: rc
DESCRIPTION:

Access the attribute name inside of comp using the convention NUOPC and purpose Instance.

Unless isPresent and isSet are provided, return with error if the attribute is not present or not set, respectively.

The arguments are:

comp
The ESMF_CplComp object to be queried.
name
The name of the queried attribute.
value
The value of the queried attribute.
[isPresent]
Set to .true. if the queried attribute is present, .false. otherwise.
[isSet]
Set to .true. if the queried attribute is set, .false. otherwise.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.6.11 NUOPC_CompAttributeGet - Get a NUOPC GridComp Attribute - string list


INTERFACE:

   ! Private name; call using NUOPC_CompAttributeGet() 
   subroutine NUOPC_GridCompAttributeGetSL(comp, name, valueList, isPresent, &
     isSet, itemCount, typekind, rc)
ARGUMENTS:
     type(ESMF_GridComp),       intent(in)            :: comp
     character(*),              intent(in)            :: name
     character(*),              intent(out), optional :: valueList(:)
     logical,                   intent(out), optional :: isPresent
     logical,                   intent(out), optional :: isSet
     integer,                   intent(out), optional :: itemCount
     type(ESMF_TypeKind_Flag),  intent(out), optional :: typekind
     integer,                   intent(out), optional :: rc
DESCRIPTION:

Access the attribute name inside of comp using the convention NUOPC and purpose Instance. Returns with error if the attribute is not present or not set.

Unless isPresent and isSet are provided, return with error if the attribute is not present or not set, respectively.

The arguments are:

comp
The ESMF_GridComp object to be queried.
name
The name of the queried attribute.
[valueList]
The list of values of the queried attribute.
[isPresent]
Set to .true. if the queried attribute is present, .false. otherwise.
[isSet]
Set to .true. if the queried attribute is set, .false. otherwise.
[itemCount]
Number of items in the attribute. Return 0 if not present or not set.
[typekind]
The typekind of the queried attribute.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.6.12 NUOPC_CompAttributeGet - Get a NUOPC CplComp Attribute - string list


INTERFACE:

   ! Private name; call using NUOPC_CompAttributeGet() 
   subroutine NUOPC_CplCompAttributeGetSL(comp, name, valueList, isPresent, &
     isSet, itemCount, typekind, rc)
ARGUMENTS:
     type(ESMF_CplComp),       intent(in)            :: comp
     character(*),             intent(in)            :: name
     character(*),             intent(out), optional :: valueList(:)
     logical,                  intent(out), optional :: isPresent
     logical,                  intent(out), optional :: isSet
     integer,                  intent(out), optional :: itemCount
     type(ESMF_TypeKind_Flag), intent(out), optional :: typekind
     integer,                  intent(out), optional :: rc
DESCRIPTION:

Access the attribute name inside of comp using the convention NUOPC and purpose Instance. Returns with error if the attribute is not present or not set.

Unless isPresent and isSet are provided, return with error if the attribute is not present or not set, respectively.

The arguments are:

comp
The ESMF_CplComp object to be queried.
name
The name of the queried attribute.
[valueList]
The list of values of the queried attribute.
[isPresent]
Set to .true. if the queried attribute is present, .false. otherwise.
[isSet]
Set to .true. if the queried attribute is set, .false. otherwise.
[itemCount]
Number of items in the attribute. Return 0 if not present or not set.
[typekind]
The typekind of the queried attribute.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.6.13 NUOPC_CompAttributeGet - Get a NUOPC GridComp Attribute - integer list


INTERFACE:

   ! Private name; call using NUOPC_CompAttributeGet() 
   subroutine NUOPC_GridCompAttributeGetIL(comp, name, valueList, isPresent, &
     isSet, itemCount, typekind, rc)
ARGUMENTS:
     type(ESMF_GridComp),       intent(in)            :: comp
     character(*),              intent(in)            :: name
     integer,                   intent(out)           :: valueList(:)
     logical,                   intent(out), optional :: isPresent
     logical,                   intent(out), optional :: isSet
     integer,                   intent(out), optional :: itemCount
     type(ESMF_TypeKind_Flag),  intent(out), optional :: typekind
     integer,                   intent(out), optional :: rc
DESCRIPTION:

Access the attribute name inside of comp using the convention NUOPC and purpose Instance. Returns with error if the attribute is not present or not set.

Unless isPresent and isSet are provided, return with error if the attribute is not present or not set, respectively.

The arguments are:

comp
The ESMF_GridComp object to be queried.
name
The name of the queried attribute.
valueList
The list of values of the queried attribute.
[isPresent]
Set to .true. if the queried attribute is present, .false. otherwise.
[isSet]
Set to .true. if the queried attribute is set, .false. otherwise.
[itemCount]
Number of items in the attribute. Return 0 if not present or not set.
[typekind]
The typekind of the queried attribute.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.6.14 NUOPC_CompAttributeGet - Get a NUOPC CplComp Attribute - integer list


INTERFACE:

   ! Private name; call using NUOPC_CompAttributeGet() 
   subroutine NUOPC_CplCompAttributeGetIL(comp, name, valueList, isPresent, &
     isSet, itemCount, typekind, rc)
ARGUMENTS:
     type(ESMF_CplComp),       intent(in)            :: comp
     character(*),             intent(in)            :: name
     integer,                  intent(out)           :: valueList(:)
     logical,                  intent(out), optional :: isPresent
     logical,                  intent(out), optional :: isSet
     integer,                  intent(out), optional :: itemCount
     type(ESMF_TypeKind_Flag), intent(out), optional :: typekind
     integer,                  intent(out), optional :: rc
DESCRIPTION:

Access the attribute name inside of comp using the convention NUOPC and purpose Instance. Returns with error if the attribute is not present or not set.

Unless isPresent and isSet are provided, return with error if the attribute is not present or not set, respectively.

The arguments are:

comp
The ESMF_CplComp object to be queried.
name
The name of the queried attribute.
valueList
The list of values of the queried attribute.
[isPresent]
Set to .true. if the queried attribute is present, .false. otherwise.
[isSet]
Set to .true. if the queried attribute is set, .false. otherwise.
[itemCount]
Number of items in the attribute. Return 0 if not present or not set.
[typekind]
The typekind of the queried attribute.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.6.15 NUOPC_CompAttributeIngest - Ingest free format NUOPC GridComp Attributes


INTERFACE:

   ! Private name; call using NUOPC_CompAttributeIngest() 
   subroutine NUOPC_GridCompAttributeIng(comp, freeFormat, addFlag, rc)
ARGUMENTS:
     type(ESMF_GridComp),    intent(in)            :: comp
     type(NUOPC_FreeFormat), intent(in)            :: freeFormat
     logical,                intent(in),  optional :: addFlag
     integer,                intent(out), optional :: rc
DESCRIPTION:

Ingest the Attributes from a FreeFormat object onto the highest level of the standard NUOPC AttPack hierarchy (convention="NUOPC", purpose="Instance").

Important: Attributes ingested by this method are stored as type character strings, and must be accessed accordingly. Conversion from string into a different data type, e.g. integer or real, is the user's responsibility.

If addFlag is .false. (default), an error will be returned if an attribute is to be ingested that was not previously added to the comp object. If addFlag is .true., all missing attributes will be added by this method automatically as needed.

Each line in freeFormat is of this format:

       attributeName = attributeValue

For example:

       Verbosity  = 0
       Profiling  = 0
       Diagnostic = 0
could directly be ingested as Attributes for any instance of the four standard NUOPC component kinds. This is because Verbosity, Profiling, and Diagnostic are pre-defined Attributes of the NUOPC component kinds according to sections 2.3.1, 2.3.2, 2.3.3, and 2.3.4.

When Attributes are specified in freeFormat that are not pre-defined for a specific component kind, they can still be ingested by a component instance using the addFlag=.true. option. For instance:

       ModelOutputChoice = 2
specifies a user-level Attribute, which is not part of the pre-defined Attributes of any of the standard NUOPC component kinds.

3.6.16 NUOPC_CompAttributeIngest - Ingest free format NUOPC CplComp Attributes


INTERFACE:

   ! Private name; call using NUOPC_CompAttributeIngest() 
   subroutine NUOPC_CplCompAttributeIng(comp, freeFormat, addFlag, rc)
ARGUMENTS:
     type(ESMF_CplComp),     intent(in)            :: comp
     type(NUOPC_FreeFormat), intent(in)            :: freeFormat
     logical,                intent(in),  optional :: addFlag
     integer,                intent(out), optional :: rc
DESCRIPTION:

Ingest the Attributes from a FreeFormat object onto the highest level of the standard NUOPC AttPack hierarchy (convention="NUOPC", purpose="Instance").

Important: Attributes ingested by this method are stored as type character strings, and must be accessed accordingly. Conversion from string into a different data type, e.g. integer or real, is the user's responsibility.

If addFlag is .false. (default), an error will be returned if an attribute is to be ingested that was not previously added to the comp object. If addFlag is .true., all missing attributes will be added by this method automatically as needed.

Each line in freeFormat is of this format:

       attributeName = attributeValue

For example:

       Verbosity  = 0
       Profiling  = 0
       Diagnostic = 0
could directly be ingested as Attributes for any instance of the four standard NUOPC component kinds. This is because Verbosity, Profiling, and Diagnostic are pre-defined Attributes of the NUOPC component kinds according to sections 2.3.1, 2.3.2, 2.3.3, and 2.3.4.

When Attributes are specified in freeFormat that are not pre-defined for a specific component kind, they can still be ingested by a component instance using the addFlag=.true. option. For instance:

       ModelOutputChoice = 2
specifies a user-level Attribute, which is not part of the pre-defined Attributes of any of the standard NUOPC component kinds.

3.6.17 NUOPC_CompAttributeIngest - Ingest NUOPC GridComp Attributes from HConfig


INTERFACE:

   ! Private name; call using NUOPC_CompAttributeIngest() 
   subroutine NUOPC_GridCompAttributeIngHC(comp, hconfig, rc)
ARGUMENTS:
     type(ESMF_GridComp),    intent(in)            :: comp
     type(ESMF_HConfig),     intent(in)            :: hconfig
     integer,                intent(out), optional :: rc
DESCRIPTION:

Ingest component attributes from a HConfig object onto the highest level of the standard NUOPC AttPack hierarchy (convention="NUOPC", purpose="Instance").

The provided hconfig is expected to be a map. An error is returned if this condition is not met. Each key-value pair held by hconfig is added as an attribute to comp. A copy of the source contents is made.

Transfers of scalar, sequence, and map values from hconfig are supported. Maps are treated recursively. Sequences are restricted to scalar elements of the same typekind.

The keys of any map provided by the hconfig object must be of scalar type. Keys are interpreted as strings when transferred as an attribute.

Existing attributes with the same key are overridden by this operation. When attributes are overridden, the typekind of the associated value element is allowed to change.

   # A simple YAML definition of standard NUOPC attributes, followed by
   # component specific attributes.
  
   Verbosity:  4609             # decimal representation of explicit bit pattern
   Profiling:  low              # pre-defined NUOPC setting
   Diagnostic: 0                # explicit 0 turns OFF feature
   CustomSeq1: [1, 2, 3, 4]     # sequence of integers
   CustomSeq2: [1., 2., 3., 4.] # sequence of floats
   CustomSeq3: [true, false]    # sequence of bools
   CustomType: {k1: [a, aa, aaa], k2: b, k3: c}  # complex structure

3.6.18 NUOPC_CompAttributeIngest - Ingest NUOPC CplComp Attributes from HConfig


INTERFACE:

   ! Private name; call using NUOPC_CompAttributeIngest() 
   subroutine NUOPC_CplCompAttributeIngHC(comp, hconfig, rc)
ARGUMENTS:
     type(ESMF_CplComp),     intent(in)            :: comp
     type(ESMF_HConfig),     intent(in)            :: hconfig
     integer,                intent(out), optional :: rc
DESCRIPTION:

Ingest component attributes from a HConfig object onto the highest level of the standard NUOPC AttPack hierarchy (convention="NUOPC", purpose="Instance").

The provided hconfig is expected to be a map. An error is returned if this condition is not met. Each key-value pair held by hconfig is added as an attribute to comp. A copy of the source contents is made.

Transfers of scalar, sequence, and map values from hconfig are supported. Maps are treated recursively. Sequences are restricted to scalar elements of the same typekind.

The keys of any map provided by the hconfig object must be of scalar type. Keys are interpreted as strings when transferred as an attribute.

Existing attributes with the same key are overridden by this operation. When attributes are overridden, the typekind of the associated value element is allowed to change.

   # A simple YAML definition of standard NUOPC attributes, followed by
   # component specific attributes.
  
   Verbosity:  4609             # decimal representation of explicit bit pattern
   Profiling:  low              # pre-defined NUOPC setting
   Diagnostic: 0                # explicit 0 turns OFF feature
   CustomSeq1: [1, 2, 3, 4]     # sequence of integers
   CustomSeq2: [1., 2., 3., 4.] # sequence of floats
   CustomSeq3: [true, false]    # sequence of bools
   CustomType: {k1: [a, aa, aaa], k2: b, k3: c}  # complex structure

3.6.19 NUOPC_CompAttributeReset - Reset NUOPC GridComp Attributes


INTERFACE:

   ! Private name; call using NUOPC_CompAttributeReset() 
   subroutine NUOPC_GridCompAttributeReset(comp, attrList, rc)
ARGUMENTS:
     type(ESMF_GridComp)                       :: comp
     character(len=*),   intent(in)            :: attrList(:)
     integer,            intent(out), optional :: rc
DESCRIPTION:

Reset Attributes on the highest level of the standard NUOPC AttPack hierarchy (convention="NUOPC", purpose="Instance").

3.6.20 NUOPC_CompAttributeReset - Reset NUOPC CplComp Attributes


INTERFACE:

   ! Private name; call using NUOPC_CompAttributeReset() 
   subroutine NUOPC_CplCompAttributeReset(comp, attrList, rc)
ARGUMENTS:
     type(ESMF_CplComp)                        :: comp
     character(len=*),   intent(in)            :: attrList(:)
     integer,            intent(out), optional :: rc
DESCRIPTION:

Reset Attributes on the highest level of the standard NUOPC AttPack hierarchy (convention="NUOPC", purpose="Instance").

3.6.21 NUOPC_CompAttributeSet - Set a NUOPC GridComp Attribute


INTERFACE:

   ! Private name; call using NUOPC_CompAttributeSet() 
   subroutine NUOPC_GridCompAttributeSetS(comp, name, value, rc)
ARGUMENTS:
     type(ESMF_GridComp)                   :: comp
     character(*), intent(in)              :: name
     character(*), intent(in)              :: value
     integer,      intent(out), optional   :: rc
DESCRIPTION:

Set the Attribute name inside of comp on the highest level of the standard NUOPC AttPack hierarchy (convention="NUOPC", purpose="Instance").

3.6.22 NUOPC_CompAttributeSet - Set a NUOPC CplComp Attribute


INTERFACE:

   ! Private name; call using NUOPC_CompAttributeSet() 
   subroutine NUOPC_CplCompAttributeSetS(comp, name, value, rc)
ARGUMENTS:
     type(ESMF_CplComp)                    :: comp
     character(*), intent(in)              :: name
     character(*), intent(in)              :: value
     integer,      intent(out), optional   :: rc
DESCRIPTION:

Set the Attribute name inside of comp on the highest level of the standard NUOPC AttPack hierarchy (convention="NUOPC", purpose="Instance").

3.6.23 NUOPC_CompAttributeSet - Set a NUOPC GridComp Attribute


INTERFACE:

   ! Private name; call using NUOPC_CompAttributeSet() 
   subroutine NUOPC_GridCompAttributeSetI(comp, name, value, rc)
ARGUMENTS:
     type(ESMF_GridComp)                   :: comp
     character(*), intent(in)              :: name
     integer,      intent(in)              :: value
     integer,      intent(out), optional   :: rc
DESCRIPTION:

Set the Attribute name inside of comp on the highest level of the standard NUOPC AttPack hierarchy (convention="NUOPC", purpose="Instance").

3.6.24 NUOPC_CompAttributeSet - Set a NUOPC CplComp Attribute


INTERFACE:

   ! Private name; call using NUOPC_CompAttributeSet() 
   subroutine NUOPC_CplCompAttributeSetI(comp, name, value, rc)
ARGUMENTS:
     type(ESMF_CplComp)                    :: comp
     character(*), intent(in)              :: name
     integer,      intent(in)              :: value
     integer,      intent(out), optional   :: rc
DESCRIPTION:

Set the Attribute name inside of comp on the highest level of the standard NUOPC AttPack hierarchy (convention="NUOPC", purpose="Instance").

3.6.25 NUOPC_CompAttributeSet - Set a NUOPC GridComp List Attribute


INTERFACE:

   ! Private name; call using NUOPC_CompAttributeSet() 
   subroutine NUOPC_GridCompAttributeSetSL(comp, name, valueList, rc)
ARGUMENTS:
     type(ESMF_GridComp)                   :: comp
     character(*), intent(in)              :: name
     character(*), intent(in)              :: valueList(:)
     integer,      intent(out), optional   :: rc
DESCRIPTION:

Set the Attribute name inside of comp on the highest level of the standard NUOPC AttPack hierarchy (convention="NUOPC", purpose="Instance").

3.6.26 NUOPC_CompAttributeSet - Set a NUOPC CplComp List Attribute


INTERFACE:

   ! Private name; call using NUOPC_CompAttributeSet() 
   subroutine NUOPC_CplCompAttributeSetSL(comp, name, valueList, rc)
ARGUMENTS:
     type(ESMF_CplComp)                    :: comp
     character(*), intent(in)              :: name
     character(*), intent(in)              :: valueList(:)
     integer,      intent(out), optional   :: rc
DESCRIPTION:

Set the Attribute name inside of comp on the highest level of the standard NUOPC AttPack hierarchy (convention="NUOPC", purpose="Instance").

3.6.27 NUOPC_CompCheckSetClock - Check Clock compatibility and set stopTime


INTERFACE:

   ! Private name; call using NUOPC_CompCheckSetClock() 
   subroutine NUOPC_GridCompCheckSetClock(comp, externalClock, checkTimeStep, &
     forceTimeStep, rc)
ARGUMENTS:
     type(ESMF_GridComp),     intent(inout)         :: comp
     type(ESMF_Clock),        intent(in)            :: externalClock
     logical,                 intent(in),  optional :: checkTimeStep
     logical,                 intent(in),  optional :: forceTimeStep
     integer,                 intent(out), optional :: rc
DESCRIPTION:

Compare externalClock to the internal clock of comp to make sure they match in their current time. Also ensure that the time step of the external clock is a multiple of the time step of the internal clock. If both conditions are satisfied then set the stop time of the internal clock so it is reached in one time step of the external clock. Otherwise leave the internal clock unchanged and return with error. The direction of the involved clocks is taking into account. Setting the forceTimeStep argument to .true. forces the timeStep of the externalClock to be used to reset the timeStep of the internal clock.

3.6.28 NUOPC_CompDerive - Derive a GridComp from a generic component


INTERFACE:

   ! Private name; call using NUOPC_CompDerive() 
   recursive subroutine NUOPC_GridCompDerive(comp, genericSetServicesRoutine, rc)
ARGUMENTS:
     type(ESMF_GridComp), intent(in)            :: comp
     interface
       subroutine genericSetServicesRoutine(gridcomp, rc)
         use ESMF
         implicit none
         type(ESMF_GridComp)        :: gridcomp ! must not be optional
         integer, intent(out)       :: rc       ! must not be optional
       end subroutine
     end interface
     integer,             intent(out), optional :: rc
DESCRIPTION:

Derive a GridComp (i.e. Model, Mediator, or Driver) from a generic component by calling into the specified SetServices() routine of the generic component. This is typically the first call in the SetServices() routine of the specializing component, and is followed by NUOPC_CompSpecialize() calls.

3.6.29 NUOPC_CompDerive - Derive a CplComp from a generic component


INTERFACE:

   ! Private name; call using NUOPC_CompDerive() 
   recursive subroutine NUOPC_CplCompDerive(comp, genericSetServicesRoutine, rc)
ARGUMENTS:
     type(ESMF_CplComp),  intent(in)            :: comp
     interface
       subroutine genericSetServicesRoutine(cplcomp, rc)
         use ESMF
         implicit none
         type(ESMF_CplComp)         :: cplcomp  ! must not be optional
         integer, intent(out)       :: rc       ! must not be optional
       end subroutine
     end interface
     integer,             intent(out), optional :: rc
DESCRIPTION:

Derive a CplComp (i.e. Connector) from a generic component by calling into the specified SetServices() routine of the generic component. This is typically the first call in the SetServices() routine of the specializing component, and is followed by NUOPC_CompSpecialize() calls.

3.6.30 NUOPC_CompFilterPhaseMap - Filter the Phase Map of a GridComp


INTERFACE:

   ! Private name; call using NUOPC_CompFilterPhaseMap()
   subroutine NUOPC_GridCompFilterPhaseMap(comp, methodflag, acceptStringList, &
     rc)
ARGUMENTS:
     type(ESMF_GridComp)                           :: comp
     type(ESMF_Method_Flag), intent(in)            :: methodflag
     character(len=*),       intent(in)            :: acceptStringList(:)
     integer,                intent(out), optional :: rc
DESCRIPTION:

Filter all PhaseMap entries in a GridComp (i.e. Model, Mediator, or Driver) that do not match any entry in the acceptStringList.

3.6.31 NUOPC_CompFilterPhaseMap - Filter the Phase Map of a CplComp


INTERFACE:

   ! Private name; call using NUOPC_CompFilterPhaseMap()
   subroutine NUOPC_CplCompFilterPhaseMap(comp, methodflag, acceptStringList, &
     rc)
ARGUMENTS:
     type(ESMF_CplComp)                            :: comp
     type(ESMF_Method_Flag), intent(in)            :: methodflag
     character(len=*),       intent(in)            :: acceptStringList(:)
     integer,                intent(out), optional :: rc
DESCRIPTION:

Filter all PhaseMap entries in a CplComp (i.e. Connector) that do not match any entry in the acceptStringList.

3.6.32 NUOPC_CompGet - Access info from GridComp


INTERFACE:

   ! Private name; call using NUOPC_CompGet()
   subroutine NUOPC_GridCompGet(comp, name, verbosity, profiling, diagnostic, rc)
ARGUMENTS:
     type(ESMF_GridComp)                       :: comp
     character(len=*),   intent(out), optional :: name
     integer,            intent(out), optional :: verbosity
     integer,            intent(out), optional :: profiling
     integer,            intent(out), optional :: diagnostic
     integer,            intent(out), optional :: rc
DESCRIPTION:

Access information from a GridComp.

3.6.33 NUOPC_CompGet - Access info from CplComp


INTERFACE:

   ! Private name; call using NUOPC_CompGet()
   subroutine NUOPC_CplCompGet(comp, name, verbosity, profiling, diagnostic, rc)
ARGUMENTS:
     type(ESMF_CplComp)                        :: comp
     character(len=*),   intent(out), optional :: name
     integer,            intent(out), optional :: verbosity
     integer,            intent(out), optional :: profiling
     integer,            intent(out), optional :: diagnostic
     integer,            intent(out), optional :: rc
DESCRIPTION:

Access information from a CplComp.

3.6.34 NUOPC_CompSearchPhaseMap - Search the Phase Map of a GridComp


INTERFACE:

   ! Private name; call using NUOPC_CompSearchPhaseMap()
   subroutine NUOPC_GridCompSearchPhaseMap(comp, methodflag, internalflag, &
     phaseLabel, phaseIndex, rc)
ARGUMENTS:
     type(ESMF_GridComp)                           :: comp
     type(ESMF_Method_Flag), intent(in)            :: methodflag
     logical,                intent(in),  optional :: internalflag
     character(len=*),       intent(in),  optional :: phaseLabel
     integer,                intent(out)           :: phaseIndex
     integer,                intent(out), optional :: rc
DESCRIPTION:

Search all PhaseMap entries in a GridComp (i.e. Model, Mediator, or Driver) to see if phaseLabel is found. Return the associated ESMF phaseIndex, or -1 if not found. If phaseLabel is not specified, set phaseIndex to the first entry in the PhaseMap, or -1 if there are no entries. The internalflag argument allows to search the internal phase maps of driver components. The default is internalflag=.false..

3.6.35 NUOPC_CompSearchPhaseMap - Search the Phase Map of a CplComp


INTERFACE:

   ! Private name; call using NUOPC_CompSearchPhaseMap()
   subroutine NUOPC_CplCompSearchPhaseMap(comp, methodflag, phaseLabel, &
     phaseIndex, rc)
ARGUMENTS:
     type(ESMF_CplComp)                            :: comp
     type(ESMF_Method_Flag), intent(in)            :: methodflag
     character(len=*),       intent(in),  optional :: phaseLabel
     integer,                intent(out)           :: phaseIndex
     integer,                intent(out), optional :: rc
DESCRIPTION:

Search all PhaseMap entries in a CplComp (i.e. Connector) to see if phaseLabel is found. Return the associated ESMF phaseIndex, or -1 if not found. If phaseLabel is not specified, set phaseIndex to the first entry in the PhaseMap, or -1 if there are no entries.

3.6.36 NUOPC_CompSearchRevPhaseMap - Reverse Search the Phase Map of a GridComp


INTERFACE:

   ! Private name; call using NUOPC_CompSearchRevPhaseMap()
   subroutine NUOPC_GridCompSearchRevPhaseMap(comp, methodflag, internalflag, &
     phaseIndex, phaseLabel, rc)
ARGUMENTS:
     type(ESMF_GridComp)                           :: comp
     type(ESMF_Method_Flag), intent(in)            :: methodflag
     logical,                intent(in),  optional :: internalflag
     integer,                intent(in),  optional :: phaseIndex
     character(len=*),       intent(out)           :: phaseLabel
     integer,                intent(out), optional :: rc
DESCRIPTION:

Search all PhaseMap entries in a GridComp (i.e. Model, Mediator, or Driver) to see if the ESMF phaseIndex is found. Return the associated phaseLabel, or an empty string if not found. If phaseIndex is not specified, set phaseLabel to the first entry in the PhaseMap, or an empty string if there are no entries. The internalflag argument allows to search the internal phase maps of driver components. The default is internalflag=.false..

3.6.37 NUOPC_CompSearchRevPhaseMap - Reverse Search the Phase Map of a CplComp


INTERFACE:

   ! Private name; call using NUOPC_CompSearchRevPhaseMap()
   subroutine NUOPC_CplCompSearchRevPhaseMap(comp, methodflag, phaseIndex, &
     phaseLabel, rc)
ARGUMENTS:
     type(ESMF_CplComp)                            :: comp
     type(ESMF_Method_Flag), intent(in)            :: methodflag
     integer,                intent(in),  optional :: phaseIndex
     character(len=*),       intent(out)           :: phaseLabel
     integer,                intent(out), optional :: rc
DESCRIPTION:

Search all PhaseMap entries in a CplComp (i.e. Connector) to see if the ESMF phaseIndex is found. Return the associated phaseLabel, or an empty string if not found. If phaseIndex is not specified, set phaseLabel to the first entry in the PhaseMap, or an empty string if there are no entries.

3.6.38 NUOPC_CompSetClock - Initialize and set the internal Clock of a GridComp


INTERFACE:

   ! Private name; call using NUOPC_CompSetClock()
   subroutine NUOPC_GridCompSetClock(comp, externalClock, stabilityTimeStep, rc)
ARGUMENTS:
     type(ESMF_GridComp),     intent(inout)         :: comp
     type(ESMF_Clock),        intent(in)            :: externalClock
     type(ESMF_TimeInterval), intent(in),  optional :: stabilityTimeStep
     integer,                 intent(out), optional :: rc
DESCRIPTION:

Set the component internal clock as a copy of externalClock, but with a timeStep that is less than or equal to the stabilityTimeStep. At the same time ensure that the timeStep of the external clock is a multiple of the timeStep of the internal clock. If the stabilityTimeStep argument is not provided then the internal clock will simply be set as a copy of the external clock.

3.6.39 NUOPC_CompSetEntryPoint - Set entry point for a GridComp (DEPRECATED!)


INTERFACE:

   ! Private name; call using NUOPC_CompSetEntryPoint()
   subroutine NUOPC_GridCompSetEntryPoint(comp, methodflag, phaseLabelList, &
     userRoutine, rc)
ARGUMENTS:
     type(ESMF_GridComp)                     :: comp
     type(ESMF_Method_Flag), intent(in)      :: methodflag
     character(len=*),       intent(in)      :: phaseLabelList(:)
     interface
       subroutine userRoutine(gridcomp, importState, exportState, clock, rc)
         use ESMF_CompMod
         use ESMF_StateMod
         use ESMF_ClockMod
         implicit none
         type(ESMF_GridComp)         :: gridcomp     ! must not be optional
         type(ESMF_State)            :: importState  ! must not be optional
         type(ESMF_State)            :: exportState  ! must not be optional
         type(ESMF_Clock)            :: clock        ! must not be optional
         integer, intent(out)        :: rc           ! must not be optional
       end subroutine
     end interface
     integer,          intent(out), optional :: rc
DESCRIPTION:

Set an entry point for a GridComp (i.e. Model, Mediator, or Driver). Publish the new entry point in the correct PhaseMap component attribute.

Starting with version 8.1.0, the use of this method is deprecated. Components should instead specialize exclusively using the NUOPC_CompSpecialize() method.

3.6.40 NUOPC_CompSetEntryPoint - Set entry point for a CplComp (DEPRECATED!)


INTERFACE:

   ! Private name; call using NUOPC_CompSetEntryPoint()
   subroutine NUOPC_CplCompSetEntryPoint(comp, methodflag, phaseLabelList, &
     userRoutine, rc)
ARGUMENTS:
     type(ESMF_CplComp)                      :: comp
     type(ESMF_Method_Flag), intent(in)      :: methodflag
     character(len=*),       intent(in)      :: phaseLabelList(:)
     interface
       subroutine userRoutine(cplcomp, importState, exportState, clock, rc)
         use ESMF_CompMod
         use ESMF_StateMod
         use ESMF_ClockMod
         implicit none
         type(ESMF_CplComp)          :: cplcomp      ! must not be optional
         type(ESMF_State)            :: importState  ! must not be optional
         type(ESMF_State)            :: exportState  ! must not be optional
         type(ESMF_Clock)            :: clock        ! must not be optional
         integer, intent(out)        :: rc           ! must not be optional
       end subroutine
     end interface
     integer,          intent(out), optional :: rc
DESCRIPTION:

Set an entry point for a CplComp (i.e. Connector). Publish the new entry point in the correct PhaseMap component attribute.

Starting with version 8.1.0, the use of this method is deprecated. Components should instead specialize exclusively using the NUOPC_CompSpecialize() method.

3.6.41 NUOPC_CompSetInternalEntryPoint - Set internal entry point for a GridComp


INTERFACE:

   ! Private name; call using NUOPC_CompSetInternalEntryPoint()
   subroutine NUOPC_GridCompSetIntEntryPoint(comp, methodflag, phaseLabelList, &
     userRoutine, rc)
ARGUMENTS:
     type(ESMF_GridComp)                     :: comp
     type(ESMF_Method_Flag), intent(in)      :: methodflag
     character(len=*),       intent(in)      :: phaseLabelList(:)
     interface
       subroutine userRoutine(gridcomp, importState, exportState, clock, rc)
         use ESMF_CompMod
         use ESMF_StateMod
         use ESMF_ClockMod
         implicit none
         type(ESMF_GridComp)         :: gridcomp     ! must not be optional
         type(ESMF_State)            :: importState  ! must not be optional
         type(ESMF_State)            :: exportState  ! must not be optional
         type(ESMF_Clock)            :: clock        ! must not be optional
         integer, intent(out)        :: rc           ! must not be optional
       end subroutine
     end interface
     integer,          intent(out), optional :: rc
DESCRIPTION:

Set an internal entry point for a GridComp (i.e. Driver). Only Drivers currently utilize internal entry points. Internal entry points allow user specialization on the driver level during initialization and run sequencing.

3.6.42 NUOPC_CompSetServices - Try to find and call SetServices in a shared object


INTERFACE:

   ! Private name; call using NUOPC_CompSetServices()
   recursive subroutine NUOPC_GridCompSetServices(comp, sharedObj, userRc, rc)
ARGUMENTS:
     type(ESMF_GridComp),     intent(inout)         :: comp
     character(len=*),        intent(in),  optional :: sharedObj
     integer,                 intent(out), optional :: userRc
     integer,                 intent(out), optional :: rc
DESCRIPTION:

Try to find a routine called "SetServices" in the sharedObj file and execute the routine. An attempt is made to find a routine that is close in name to "SetServices", allowing for compiler name mangling, i.e. upper and lower case, as well as trailing underscores.

3.6.43 NUOPC_CompSetVM - Try to find and call SetVM in a shared object


INTERFACE:

   ! Private name; call using NUOPC_CompSetVM()
   recursive subroutine NUOPC_GridCompSetVM(comp, sharedObj, userRc, rc)
ARGUMENTS:
     type(ESMF_GridComp),     intent(inout)         :: comp
     character(len=*),        intent(in),  optional :: sharedObj
     integer,                 intent(out), optional :: userRc
     integer,                 intent(out), optional :: rc
DESCRIPTION:

Try to find a routine called "SetVM" in the sharedObj file and execute the routine. An attempt is made to find a routine that is close in name to "SetVM", allowing for compiler name mangling, i.e. upper and lower case, as well as trailing underscores.

3.6.44 NUOPC_CompSpecialize - Specialize a derived GridComp


INTERFACE:

   ! Private name; call using NUOPC_CompSpecialize()
   subroutine NUOPC_GridCompSpecialize(comp, specLabel, specPhaseLabel, &
     specRoutine, rc)
ARGUMENTS:
     type(ESMF_GridComp)                     :: comp
     character(len=*), intent(in)            :: specLabel
     character(len=*), intent(in),  optional :: specPhaseLabel
     interface
       subroutine specRoutine(gridcomp, rc)
         use ESMF
         implicit none
         type(ESMF_GridComp)        :: gridcomp ! must not be optional
         integer, intent(out)       :: rc       ! must not be optional
       end subroutine
     end interface
     integer,          intent(out), optional :: rc
DESCRIPTION:

Specialize a derived GridComp (i.e. Model, Mediator, or Driver). If specPhaseLabel is specified, the specialization only applies to the associated phase. Otherwise the specialization applies to all phases.

3.6.45 NUOPC_CompSpecialize - Specialize a derived CplComp


INTERFACE:

   ! Private name; call using NUOPC_CompSpecialize()
   subroutine NUOPC_CplCompSpecialize(comp, specLabel, specPhaseLabel, &
     specRoutine, rc)
ARGUMENTS:
     type(ESMF_CplComp)                      :: comp
     character(len=*), intent(in)            :: specLabel
     character(len=*), intent(in),  optional :: specPhaseLabel
     interface
       subroutine specRoutine(cplcomp, rc)
         use ESMF
         implicit none
         type(ESMF_CplComp)         :: cplcomp  ! must not be optional
         integer, intent(out)       :: rc       ! must not be optional
       end subroutine
     end interface
     integer,          intent(out), optional :: rc
DESCRIPTION:

Specialize a derived CplComp (i.e. Connector). If specPhaseLabel is specified, the specialization only applies to the associated phase. Otherwise the specialization applies to all phases.

3.7 Field Dictionary Methods

3.7.1 NUOPC_FieldDictionaryAddEntry - Add an entry to the NUOPC Field dictionary


INTERFACE:

   subroutine NUOPC_FieldDictionaryAddEntry(standardName, canonicalUnits, rc)
ARGUMENTS:
     character(*),                 intent(in)            :: standardName
     character(*),                 intent(in)            :: canonicalUnits
     integer,                      intent(out), optional :: rc
DESCRIPTION:

Add an entry to the NUOPC Field dictionary. If necessary the dictionary is first set up.

3.7.2 NUOPC_FieldDictionaryEgest - Egest NUOPC Field dictionary into FreeFormat


INTERFACE:

   subroutine NUOPC_FieldDictionaryEgest(freeFormat, iofmt, rc)
ARGUMENTS:
     type(NUOPC_FreeFormat), intent(out)           :: freeFormat
     type(ESMF_IOFmt_Flag),  intent(in),  optional :: iofmt
     integer,                intent(out), optional :: rc
DESCRIPTION:

Egest the contents of the NUOPC Field dictionary into a FreeFormat object. If I/O format option iofmt is provided and equal to ESMF_IOFMT_YAML, the FreeFormat object will contain the NUOPC Field dictionary expressed in YAML format. Other values for iofmt are ignored and this method behaves as if the optional iofmt argument were missing. In such a case, freeFormat will contain NUOPC Field dictionary entries in the traditional format. It is the caller's responsibility to destroy the created freeFormat object.

3.7.3 NUOPC_FieldDictionaryGetEntry - Get information about a NUOPC Field dictionary entry


INTERFACE:

   subroutine NUOPC_FieldDictionaryGetEntry(standardName, canonicalUnits, rc)
ARGUMENTS:
     character(*),                 intent(in)            :: standardName
     character(*),                 intent(out), optional :: canonicalUnits
     integer,                      intent(out), optional :: rc
DESCRIPTION:

Return the canonical units that the NUOPC Field dictionary associates with the standardName.

3.7.4 NUOPC_FieldDictionaryHasEntry - Check whether the NUOPC Field dictionary has a specific entry


INTERFACE:

   function NUOPC_FieldDictionaryHasEntry(standardName, rc)
RETURN VALUE:
     logical :: NUOPC_FieldDictionaryHasEntry
ARGUMENTS:
     character(*),                 intent(in)            :: standardName
     integer,                      intent(out), optional :: rc
DESCRIPTION:

Return .true. if the NUOPC Field dictionary has an entry with the specified standardName, .false. otherwise.

3.7.5 NUOPC_FieldDictionaryMatchSyno - Check whether the NUOPC Field dictionary considers the standard names synonyms


INTERFACE:

   function NUOPC_FieldDictionaryMatchSyno(standardName1, standardName2, rc)
RETURN VALUE:
     logical :: NUOPC_FieldDictionaryMatchSyno
ARGUMENTS:
     character(*),                 intent(in)            :: standardName1
     character(*),                 intent(in)            :: standardName2
     integer,                      intent(out), optional :: rc
DESCRIPTION:

Return .true. if the NUOPC Field dictionary considers standardName1 and standardName2 synonyms, .false. otherwise. Also, if standardName1 and/or standardName2 do not correspond to an existing dictionary entry, .false. will be returned.

3.7.6 NUOPC_FieldDictionarySetSyno - Set synonyms in the NUOPC Field dictionary


INTERFACE:

   subroutine NUOPC_FieldDictionarySetSyno(standardNames, rc)
ARGUMENTS:
     character(*),                 intent(in)            :: standardNames(:)
     integer,                      intent(out), optional :: rc
DESCRIPTION:

Set all of the elements of the standardNames argument to be considered synonyms by the field dictionary. Every element in standardNames must correspond to the standard name of already existing entries in the field dictionary, or else an error will be returned.

3.7.7 NUOPC_FieldDictionarySetup - Setup the default NUOPC Field dictionary


INTERFACE:

   ! Private name; call using NUOPC_FieldDictionarySetup()
   subroutine NUOPC_FieldDictionarySetupDefault(rc)
ARGUMENTS:
     integer,      intent(out), optional   :: rc
DESCRIPTION:

Setup the default NUOPC Field dictionary.

3.7.8 NUOPC_FieldDictionarySetup - Setup the NUOPC Field dictionary from YAML file


INTERFACE:

   ! Private name; call using NUOPC_FieldDictionarySetup()
   subroutine NUOPC_FieldDictionarySetupFile(fileName, rc)
ARGUMENTS:
     character(len=*),      intent(in)              :: fileName
     integer,               intent(out), optional   :: rc
DESCRIPTION:

Setup the NUOPC Field dictionary by reading its content from YAML file. If the NUOPC Field dictionary already exists, remove it and create a new one. This feature requires ESMF built with YAML support. Please see the ESMF User's Guide for details.

3.8 Free Format Methods

3.8.1 NUOPC_FreeFormatAdd - Add lines to a FreeFormat object


INTERFACE:

   subroutine NUOPC_FreeFormatAdd(freeFormat, stringList, line, rc)
ARGUMENTS:
     type(NUOPC_FreeFormat),           intent(inout) :: freeFormat
     character(len=*),                 intent(in)    :: stringList(:)
     integer,                optional, intent(in)    :: line
     integer,                optional, intent(out)   :: rc
DESCRIPTION:

Add lines to a FreeFormat object. The capacity of freeFormat may increase during this operation. The new lines provided in stringList are added starting at position line. If line is greater than the current lineCount of freeFormat, blank lines are inserted to fill the gap. By default, i.e. without specifying the line argument, the elements in stringList are added to the end of the freeFormat object.

3.8.2 NUOPC_FreeFormatCreate - Create a FreeFormat object


INTERFACE:

   ! Private name; call using NUOPC_FreeFormatCreate()
   function NUOPC_FreeFormatCreateDefault(freeFormat, stringList, capacity, rc)
RETURN VALUE:
     type(NUOPC_FreeFormat) :: NUOPC_FreeFormatCreateDefault
ARGUMENTS:
     type(NUOPC_FreeFormat), optional, intent(in)  :: freeFormat
     character(len=*),       optional, intent(in)  :: stringList(:)
     integer,                optional, intent(in)  :: capacity
     integer,                optional, intent(out) :: rc
DESCRIPTION:

Create a new FreeFormat object, which by default is empty. If freeFormat is provided, then the newly created object starts as a copy of freeFormat. If stringList is provided, then it is added to the end of the newly created object. If capacity is provided, it is used for the initial creation of the newly created FreeFormat object. However, if the freeFormat or stringList arguments are present, the final capacity may be larger than specified by capacity.

3.8.3 NUOPC_FreeFormatCreate - Create a FreeFormat object from Config


INTERFACE:

   ! Private name; call using NUOPC_FreeFormatCreate()
   function NUOPC_FreeFormatCreateRead(config, label, relaxedflag, rc)
RETURN VALUE:
     type(NUOPC_FreeFormat) :: NUOPC_FreeFormatCreateRead
ARGUMENTS:
     type(ESMF_Config)                            :: config
     character(len=*),      intent(in)            :: label
     logical,               intent(in),  optional :: relaxedflag
     integer,               intent(out), optional :: rc
DESCRIPTION:

Create a new FreeFormat object from ESMF_Config object. The config object must exist, and label must reference either a single line or a table attribute within config. The content of the attribute is read and returned in the newly created FreeFormat object.

By default an error is returned if label is not found in config. This error can be suppressed by setting relaxedflag=.true., in which case an empty FreeFormat object is returned.

3.8.4 NUOPC_FreeFormatDestroy - Destroy a FreeFormat object


INTERFACE:

   subroutine NUOPC_FreeFormatDestroy(freeFormat, rc)
ARGUMENTS:
     type(NUOPC_FreeFormat),           intent(inout) :: freeFormat
     integer,                optional, intent(out)   :: rc
DESCRIPTION:

Destroy a FreeFormat object. All internal memory is deallocated.

3.8.5 NUOPC_FreeFormatGet - Get information from a FreeFormat object


INTERFACE:

   subroutine NUOPC_FreeFormatGet(freeFormat, lineCount, capacity, stringList, rc)
ARGUMENTS:
     type(NUOPC_FreeFormat),                       intent(in)  :: freeFormat
     integer,                            optional, intent(out) :: lineCount
     integer,                            optional, intent(out) :: capacity
     character(len=NUOPC_FreeFormatLen), optional, pointer     :: stringList(:)
     integer,                            optional, intent(out) :: rc
DESCRIPTION:

Get information from a FreeFormat object.

3.8.6 NUOPC_FreeFormatGetLine - Get line info from a FreeFormat object


INTERFACE:

   subroutine NUOPC_FreeFormatGetLine(freeFormat, line, commentChar, lineString, &
     tokenCount, tokenList, rc)
ARGUMENTS:
     type(NUOPC_FreeFormat),                       intent(in)  :: freeFormat
     integer,                                      intent(in)  :: line
     character,                          optional, intent(in)  :: commentChar
     character(len=NUOPC_FreeFormatLen), optional, intent(out) :: lineString
     integer,                            optional, intent(out) :: tokenCount
     character(len=NUOPC_FreeFormatLen), optional, intent(out) :: tokenList(:)
     integer,                            optional, intent(out) :: rc
DESCRIPTION:

Get information about a specific line in a FreeFormat object. If commentChar is specified, anything on the line, starting with commentChar is considered a comment, and subsequently ignored.

3.8.7 NUOPC_FreeFormatLog - Write a FreeFormat object to the default Log


INTERFACE:

   subroutine NUOPC_FreeFormatLog(freeFormat, rc)
ARGUMENTS:
     type(NUOPC_FreeFormat),           intent(in)    :: freeFormat
     integer,                optional, intent(out)   :: rc
DESCRIPTION:

Write a FreeFormat object to the default Log.

3.8.8 NUOPC_FreeFormatPrint - Print a FreeFormat object


INTERFACE:

   subroutine NUOPC_FreeFormatPrint(freeFormat, rc)
ARGUMENTS:
     type(NUOPC_FreeFormat),           intent(in)    :: freeFormat
     integer,                optional, intent(out)   :: rc
DESCRIPTION:

Print a FreeFormat object.

3.9 Utility Routines

3.9.1 NUOPC_AddNamespace - Add a nested state with Namespace to a State


INTERFACE:

   subroutine NUOPC_AddNamespace(state, Namespace, nestedStateName, &
     nestedState, rc)
ARGUMENTS:
     type(ESMF_State), intent(inout)         :: state
     character(len=*), intent(in)            :: Namespace
     character(len=*), intent(in),  optional :: nestedStateName
     type(ESMF_State), intent(out), optional :: nestedState
     integer,          intent(out), optional :: rc
DESCRIPTION:

Add a Namespace to state. Namespaces are implemented via nested states. This creates a nested state inside of state. The nested state is returned as nestedState. If provided, nestedStateName will be used to name the newly created nested state. The default name of the nested state is equal to Namespace.

The arguments are:

state
The ESMF_State object to which the Namespace is added.
Namespace
The Namespace string.
[nestedStateName]
Name of the nested state. Defaults to Namespace.
[nestedState]
Optional return of the newly created nested state.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.2 NUOPC_AddNestedState - Add a nested state to a state with NUOPC attributes


INTERFACE:

   subroutine NUOPC_AddNestedState(state, Namespace, CplSet, nestedStateName, &
     vm, nestedState, rc)
ARGUMENTS:
     type(ESMF_State), intent(inout)         :: state
     character(len=*), intent(in),  optional :: Namespace
     character(len=*), intent(in),  optional :: CplSet
     character(len=*), intent(in),  optional :: nestedStateName
     type(ESMF_VM),    intent(in),  optional :: vm
     type(ESMF_State), intent(out), optional :: nestedState
     integer,          intent(out), optional :: rc
DESCRIPTION:

Create a nested state inside of state. The arguments Namespace and tt CplSet are used to set NUOPC attributes on the newly created state. The nested state is returned as nestedState. If provided, nestedStateName will be used to name the newly created nested state. The default name of the nested state is equal to Namespace_CplSet, Namespace, or CplSet if the arguments are provided.

The arguments are:

state
The ESMF_State object to which the namespace is added.
[Namespace]
Optional The Namespace string. Defaults to "__UNSPECIFIED__".
[CplSet]
Optional The CplSet string. Defaults to "__UNSPECIFIED__".
[nestedStateName]
Name of the nested state. Defaults to Namespace_CplSet, Namespace, or CplSet if arguments are provided.
[vm]
If present, the nested state object is created on the specified ESMF_VM object. The default is to create the nested state object on the VM of the current component context.
[nestedState]
Optional return of the newly created nested state.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.3 NUOPC_Advertise - Advertise a single Field in a State


INTERFACE:

   ! Private name; call using NUOPC_Advertise() 
   subroutine NUOPC_AdvertiseField(state, StandardName, Units, &
     LongName, ShortName, name, TransferOfferGeomObject, SharePolicyField, &
     SharePolicyGeomObject, vm, field, rc)
ARGUMENTS:
     type(ESMF_State), intent(inout)         :: state
     character(*),     intent(in)            :: StandardName
     character(*),     intent(in),  optional :: Units
     character(*),     intent(in),  optional :: LongName
     character(*),     intent(in),  optional :: ShortName
     character(*),     intent(in),  optional :: name
     character(*),     intent(in),  optional :: TransferOfferGeomObject
     character(*),     intent(in),  optional :: SharePolicyField
     character(*),     intent(in),  optional :: SharePolicyGeomObject
     type(ESMF_VM),    intent(in),  optional :: vm
     type(ESMF_Field), intent(out), optional :: field
     integer,          intent(out), optional :: rc
DESCRIPTION:

Advertise a field in a state. This creates an empty field and adds it to state. The "StandardName", "Units", "LongName", "ShortName", and "TransferOfferGeomObject" attributes of the field are set according to the provided input..

The call checks the provided information against the NUOPC Field Dictionary to ensure correctness. Defaults are set according to the NUOPC Field Dictionary.

The arguments are:

state
The ESMF_State object through which the field is advertised.
StandardName
The "StandardName" attribute of the advertised field. Must be a StandardName found in the NUOPC Field Dictionary.
NOTE that if by below default rules, StandardName is also used as the input for name, then it must not contain the slash ("/") character.
[Units]
The "Units" attribute of the advertised field. Must be convertible to the canonical units specified in the NUOPC Field Dictionary for the specified StandardName. (Currently this is restricted to be identical to the canonical untis specified in the NUOPC Field Dictionary.) If omitted, the default is to use the canonical units associated with the StandardName in the NUOPC Field Dictionary.
[LongName]
The "LongName" attribute of the advertised field. NUOPC does not restrict the value of this attribute. If omitted, the default is to use the StandardName.
[ShortName]
The "ShortName" attribute of the advertised field. NUOPC does not restrict the value of this attribute. If omitted, the default is to use the StandardName.
NOTE that if by below default rules, ShortName is also used as the input for name, then it must not contain the slash ("/") character.
[name]
The actual name of the advertised field by which it is accessed in the state object. The string provided for name must not contain the slash ("/") character. If omitted, the default is to use the value of the ShortName.
[TransferOfferGeomObject]
If the state intent of state is ESMF_STATEINTENT_EXPORT, the "ProducerTransferOffer" attribute of the advertised field is set. If the state intent of state is ESMF_STATEINTENT_IMPORT, the "ConsumerTransferOffer" attribute of the advertised field is set. NUOPC controls the vocabulary of this attribute. Valid options are "will provide", "can provide", "cannot provide". If omitted, the default is "will provide".
[SharePolicyField]
The "SharePolicyField" attribute of the advertised field. NUOPC controls the vocabulary of this attribute. Valid options are "share", and "not share". If omitted, the default is "not share".
[SharePolicyGeomObject]
The "SharePolicyGeomObject" attribute of the advertised field. NUOPC controls the vocabulary of this attribute. Valid options are "share", and "not share". If omitted, the default is equal to SharePolicyField.
[vm]
If present, the Field object used during advertising is created on the specified ESMF_VM object. The default is to create the Field object on the VM of the current component context.
[field]
Returns the empty field object that was used to advertise.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.4 NUOPC_Advertise - Advertise a list of Fields in a State


INTERFACE:

   ! Private name; call using NUOPC_Advertise() 
   subroutine NUOPC_AdvertiseFields(state, StandardNames, &
     TransferOfferGeomObject, SharePolicyField, SharePolicyGeomObject, rc)
ARGUMENTS:
     type(ESMF_State), intent(inout)         :: state
     character(*),     intent(in)            :: StandardNames(:)
     character(*),     intent(in),  optional :: TransferOfferGeomObject
     character(*),     intent(in),  optional :: SharePolicyField
     character(*),     intent(in),  optional :: SharePolicyGeomObject
     integer,          intent(out), optional :: rc
DESCRIPTION:

Advertise a list of fields in a state. This creates a list of empty fields and adds it to the state. The "StandardName", "TransferOfferGeomObject", "SharePolicyField", and "SharePolicyGeomObject" attributes of all the fields are set according to the provided input. The "Units", "LongName", and "ShortName" attributes for each field are set according to the defaults documented under method 3.9.3

The call checks the provided information against the NUOPC Field Dictionary to ensure correctness.

The arguments are:

state
The ESMF_State object through which the fields are advertised.
StandardNames
A list of "StandardName" attributes of the advertised fields. Must be StandardNames found in the NUOPC Field Dictionary.
[TransferOfferGeomObject]
The "TransferOfferGeomObject" attribute of the advertised fields. This setting applies to all the fields advertised in this call. NUOPC controls the vocabulary of this attribute. Valid options are "will provide", "can provide", "cannot provide". If omitted, the default is "will provide".
[SharePolicyField]
The "SharePolicyField" attribute of the advertised fields. This setting applies to all the fields advertised in this call. NUOPC controls the vocabulary of this attribute. Valid options are "share", and "not share". If omitted, the default is "not share".
[SharePolicyGeomObject]
The "SharePolicyGeomObject" attribute of the advertised fields. This setting applies to all the fields advertised in this call. NUOPC controls the vocabulary of this attribute. Valid options are "share", and "not share". If omitted, the default is equal to SharePolicyField.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.5 NUOPC_AdjustClock - Adjust the timestep in a clock


INTERFACE:

   subroutine NUOPC_AdjustClock(clock, maxTimestep, rc)
ARGUMENTS:
     type(ESMF_Clock)                               :: clock
     type(ESMF_TimeInterval), intent(in),  optional :: maxTimestep
     integer,                 intent(out), optional :: rc
DESCRIPTION:

Adjust the clock to have a potentially smaller timestep. The timestep on the incoming clock object is compared to the maxTimestep, and reset to the smaller of the two.

The arguments are:

clock
The clock to be adjusted.
[maxTimestep]
Upper bound of the timestep allowed in clock.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.6 NUOPC_CheckSetClock - Check a Clock for compatibility and set its values


INTERFACE:

   subroutine NUOPC_CheckSetClock(setClock, checkClock, setStartTimeToCurrent, &
     currTime, forceCurrTime, checkTimeStep, forceTimeStep, rc)
ARGUMENTS:
     type(ESMF_Clock),        intent(inout)         :: setClock
     type(ESMF_Clock),        intent(in)            :: checkClock
     logical,                 intent(in),  optional :: setStartTimeToCurrent
     type(ESMF_Time),         intent(in),  optional :: currTime
     logical,                 intent(in),  optional :: forceCurrTime
     logical,                 intent(in),  optional :: checkTimeStep
     logical,                 intent(in),  optional :: forceTimeStep
     integer,                 intent(out), optional :: rc
DESCRIPTION:

By default compare setClock to checkClock to ensure they match in their current time. Further ensure that the timeStep of checkClock is a multiple of the timeStep of setClock. If both conditions are satisfied then the stopTime of the setClock is set one checkClock timeStep, or setClock runDuration, ahead of the current time, which ever is shorter. The direction of checkClock is considered when setting the stopTime.

By default the startTime of the setClock is not modified. However, if setStartTimeToCurrent == .true. the startTime of setClock is set to the currentTime of checkClock.

The arguments are:

setClock
The ESMF_Clock object to be checked and set.
checkClock
The reference clock object.
[setStartTimeToCurrent]
If .true. then also set the startTime in setClock according to the startTime in checkClock. The default is .false..
[currTime]
If provided, use currTime instead of checkClock when checking or setting the current time of setClock.
[forceCurrTime]
If .true. then do not check the current time of the setClock, but instead force it to align with the checkClock, or currTime, if it was provided. The default is .false..
[checkTimeStep]
If .true. then check that timeStep of the setClock can reach the next increment on the checkClock by an integer number of steps. For .false. do not check this condition. The default is .true..
[forceTimeStep]
If .true. then do not use the timeStep of the setClock to check if the next increment on the checkClock can be reached in an integer number of steps. Instead set the timeStep of the setClock to the timeStep of the checkClock. The default is .false..
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.7 NUOPC_GetAttribute - Get the value of a NUOPC Field Attribute


INTERFACE:

   ! Private name; call using NUOPC_GetAttribute()
   subroutine NUOPC_GetAttributeFieldVal(field, name, value, isPresent, isSet, rc)
ARGUMENTS:
     type(ESMF_Field), intent(in)            :: field
     character(*),     intent(in)            :: name
     character(*),     intent(out)           :: value
     logical,          intent(out), optional :: isPresent
     logical,          intent(out), optional :: isSet
     integer,          intent(out), optional :: rc
DESCRIPTION:

Access the attribute name inside of field using the convention NUOPC and purpose Instance.

Unless isPresent and isSet are provided, return with error if the Attribute is not present or not set, respectively.

The arguments are:

field
The ESMF_Field object to be queried.
name
The name of the queried attribute.
value
The value of the queried attribute.
[isPresent]
Set to .true. if the queried attribute is present, .false. otherwise.
[isSet]
Set to .true. if the queried attribute is set, .false. otherwise.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.8 NUOPC_GetAttribute - Get the typekind of a NUOPC Field Attribute


INTERFACE:

   ! Private name; call using NUOPC_GetAttribute()
   subroutine NUOPC_GetAttributeFieldTK(field, name, isPresent, isSet, &
     itemCount, typekind, rc)
ARGUMENTS:
     type(ESMF_Field),         intent(in)            :: field
     character(*),             intent(in)            :: name
     logical,                  intent(out), optional :: isPresent
     logical,                  intent(out), optional :: isSet
     integer,                  intent(out), optional :: itemCount
     type(ESMF_TypeKind_Flag), intent(out), optional :: typekind
     integer,                  intent(out), optional :: rc
DESCRIPTION:

Query the typekind of the attribute name inside of field using the convention NUOPC and purpose Instance.

Unless isPresent and isSet are provided, return with error if the Attribute is not present or not set, respectively.

The arguments are:

field
The ESMF_Field object to be queried.
name
The name of the queried attribute.
[isPresent]
Set to .true. if the queried attribute is present, .false. otherwise.
[isSet]
Set to .true. if the queried attribute is set, .false. otherwise.
[itemCount]
Number of items in the attribute. Return 0 if not present or not set.
[typekind]
The typekind of the queried attribute.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.9 NUOPC_GetAttribute - Get the value of a NUOPC State Attribute


INTERFACE:

   ! Private name; call using NUOPC_GetAttribute()
   subroutine NUOPC_GetAttributeState(state, name, value, isPresent, isSet, &
     itemCount, typekind, rc)
ARGUMENTS:
     type(ESMF_State),         intent(in)            :: state
     character(*),             intent(in)            :: name
     character(*),             intent(out), optional :: value
     logical,                  intent(out), optional :: isPresent
     logical,                  intent(out), optional :: isSet
     integer,                  intent(out), optional :: itemCount
     type(ESMF_TypeKind_Flag), intent(out), optional :: typekind
     integer,                  intent(out), optional :: rc
DESCRIPTION:

Access the attribute name inside of state using the convention NUOPC and purpose Instance. Returns with error if the attribute is not present or not set.

Unless isPresent and isSet are provided, return with error if the Attribute is not present or not set, respectively.

The arguments are:

state
The ESMF_State object to be queried.
name
The name of the queried attribute.
[value]
The value of the queried attribute.
[isPresent]
Set to .true. if the queried attribute is present, .false. otherwise.
[isSet]
Set to .true. if the queried attribute is set, .false. otherwise.
[itemCount]
Number of items in the attribute. Return 0 if not present or not set.
[typekind]
The typekind of the queried attribute.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.10 NUOPC_GetStateMemberLists - Build lists of information of State members


INTERFACE:

   subroutine NUOPC_GetStateMemberLists(state, StandardNameList, &
     ConnectedList, NamespaceList, CplSetList, itemNameList, fieldList, &
     stateList, nestedFlag, rc)
ARGUMENTS:
     type(ESMF_State),       intent(in)            :: state
     character(ESMF_MAXSTR), pointer, optional     :: StandardNameList(:)
     character(ESMF_MAXSTR), pointer, optional     :: ConnectedList(:)
     character(ESMF_MAXSTR), pointer, optional     :: NamespaceList(:)
     character(ESMF_MAXSTR), pointer, optional     :: CplSetList(:)
     character(ESMF_MAXSTR), pointer, optional     :: itemNameList(:)
     type(ESMF_Field),       pointer, optional     :: fieldList(:)
     type(ESMF_State),       pointer, optional     :: stateList(:)
     logical,                intent(in), optional  :: nestedFlag
     integer,                intent(out), optional :: rc
DESCRIPTION:

Construct lists containing the StandardNames, field names, and connected status of the fields in state. Return this information in the list arguments. Recursively parse through nested States.

All pointer arguments present must enter this method unassociated. This means that the user code must explicitly call nullify() or use the => null() syntax on the variables passed in as any of the pointer arguments. On return, the pointer arguments may either be unassociated or associated. Consequently the user code must first check the status of any of the returned pointer arguments via the associated() intrinsic before accessing the argument. The responsibility for deallocation of associated pointer arguments transfers to the caller.

The arguments are:

state
The ESMF_State object to be queried.
[StandardNameList]
If present, return a list of the "StandardName" attribute of each member. See the note about pointer arguments in the description section above for correct usage.
[ConnectedList]
If present, return a list of the "Connected" attribute of each member. See the note about pointer arguments in the description section above for correct usage.
[NamespaceList]
If present, return a list of the "Namespace" attribute of each member. See the note about pointer arguments in the description section above for correct usage.
[CplSetList]
If present, return a list of the "CplSet" attribute of each member. See the note about pointer arguments in the description section above for correct usage.
[itemNameList]
If present, return a list of each member name. See the note about pointer arguments in the description section above for correct usage.
[fieldList]
If present, return a list of the member fields. See the note about pointer arguments in the description section above for correct usage.
[stateList]
If present, return a list of the states corresonding to the owner of the fields returned under fieldList. See the note about pointer arguments in the description section above for correct usage.
[nestedFlag]
When set to .true., returns information from nested States (default). When set to .false., returns information at the current State level only.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.11 NUOPC_GetStateMemberCount - Determing number of State members


INTERFACE:

   subroutine NUOPC_GetStateMemberCount(state, fieldCount, nestedFlag, rc)
ARGUMENTS:
     type(ESMF_State),       intent(in)            :: state
     integer,                intent(out), optional :: fieldCount
     logical,                intent(in),  optional :: nestedFlag
     integer,                intent(out), optional :: rc
DESCRIPTION:

Determine the number of fields in state. By default recursively parse through nested States.

The arguments are:

state
The ESMF_State object to be queried.
[fieldCount]
Number of fields.
[nestedFlag]
When set to .true., returns information from nested States (default). When set to .false., returns information at the current State level only.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.12 NUOPC_GetTimestamp - Get the timestamp of a Field


INTERFACE:

   subroutine NUOPC_GetTimestamp(field, isValid, time, rc)
ARGUMENTS:
     type(ESMF_Field), intent(in)            :: field
     logical,          intent(out), optional :: isValid
     type(ESMF_Time),  intent(out), optional :: time
     integer,          intent(out), optional :: rc
DESCRIPTION:

Access the timestamp on field in form of an ESMF_Time object.

The arguments are:

field
The ESMF_Field object to be checked.
[isValid]
Set to .true. if the timestamp is valid, .false. otherwise.
[time]
The timestamp as ESMF_Time object.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.13 NUOPC_IngestPetList - Ingest a petList from FreeFormat


INTERFACE:

   ! Private name; call using NUOPC_IngestPetList()
   subroutine NUOPC_IngestPetListFF(petList, freeFormat, rc)
ARGUMENTS:
     integer, allocatable,   intent(out)           :: petList(:)
     type(NUOPC_FreeFormat), intent(in),  target   :: freeFormat
     integer,                intent(out), optional :: rc
DESCRIPTION:

Construct a petList from a FreeFormat object.

The arguments are:

petList
The constructed petList. The size and content is set by this method.
freeFormat
The incoming petList information in free format. The format supports two types of elements: Any number of elements may be listed in the free format. The idividual elements are separated by white spaces.

For an example, the free format petList definition

       "2-5 12 0 15-23"
would translate into a petList output of
       (/2, 3, 4, 5, 12, 0, 15, 16, 17, 18, 19, 20, 21, 22, 23/)
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.14 NUOPC_IngestPetList - Ingest a petList from HConfig


INTERFACE:

   ! Private name; call using NUOPC_IngestPetList()
   subroutine NUOPC_IngestPetListHC(petList, hconfig, rc)
ARGUMENTS:
     integer, allocatable,   intent(out)           :: petList(:)
     type(ESMF_HConfig),     intent(in)            :: hconfig
     integer,                intent(out), optional :: rc
DESCRIPTION:

Construct a petList from a HConfig object.

The arguments are:

petList
The constructed petList. The size and content is set by this method.
hconfig
The incoming petList information as HConfig. The provided hconfig must be a scalar, or a list of lists and scalars. The input is recursively processed, and each scalar fed into the FreeFormat version of the NUOPC_IngestPetList() interface as a single string. The resulting petList is the union of all PETs determined by all of the elements contained in hconfig.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.15 NUOPC_IsAtTime - Check if a Field is at the given Time


INTERFACE:

   ! Private name; call using NUOPC_IsAtTime()
   function NUOPC_IsAtTimeField(field, time, rc)
RETURN VALUE:
     logical :: NUOPC_IsAtTimeField
ARGUMENTS:
     type(ESMF_Field), intent(in)            :: field
     type(ESMF_Time),  intent(in)            :: time
     integer,          intent(out), optional :: rc
DESCRIPTION:

Returns .true. if field has a timestamp that matches time. Otherwise returns .false.. On PETs with only a proxy instance of the field, .true. is returned regardless of the actual timestamp.

The arguments are:

field
The ESMF_Field object to be checked.
time
The time to compare against.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.16 NUOPC_IsAtTime - Check if Field(s) in a State are at the given Time


INTERFACE:

   ! Private name; call using NUOPC_IsAtTime()
   function NUOPC_IsAtTimeState(state, time, fieldName, count, fieldList, rc)
RETURN VALUE:
     logical :: NUOPC_IsAtTimeState
ARGUMENTS:
     type(ESMF_State),              intent(in)            :: state
     type(ESMF_Time),               intent(in)            :: time
     character(*),                  intent(in),  optional :: fieldName
     integer,                       intent(out), optional :: count
     type(ESMF_Field), allocatable, intent(out), optional :: fieldList(:)
     integer,                       intent(out), optional :: rc
DESCRIPTION:

Return .true. if the field(s) in state have a timestamp that matches time. Otherwise return .false..

The arguments are:

state
The ESMF_State object to be checked.
time
The time to compare against.
[fieldName]
The name of the field in state to be checked. If provided, and the state does not contain a field with fieldName, return an error in rc. If not provided, check all the fields contained in state and return .true. if all the fields are at the correct time.
[count]
If provided, the number of fields that are at time are returned. If fieldName is present then count cannot be greater than 1.
[fieldList]
If provided, the fields that are not at time are returned. If fieldName is present then fieldList can contain a maximum of 1 field.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.17 NUOPC_IsConnected - Check if a Field is connected


INTERFACE:

   ! Private name; call using NUOPC_IsConnected()
   function NUOPC_IsConnectedField(field, rc)
RETURN VALUE:
     logical :: NUOPC_IsConnectedField
ARGUMENTS:
     type(ESMF_Field), intent(in)            :: field
     integer,          intent(out), optional :: rc
DESCRIPTION:

Return .true. if the field is connected. Otherwise return .false..

The arguments are:

field
The ESMF_Field object to be checked.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.18 NUOPC_IsConnected - Check if Field(s) in a State are connected


INTERFACE:

   ! Private name; call using NUOPC_IsConnected()
   function NUOPC_IsConnectedState(state, fieldName, count, rc)
RETURN VALUE:
     logical :: NUOPC_IsConnectedState
ARGUMENTS:
     type(ESMF_State), intent(in)            :: state
     character(*),     intent(in),  optional :: fieldName
     integer,          intent(out), optional :: count
     integer,          intent(out), optional :: rc
DESCRIPTION:

Return .true. if the field(s) in state are connected. Otherwise return .false..

The arguments are:

state
The ESMF_State object to be checked.
[fieldName]
The name of the field in state to be checked. If provided, and the state does not contain a field with fieldName, return an error in rc. If not provided, check all the fields contained in state and return .true. if all the fields are connected.
[count]
If provided, the number of fields that are connected are returned. If fieldName is present then count cannot be greater than 1.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.19 NUOPC_IsUpdated - Check if a Field is marked as updated


INTERFACE:

   ! Private name; call using NUOPC_IsUpdated()
   function NUOPC_IsUpdatedField(field, rc)
RETURN VALUE:
     logical :: NUOPC_IsUpdatedField
ARGUMENTS:
     type(ESMF_Field), intent(in)            :: field
     integer,          intent(out), optional :: rc
DESCRIPTION:

Return .true. if the field has its "Updated" attribute set to "true". Otherwise return .false..

The arguments are:

field
The ESMF_Field object to be checked.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.20 NUOPC_IsUpdated - Check if Field(s) in a State are marked as updated


INTERFACE:

   ! Private name; call using NUOPC_IsUpdated()
   function NUOPC_IsUpdatedState(state, fieldName, count, rc)
RETURN VALUE:
     logical :: NUOPC_IsUpdatedState
ARGUMENTS:
     type(ESMF_State), intent(in)            :: state
     character(*),     intent(in),  optional :: fieldName
     integer,          intent(out), optional :: count
     integer,          intent(out), optional :: rc
DESCRIPTION:

Return .true. if the field(s) in state have the "Updated" attribute set to "true". Otherwise return .false..

The arguments are:

state
The ESMF_State object to be checked.
[fieldName]
The name of the field in state to be checked. If provided, and the state does not contain a field with fieldName, return an error in rc. If not provided, check all the fields contained in state and return .true. if all the fields are updated.
[count]
If provided, the number of fields that are updated are returned. If fieldName is present then count cannot be greater than 1.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.21 NUOPC_NoOp - No-Operation attachable method for GridComp


INTERFACE:

   subroutine NUOPC_NoOp(gcomp, rc)
ARGUMENTS:
     type(ESMF_GridComp)   :: gcomp
     integer, intent(out)  :: rc
DESCRIPTION:

No-Op method with an interface that matches the requirements for a attachable method for ESMF_GridComp objects.

The arguments are:

gcomp
The ESMF_GridComp object to which this method is attached.
rc
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.22 NUOPC_Realize - Realize previously advertised Fields inside a State on a single Grid with internal allocation


INTERFACE:

   ! Private name; call using NUOPC_Realize()
   subroutine NUOPC_RealizeCompleteG(state, grid, fieldName, typekind, &
     staggerloc, selection, dataFillScheme, field, rc)
ARGUMENTS:
     type(ESMF_State)                                :: state
     type(ESMF_Grid),          intent(in)            :: grid
     character(*),             intent(in),  optional :: fieldName
     type(ESMF_TypeKind_Flag), intent(in),  optional :: typekind
     type(ESMF_StaggerLoc),    intent(in),  optional :: staggerloc
     character(len=*),         intent(in),  optional :: selection
     character(len=*),         intent(in),  optional :: dataFillScheme    
     type(ESMF_Field),         intent(out), optional :: field
     integer,                  intent(out), optional :: rc
DESCRIPTION:

Realize or remove fields inside of state according to selection. All of the fields that are realized are created internally on the same grid object, allocating memory for as many field dimensions as there are grid dimensions.

The type and kind of the created fields is according to argument typekind.

Realized fields are filled with data according to the dataFillScheme argument.

The arguments are:

state
The ESMF_State object in which the fields are realized.
grid
The ESMF_Grid object on which to realize the fields.
[fieldName]
The name of the field in state to be realized, or removed, according to selection. If provided, and the state does not contain a field with name fieldName, return an error in rc. If not provided, realize all the fields contained in state according to selection.
[typekind]
The typekind of the internally created field(s). The valid options are ESMF_TYPEKIND_I4, ESMF_TYPEKIND_I8, ESMF_TYPEKIND_R4, and ESMF_TYPEKIND_R8. By default use the typekind of the partially created field used during advertise, or ESMF_TYPEKIND_R8, if the advertised field did not have a typekind defined.
[staggerloc]
Stagger location of data in grid cells. By default use the same stagger location as the advertising field, or ESMF_STAGGERLOC_CENTER if the advertising field was created empty.
[selection]
Selection of mode of operation:
[dataFillScheme]
Realized fields will be filled according to the selected fill scheme. See ESMF_FieldFill() for fill schemes. Default is to leave the data in realized fields uninitialized.
[field]
Returns the completed field that was realized by this method. This option is only supported if also argument fieldName was specified.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.23 NUOPC_Realize - Realize previously advertised Fields inside a State on a single LocStream with internal allocation


INTERFACE:

   ! Private name; call using NUOPC_Realize()
   subroutine NUOPC_RealizeCompleteLS(state, locstream, fieldName, typekind, selection,&
     dataFillScheme, field, rc)
ARGUMENTS:
     type(ESMF_State)                                :: state
     type(ESMF_LocStream),     intent(in)            :: locstream
     character(*),             intent(in),  optional :: fieldName
     type(ESMF_TypeKind_Flag), intent(in),  optional :: typekind
     character(len=*),         intent(in),  optional :: selection
     character(len=*),         intent(in),  optional :: dataFillScheme    
     type(ESMF_Field),         intent(out), optional :: field
     integer,                  intent(out), optional :: rc
DESCRIPTION:

Realize or remove fields inside of state according to selection. All of the fields that are realized are created internally on the same locstream object, allocating memory accordingly.

The type and kind of the created fields is according to argument typekind.

Realized fields are filled with data according to the dataFillScheme argument.

The arguments are:

state
The ESMF_State object in which the fields are realized.
locstream
The ESMF_LocStream object on which to realize the fields.
[fieldName]
The name of the field in state to be realized, or removed, according to selection. If provided, and the state does not contain a field with name fieldName, return an error in rc. If not provided, realize all the fields contained in state according to selection.
[typekind]
The typekind of the internally created field(s). The valid options are ESMF_TYPEKIND_I4, ESMF_TYPEKIND_I8, ESMF_TYPEKIND_R4, and ESMF_TYPEKIND_R8. By default use the typekind of the partially created field used during advertise, or ESMF_TYPEKIND_R8, if the advertised field did not have a typekind defined.
[selection]
Selection of mode of operation:
[dataFillScheme]
Realized fields will be filled according to the selected fill scheme. See ESMF_FieldFill() for fill schemes. Default is to leave the data in realized fields uninitialized.
[field]
Returns the completed field that was realized by this method. This option is only supported if also argument fieldName was specified.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.24 NUOPC_Realize - Realize previously advertised Fields inside a State on a single Mesh with internal allocation


INTERFACE:

   ! Private name; call using NUOPC_Realize()
   subroutine NUOPC_RealizeCompleteM(state, mesh, fieldName, typekind, &
     meshloc, selection, dataFillScheme, field, rc)
ARGUMENTS:
     type(ESMF_State)                                :: state
     type(ESMF_Mesh),          intent(in)            :: mesh
     character(*),             intent(in),  optional :: fieldName
     type(ESMF_TypeKind_Flag), intent(in),  optional :: typekind
     type(ESMF_MeshLoc),       intent(in),  optional :: meshloc
     character(len=*),         intent(in),  optional :: selection
     character(len=*),         intent(in),  optional :: dataFillScheme
     type(ESMF_Field),         intent(out), optional :: field
     integer,                  intent(out), optional :: rc
DESCRIPTION:

Realize or remove fields inside of state according to selection. All of the fields that are realized are created internally on the same mesh object, allocating memory accordingly.

The type and kind of the created fields is according to argument typekind.

Realized fields are filled with data according to the dataFillScheme argument.

The arguments are:

state
The ESMF_State object in which the fields are realized.
mesh
The ESMF_Mesh object on which to realize the fields.
[fieldName]
The name of the field in state to be realized, or removed, according to selection. If provided, and the state does not contain a field with name fieldName, return an error in rc. If not provided, realize all the fields contained in state according to selection.
[typekind]
The typekind of the internally created field(s). The valid options are ESMF_TYPEKIND_I4, ESMF_TYPEKIND_I8, ESMF_TYPEKIND_R4, and ESMF_TYPEKIND_R8. By default use the typekind of the partially created field used during advertise, or ESMF_TYPEKIND_R8, if the advertised field did not have a typekind defined.
[meshloc]
Location of data in the mesh cell. By default use the same mesh location as the advertising field, or ESMF_STAGGERLOC_NODE if the advertising field was created empty.
[selection]
Selection of mode of operation:
[dataFillScheme]
Realized fields will be filled according to the selected fill scheme. See ESMF_FieldFill() for fill schemes. Default is to leave the data in realized fields uninitialized.
[field]
Returns the completed field that was realized by this method. This option is only supported if also argument fieldName was specified.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.25 NUOPC_Realize - Realize a previously advertised Field in a State


INTERFACE:

   ! Private name; call using NUOPC_Realize()
   subroutine NUOPC_RealizeField(state, field, rc)
ARGUMENTS:
     type(ESMF_State), intent(inout)         :: state
     type(ESMF_Field), intent(in)            :: field
     integer,          intent(out), optional :: rc
DESCRIPTION:

Realize a previously advertised field in state by replacing the advertised field with field of the same name.

The arguments are:

state
The ESMF_State object in which the fields are realized.
field
The new field to put in place of the previously advertised (empty) field.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.26 NUOPC_Realize - Realize a previously advertised Field in a State after Transfer of GeomObject


INTERFACE:

   ! Private name; call using NUOPC_Realize()
   subroutine NUOPC_RealizeTransfer(state, fieldName, typekind, gridToFieldMap, &
     ungriddedLBound, ungriddedUBound, totalLWidth, totalUWidth, &
     realizeOnlyConnected, removeNotConnected, realizeOnlyNotShared, &
     realizeOnlyNotComplete, field, rc)
ARGUMENTS:
     type(ESMF_State)                                :: state
     character(*),             intent(in)            :: fieldName
     type(ESMF_TypeKind_Flag), intent(in),  optional :: typekind
     integer, target,          intent(in),  optional :: gridToFieldMap(:)
     integer, target,          intent(in),  optional :: ungriddedLBound(:)
     integer, target,          intent(in),  optional :: ungriddedUBound(:)
     integer,                  intent(in),  optional :: totalLWidth(:)
     integer,                  intent(in),  optional :: totalUWidth(:)
     logical,                  intent(in),  optional :: realizeOnlyConnected
     logical,                  intent(in),  optional :: removeNotConnected
     logical,                  intent(in),  optional :: realizeOnlyNotShared
     logical,                  intent(in),  optional :: realizeOnlyNotComplete
     type(ESMF_Field),         intent(out), optional :: field
     integer,                  intent(out), optional :: rc
DESCRIPTION:

Realize a field where GeomObject has been set by the NUOPC GeomObject transfer protocol.

The data of the realized field is left uninitialized by this method.

The arguments are:

state
The ESMF_State object in which the field is realized.
fieldName
The name of the field in state to be realized. If state does not contain a field with name fieldName, return an error in rc.
[typekind]
The typekind of the internally created field(s). The valid options are ESMF_TYPEKIND_I4, ESMF_TYPEKIND_I8, ESMF_TYPEKIND_R4, and ESMF_TYPEKIND_R8. By default use the typekind of the connected provider field.
[gridToFieldMap]
The mapping of grid/mesh dimensions against field dimensions. The argument is of rank 1 and with a size of dimCount. The elements correspond to the grid/mesh elements in order, and associates it with the indicated field dimension. Only entries between 1 and the field rank are allowed. There must be no duplicate entries in gridToFieldMap. By default use the gridToFieldMap of the connected provider field.
[ungriddedLBound]
Lower bounds of the ungridded dimensions of the field. The number of elements defines the number of ungridded dimensions of the field and must be consistent with ungriddedUBound. By default use the ungriddedLBound of the connected provider field.
[ungriddedUBound]
Upper bounds of the ungridded dimensions of the field. The number of elements defines the number of ungridded dimensions of the field and must be consistent with ungriddedLBound. By default use the ungriddedLBound of the connected provider field.
[totalLWidth]
This argument is only supported for fields defined on ESMF_Grid. The number elements outside the lower bound of the exclusive region. The argument is of rank 1 and with a size of dimCount, the number of gridded dimensions of the field. The ordering of the dimensions is that of the field (considering gridToFieldMap). By default a zero vector is used, resulting in no elements outside the exclusive region.
[totalUWidth]
This argument is only supported for fields defined on ESMF_Grid. The number elements outside the upper bound of the exclusive region. The argument is of rank 1 and with a size of dimCount, the number of gridded dimensions of the field. The ordering of the dimensions is that of the field (considering gridToFieldMap). By default a zero vector is used, resulting in no elements outside the exclusive region.
[realizeOnlyConnected]
If set to .false., realize the specified field irregardless of the connected status. If set to .true., only a connected field will be realized. The default is .true..
[removeNotConnected]
If set to .false., do not remove a field from the state due to its connected status. If set to .true., remove the field if it is not connected. This requires realizeOnlyConnected to be .true., and a runtime error will be returned otherwise. The default is .true..
[realizeOnlyNotShared]
If set to .false., realize the specified field irregardless of its shared status. If set to .true., only a field that has "ShareStatusField" set to "not shared" will be realized. The default is .true..
[realizeOnlyNotComplete]
If set to .false., realize the specified field irregardless of its complete status. If set to .true., only a field that has not yet been completed will be realized. The default is .true..
[field]
Returns the completed field that was realized by this method. An invalid field object will be returned if the conditions were such that the field was not realized.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.27 NUOPC_SetAttribute - Set the value of a NUOPC Field Attribute


INTERFACE:

   ! Private name; call using NUOPC_SetAttribute()
   subroutine NUOPC_SetAttributeField(field, name, value, rc)
ARGUMENTS:
     type(ESMF_Field)                      :: field
     character(*), intent(in)              :: name
     character(*), intent(in)              :: value
     integer,      intent(out), optional   :: rc
DESCRIPTION:

Set the attribute name inside of field using the convention NUOPC and purpose Instance.

The arguments are:

field
The ESMF_Field object on which to set the attribute.
name
The name of the set attribute.
value
The value of the set attribute.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.28 NUOPC_SetAttribute - Set the value of a NUOPC State Attribute


INTERFACE:

   ! Private name; call using NUOPC_SetAttribute()
   subroutine NUOPC_SetAttributeState(state, name, value, rc)
ARGUMENTS:
     type(ESMF_State)                      :: state
     character(*), intent(in)              :: name
     character(*), intent(in)              :: value
     integer,      intent(out), optional   :: rc
DESCRIPTION:

Set the attribute name inside of state using the convention NUOPC and purpose Instance.

The arguments are:

state
The ESMF_State object on which to set the attribute.
name
The name of the set attribute.
value
The value of the set attribute.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.29 NUOPC_SetTimestamp - Set the TimeStamp on a Field


INTERFACE:

   ! Private name; call using NUOPC_SetTimestamp()
   subroutine NUOPC_SetTimestampField(field, time, rc)
ARGUMENTS:
     type(ESMF_Field), intent(inout)         :: field
     type(ESMF_Time),  intent(in)            :: time
     integer,          intent(out), optional :: rc
DESCRIPTION:

Set the TimeStamp according to time on field.

This call should rarely be needed in user written code.

The arguments are:

field
The ESMF_Field object to be time stampped.
time
The ESMF_Time object defining the TimeStamp.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.30 NUOPC_SetTimestamp - Set the TimeStamp on Fields in a list


INTERFACE:

   ! Private name; call using NUOPC_SetTimestamp()
   subroutine NUOPC_SetTimestampFieldList(fieldList, time, selective, rc)
ARGUMENTS:
     type(ESMF_Field), intent(inout)         :: fieldList(:)
     type(ESMF_Time),  intent(in)            :: time
     logical,          intent(in),  optional :: selective
     integer,          intent(out), optional :: rc
DESCRIPTION:

Set the TimeStamp according to time on field.

This call should rarely be needed in user written code.

The arguments are:

fieldList
The list of ESMF_Field objects to be time stampped.
time
The ESMF_Time object defining the TimeStamp.
[selective]
If .true., then only set the TimeStamp on those fields for which the "Updated" attribute is equal to "true". Otherwise set the TimeStamp on all the fields. Default is .false..
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.31 NUOPC_SetTimestamp - Set the TimeStamp on Fields in a list from Clock


INTERFACE:

   ! Private name; call using NUOPC_SetTimestamp()
   subroutine NUOPC_SetTimestampFieldListClk(fieldList, clock, selective, rc)
ARGUMENTS:
     type(ESMF_Field), intent(inout)         :: fieldList(:)
     type(ESMF_Clock), intent(in)            :: clock
     logical,          intent(in),  optional :: selective
     integer,          intent(out), optional :: rc
DESCRIPTION:

Set the TimeStamp according to time on field.

This call should rarely be needed in user written code.

The arguments are:

fieldList
The list of ESMF_Field objects to be time stampped.
clock
The ESMF_Clock object defining the TimeStamp by its current time.
[selective]
If .true., then only set the TimeStamp on those fields for which the "Updated" attribute is equal to "true". Otherwise set the TimeStamp on all the fields. Default is .false..
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.32 NUOPC_SetTimestamp - Set the TimeStamp on all the Fields in a State


INTERFACE:

   ! Private name; call using NUOPC_SetTimestamp()
   subroutine NUOPC_SetTimestampState(state, time, selective, rc)
ARGUMENTS:
     type(ESMF_State), intent(inout)         :: state
     type(ESMF_Time),  intent(in)            :: time
     logical,          intent(in),  optional :: selective
     integer,          intent(out), optional :: rc
DESCRIPTION:

Set the TimeStamp according to clock on all the fields in state. Depending on selective, all or only some fields may be updated.

This call should rarely be needed in user written code. It is used by the generic Connector.

The arguments are:

state
The ESMF_State object holding the fields to be time stampped.
time
The ESMF_Time object defining the TimeStamp.
[selective]
If .true., then only set the TimeStamp on those fields for which the "Updated" attribute is equal to "true". Otherwise set the TimeStamp on all the fields. Default is .false..
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.9.33 NUOPC_SetTimestamp - Set the TimeStamp on all the Fields in a State from Clock


INTERFACE:

   ! Private name; call using NUOPC_SetTimestamp()
   subroutine NUOPC_SetTimestampStateClk(state, clock, selective, rc)
ARGUMENTS:
     type(ESMF_State), intent(inout)         :: state
     type(ESMF_Clock), intent(in)            :: clock
     logical,          intent(in),  optional :: selective
     integer,          intent(out), optional :: rc
DESCRIPTION:

Set the TimeStamp according to clock on all the fields in state. Depending on selective, all or only some fields may be updated.

This call should rarely be needed in user written code. It is used by the generic Connector.

The arguments are:

state
The ESMF_State object holding the fields to be time stampped.
clock
The ESMF_Clock object defining the TimeStamp by its current time.
[selective]
If .true., then only set the TimeStamp on those fields for which the "Updated" attribute is equal to "true". Otherwise set the TimeStamp on all the fields. Default is .false..
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.10 Auxiliary Routines

Auxiliary routines are provided with the NUOPC Layer as a convenience to the user. Typically more work is needed on these methods before considering them NUOPC core functionality.

3.10.1 NUOPC_Write - Write a distributed interpolation matrix to file in SCRIP format


INTERFACE:

   ! Private name; call using NUOPC_Write()
   subroutine NUOPC_SCRIPWrite(factorList, factorIndexList, fileName, &
     relaxedflag, rc)
ARGUMENTS:
     real(ESMF_KIND_R8), intent(in), target    :: factorList(:)
     integer,            intent(in), target    :: factorIndexList(:,:) 
     character(*),       intent(in)            :: fileName
     logical,            intent(in),  optional :: relaxedflag
     integer,            intent(out), optional :: rc
DESCRIPTION:

Write the destributed interpolaton matrix provided by factorList and factorIndexList to a SCRIP formatted NetCDF file. Each PET calls with its local list of factors and indices. The call then writes the distributed factors into a single file. If the file already exists, the contents is replaced by this call.

The arguments are:

factorList
The distributed factor list.
factorIndexList
The distributed list of source and destination indices.
fileName
The name of the file to be written to.
[relaxedflag]
If .true., then no error is returned even if the call cannot write the file due to library limitations. Default is .false..
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.10.2 NUOPC_Write - Write a distributed factorList to file


INTERFACE:

   ! Private name; call using NUOPC_Write()
   subroutine NUOPC_FactorsWrite(factorList, fileName, rc)
ARGUMENTS:
     real(ESMF_KIND_R8), pointer               :: factorList(:)
     character(*),       intent(in)            :: fileName
     integer,            intent(out), optional :: rc
DESCRIPTION:

THIS METHOD IS DEPRECATED. Use 3.10.1 instead.

Write the destributed factorList to file. Each PET calls with its local list of factors. The call then writes the distributed factors into a single file. The order of the factors in the file is first by PET, and within each PET the PET-local order is preserved. Changing the number of PETs for the same regrid operation will likely change the order of factors across PETs, and therefore files written will differ.

The arguments are:

factorList
The distributed factor list.
fileName
The name of the file to be written to.
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.10.3 NUOPC_Write - Write Field data to file


INTERFACE:

   ! Private name; call using NUOPC_Write()
   subroutine NUOPC_FieldWrite(field, fileName, overwrite, status, timeslice, &
     iofmt, relaxedflag, rc)
ARGUMENTS:
     type(ESMF_Field),           intent(in)            :: field
     character(*),               intent(in)            :: fileName
     logical,                    intent(in),  optional :: overwrite
     type(ESMF_FileStatus_Flag), intent(in),  optional :: status
     integer,                    intent(in),  optional :: timeslice
     type(ESMF_IOFmt_Flag),      intent(in),  optional :: iofmt
     logical,                    intent(in),  optional :: relaxedflag
     integer,                    intent(out), optional :: rc
DESCRIPTION:

Write the data in field to file under the field's "StandardName" attribute if supported by the iofmt.

The arguments are:

field
The ESMF_Field object whose data is to be written.
fileName
The name of the file to write to.
[overwrite]
A logical flag, the default is .false., i.e., existing Field data may not be overwritten. If .true., the data corresponding to each field's name will be be overwritten. If the timeslice option is given, only data for the given timeslice may be overwritten. Note that it is always an error to attempt to overwrite a NetCDF variable with data which has a different shape.
[status]
The file status. Valid options are ESMF_FILESTATUS_NEW, ESMF_FILESTATUS_OLD, ESMF_FILESTATUS_REPLACE, and ESMF_FILESTATUS_UNKNOWN (default).
[timeslice]
Time slice counter. Must be positive. The behavior of this option may depend on the setting of the overwrite flag:
overwrite = .false.:
If the timeslice value is less than the maximum time already in the file, the write will fail.
overwrite = .true.:
Any positive timeslice value is valid.
By default, i.e. by omitting the timeslice argument, no provisions for time slicing are made in the output file, however, if the file already contains a time axis for the variable, a timeslice one greater than the maximum will be written.
[iofmt]
The I/O format. Supported options are ESMF_IOFMT_NETCDF, ESMF_IOFMT_NETCDF4P, and ESMF_IOFMT_NETCDF4C. If not present, defaults to ESMF_IOFMT_NETCDF.
[relaxedflag]
If .true., then no error is returned even if the call cannot write the file due to library limitations, or because field does not contain any data. Default is .false..
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.10.4 NUOPC_Write - Write the Fields within a State to NetCDF files


INTERFACE:

   ! Private name; call using NUOPC_Write()
   subroutine NUOPC_StateWrite(state, fieldNameList, fileNamePrefix, overwrite, &
     status, timeslice, iofmt, relaxedflag, rc)
ARGUMENTS:
     type(ESMF_State),           intent(in)            :: state
     character(len=*),           intent(in),  optional :: fieldNameList(:)
     character(len=*),           intent(in),  optional :: fileNamePrefix
     logical,                    intent(in),  optional :: overwrite
     type(ESMF_FileStatus_Flag), intent(in),  optional :: status
     integer,                    intent(in),  optional :: timeslice
     type(ESMF_IOFmt_Flag),      intent(in),  optional :: iofmt
     logical,                    intent(in),  optional :: relaxedflag
     integer,                    intent(out), optional :: rc
DESCRIPTION:

Write the data of the fields within a state to NetCDF files. Each field is written to an individual file using the "StandardName" attribute as NetCDF attribute.

The arguments are:

state
The ESMF_State object containing the fields.
[fieldNameList]
List of names of the fields to be written. By default write all the fields in state.
[fileNamePrefix]
File name prefix, common to all the files written.
[overwrite]
A logical flag, the default is .false., i.e., existing Field data may not be overwritten. If .true., the data corresponding to each field's name will be be overwritten. If the timeslice option is given, only data for the given timeslice may be overwritten. Note that it is always an error to attempt to overwrite a NetCDF variable with data which has a different shape.
[status]
The file status. Valid options are ESMF_FILESTATUS_NEW, ESMF_FILESTATUS_OLD, ESMF_FILESTATUS_REPLACE, and ESMF_FILESTATUS_UNKNOWN (default).
[timeslice]
Time slice counter. Must be positive. The behavior of this option may depend on the setting of the overwrite flag:
overwrite = .false.:
If the timeslice value is less than the maximum time already in the file, the write will fail.
overwrite = .true.:
Any positive timeslice value is valid.
By default, i.e. by omitting the timeslice argument, no provisions for time slicing are made in the output file, however, if the file already contains a time axis for the variable, a timeslice one greater than the maximum will be written.
[iofmt]
The I/O format. Supported options are ESMF_IOFMT_NETCDF, ESMF_IOFMT_NETCDF4P, and ESMF_IOFMT_NETCDF4C. If not present, defaults to ESMF_IOFMT_NETCDF.
[relaxedflag]
If .true., then no error is returned even if the call cannot write the file due to library limitations. Default is .false..
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

3.10.5 NUOPC_Write - Write the Fields within a FieldBundle to NetCDF files


INTERFACE:

   ! Private name; call using NUOPC_Write()
   subroutine NUOPC_FieldBundleWrite(fieldbundle, fieldNameList, fileNamePrefix, overwrite, &
     status, timeslice, iofmt, relaxedflag, rc)
ARGUMENTS:
     type(ESMF_FieldBundle),     intent(in)            :: fieldbundle
     character(len=*),           intent(in),  optional :: fieldNameList(:)
     character(len=*),           intent(in),  optional :: fileNamePrefix
     logical,                    intent(in),  optional :: overwrite
     type(ESMF_FileStatus_Flag), intent(in),  optional :: status
     integer,                    intent(in),  optional :: timeslice
     type(ESMF_IOFmt_Flag),      intent(in),  optional :: iofmt
     logical,                    intent(in),  optional :: relaxedflag
     integer,                    intent(out), optional :: rc
DESCRIPTION:

Write the data of the fields within a fieldbundle to NetCDF files. Each field is written to an individual file using the "StandardName" attribute as NetCDF attribute.

The arguments are:

fieldbundle
The ESMF_FieldBundle object containing the fields.
[fieldNameList]
List of names of the fields to be written. By default write all the fields in fieldbundle.
[fileNamePrefix]
File name prefix, common to all the files written.
[overwrite]
A logical flag, the default is .false., i.e., existing Field data may not be overwritten. If .true., the data corresponding to each field's name will be be overwritten. If the timeslice option is given, only data for the given timeslice may be overwritten. Note that it is always an error to attempt to overwrite a NetCDF variable with data which has a different shape.
[status]
The file status. Valid options are ESMF_FILESTATUS_NEW, ESMF_FILESTATUS_OLD, ESMF_FILESTATUS_REPLACE, and ESMF_FILESTATUS_UNKNOWN (default).
[timeslice]
Time slice counter. Must be positive. The behavior of this option may depend on the setting of the overwrite flag:
overwrite = .false.:
If the timeslice value is less than the maximum time already in the file, the write will fail.
overwrite = .true.:
Any positive timeslice value is valid.
By default, i.e. by omitting the timeslice argument, no provisions for time slicing are made in the output file, however, if the file already contains a time axis for the variable, a timeslice one greater than the maximum will be written.
[iofmt]
The I/O format. Supported options are ESMF_IOFMT_NETCDF, ESMF_IOFMT_NETCDF4P, and ESMF_IOFMT_NETCDF4C. If not present, defaults to ESMF_IOFMT_NETCDF.
[relaxedflag]
If .true., then no error is returned even if the call cannot write the file due to library limitations. Default is .false..
[rc]
Return code; equals ESMF_SUCCESS if there are no errors.

esmf_support@ucar.edu