ESMF superstructure classes define an architecture for assembling Earth system applications from modeling components. A component may be defined in terms of the physical domain that it represents, such as an atmosphere or sea ice model. It may also be defined in terms of a computational function, such as a data assimilation system. Earth system research often requires that such components be coupled together to create an application. By coupling we mean the data transformations and, on parallel computing systems, data transfers, that are necessary to allow data from one component to be utilized by another. ESMF offers regridding methods and other tools to simplify the organization and execution of inter-component data exchanges.
In addition to components defined at the level of major physical domains and computational functions, components may be defined that represent smaller computational functions within larger components, such as the transformation of data between the physics and dynamics in a spectral atmosphere model, or the creation of nested higher resolution regions within a coarser grid. The objective is to couple components at varying scales both flexibly and efficiently. ESMF encourages a hierachical application structure, in which large components branch into smaller sub-components (see Figure 2). ESMF also makes it easier for the same component to be used in multiple contexts without changes to its source code.
Key Features |
Modular, component-based architecture. |
Hierarchical assembly of components into applications. |
Use of components in multiple contexts without modification. |
Sequential or concurrent component execution. |
Single program, multiple datastream (SPMD) applications for maximum portability and reconfigurability. |
Multiple program, multiple datastream (MPMD) option for flexibility. |
There are a small number of classes in the ESMF superstructure:
The second part of an ESMF Component is user code, such as a model or data assimilation system. Users set entry points within their code so that it is callable by the framework. In practice, setting entry points means that within user code there are calls to ESMF methods that associate the name of a Fortran subroutine with a corresponding standard ESMF operation. For example, a user-written initialization routine called myOceanInit might be associated with the standard initialize routine of an ESMF Gridded Component named ``myOcean'' that represents an ocean model.
An ESMF coupled application typically involves an AppDriver, a parent Gridded Component, two or more child Gridded Components that require an inter-component data exchange, and one or more Coupler Components.
The parent Gridded Component is responsible for creating the child Gridded Components that are exchanging data, for creating the Coupler, for creating the necessary Import and Export States, and for setting up the desired sequencing. The AppDriver ``main'' routine calls the parent Gridded Component's initialize, run, and finalize methods in order to execute the application. For each of these standard methods, the parent Gridded Component in turn calls the corresponding methods in the child Gridded Components and the Coupler Component. For example, consider a simple coupled ocean/atmosphere simulation. When the initialize method of the parent Gridded Component is called by the AppDriver, it in turn calls the initialize methods of its child atmosphere and ocean Gridded Components, and the initialize method of an ocean-to-atmosphere Coupler Component. Figure 3 shows this schematically.
Components are allocated computational resources in the form of Persistent Execution Threads, or PETs. A list of a Component's PETs is contained in a structure called a Virtual Machine, or VM. The VM also contains information about the topology and characteristics of the underlying computer. Components are created hierarchically, with parent Components creating child Components and allocating some or all of their PETs to each one. By default ESMF creates a new VM for each child Component, which allows Components to tailor their VM resources to match their needs. In some cases a child may want to share its parent's VM - ESMF supports this too.
A Gridded Component may exist across all the PETs in an application. A Gridded Component may also reside on a subset of PETs in an application. These PETs may wholly coincide with, be wholly contained within, or wholly contain another Component.
When a set of Gridded Components and a Coupler runs in sequence on the same set of PETs the application is executing in a sequential mode. When Gridded Components are created and run on mutually exclusive sets of PETs, and are coupled by a Coupler Component that extends over the union of these sets, the mode of execution is concurrent.
Figure 4 illustrates a typical configuration for a simple coupled sequential application, and Figure 5 shows a possible configuration for the same application running in a concurrent mode.
Parent Components can select if and when to wait for concurrently executing child Components, synchronizing only when required.
It is possible for ESMF applications to contain some Component sets that are executing sequentially and others that are executing concurrently. We might have, for example, atmosphere and land Components created on the same subset of PETs, ocean and sea ice Components created on the remainder of PETs, and a Coupler created across all the PETs in the application.
All data transfers within an ESMF application occur within a component. For example, a Gridded Component may contain halo updates. Another example is that a Coupler Component may redistribute data between two Gridded Components. As a result, the architecture of ESMF does not depend on any particular data communication mechanism, and new communication schemes can be introduced without affecting the overall structure of the application.
Since all data communication happens within a component, a Coupler Component must be created on the union of the PETs of all the Gridded Components that it couples.
The scope of distributed objects is the VM of the currently executing Component. For this reason, all PETs in the current VM must make the same distributed object creation calls. When a Coupler Component running on a superset of a Gridded Component's PETs needs to make communication calls involving objects created by the Gridded Component, an ESMF-supplied function called ESMF_StateReconcile() creates proxy objects for those PETs that had no previous information about the distributed objects. Proxy objects contain no local data but can be used in communication calls (such as regrid or redistribute) to describe the remote source for data being moved to the current PET, or to describe the remote destination for data being moved from the local PET. Figure 6 is a simple schematic that shows the sequence of events in a reconcile call.
The ESMF design enables the user to configure ESMF applications so that data is transferred directly from one component to another, without requiring that it be copied or sent to a different data buffer as an interim step. This is likely to be the most efficient way of performing inter-component coupling. However, if desired, an application can also be configured so that data from a source component is sent to a distinct set of Coupler Component PETs for processing before being sent to its destination.
The ability to overlap computation with communication is essential for performance. When running with ESMF the user can initiate data sends during Gridded Component execution, as soon as the data is ready. Computations can then proceed simultaneously with the data transfer.
The following is a simplified UML diagram showing the relationships among ESMF superstructure classes. See Appendix A, A Brief Introduction to UML, for a translation table that lists the symbols in the diagram and their meaning.
The ESMF Application Driver (ESMF_AppDriver), is a generic ESMF driver program that contains a ``main.'' Simpler applications may be able to use an Application Driver without modification; for more complex applications, an Application Driver can be used as an extendable template.
ESMF provides a number of different Application Drivers in the $ESMF_DIR/src/Superstructure/AppDriver directory. An appropriate one can be chosen depending on how the application is to be structured. Options when deciding how to structure an application include choices about:
In a sequential execution model every Component executes on all PETs, with each Component completing execution before the next Component begins. This has the appeal of simplicity of data consumption and production: when a Gridded Component starts all required data is available for use, and when a Gridded Component finishes all data produced is ready for consumption by the next Gridded Component. This approach also has the possibility of less data movement if the grid and data decomposition is done such that each processor's memory contains the data needed by the next Component.
In a concurrent execution model subgroups of PETs run Gridded Components and multiple Gridded Components are active at the same time. Data exchange must be coordinated between Gridded Components so that data deadlock does not occur. This strategy has the advantage of allowing coupling to other Gridded Components at any time during the computational process, including not having to return to the calling level of code before making data available.
Coupler Components are responsible for taking data from one Gridded Component and putting it into the form expected by another Gridded Component. This might include regridding, change of units, averaging, or binning.
Coupler Components can be written for pairwise data exchange: the Coupler Component takes data from a single Component and transforms it for use by another single Gridded Component. This simplifies the structure of the Coupler Component code.
Couplers can also be written using a hub and spoke model where a single Coupler accepts data from all other Components, can do data merging or splitting, and formats data for all other Components.
Multiple Couplers, using either of the above two models or some mixture of these approaches, are also possible.
The ESMF framework currently has Fortran interfaces for all public functions. Some functions also have C interfaces, and the number of these is expected to increase over time.
The simplest way to run an application is to run the same executable program on all PETs. Different Components can still be run on mutually exclusive PETs by using branching (e.g., if this is PET 1, 2, or 3, run Component A, if it is PET 4, 5, or 6 run Component B). This is a SPMD model, Single Program Multiple Data.
The alternative is to start a different executable program on different PETs. This is a MPMD model, Multiple Program Multiple Data. There are complications with many job control systems on multiprocessor machines in getting the different executables started, and getting inter-process communcations established. ESMF currently has some support for MPMD: different Components can run as separate executables, but the Coupler that transfers data between the Components must still run on the union of their PETs. This means that the Coupler Component must be linked into all of the executables.
DESCRIPTION:
The ESMF_TerminationFlag determines how an ESMF application is shut down.
Valid values are:
ESMF encourages application organization in which there is a single top-level Gridded Component. This provides a simple, clear sequence of operations at the highest level, and also enables the entire application to be treated as a sub-Component of another, larger application if desired. When a simple application is organized in this fashion the standard AppDriver can probably be used without much modification.
Examples of program organization using the AppDriver can be found in the src/Superstructure/AppDriver directory. A set of subdirectories within the AppDriver directory follows the naming convention:
<seq|concur>_<pairwise|hub>_<f|c>driver_<spmd|mpmd>
The example that is currently implemented is seq_pairwise_fdriver_spmd, which has sequential component execution, a pairwise coupler, a main program in Fortran, and all processors launching the same executable. It is also copied automatically into a top-level quick_start directory at compilation time.
The user can copy the AppDriver files into their own local directory. Some of the files can be used unchanged. Others are template files which have the rough outline of the code but need additional application-specific code added in order to perform a meaningful function. The README file in the AppDriver subdirectory or quick_start directory contains instructions about which files to change.
Examples of concurrent component execution can be found in the system tests that are bundled with the ESMF distribution.
------------------------------------------------------------------------------ ------------------------------------------------------------------------------ EXAMPLE: This is an AppDriver.F90 file for a sequential ESMF application. ------------------------------------------------------------------------------ ------------------------------------------------------------------------------ The ChangeMe.F90 file that's included below contains a number of definitions that are used by the AppDriver, such as the name of the application's main configuration file and the name of the application's SetServices routine. This file is in the same directory as the AppDriver.F90 file. ------------------------------------------------------------------------------ #include "ChangeMe.F90" program ESMF_AppDriver #define ESMF_METHOD "program ESMF_AppDriver" #include "ESMF.h" ! ESMF module, defines all ESMF data types and procedures use ESMF_Mod ! Gridded Component registration routines. Defined in "ChangeMe.F90" use USER_APP_Mod, only : SetServices => USER_APP_SetServices implicit none ------------------------------------------------------------------------------ Define local variables ------------------------------------------------------------------------------ ! Components and States type(ESMF_GridComp) :: compGridded type(ESMF_State) :: defaultstate ! Configuration information type(ESMF_Config) :: config ! A common Grid type(ESMF_Grid) :: grid ! A Clock, a Calendar, and timesteps type(ESMF_Clock) :: clock type(ESMF_TimeInterval) :: timeStep type(ESMF_Time) :: startTime type(ESMF_Time) :: stopTime ! Variables related to the Grid integer :: i_max, j_max ! Return codes for error checks integer :: rc, localrc ------------------------------------------------------------------------------ Initialize ESMF. Note that an output Log is created by default. ------------------------------------------------------------------------------ call ESMF_Initialize(defaultCalendar=ESMF_CAL_GREGORIAN, rc=localrc) if (ESMF_LogMsgFoundError(localrc, ESMF_ERR_PASSTHRU, & ESMF_CONTEXT, rcToReturn=rc)) & call ESMF_Finalize(rc=localrc, terminationflag=ESMF_ABORT) call ESMF_LogWrite("ESMF AppDriver start", ESMF_LOG_INFO) ------------------------------------------------------------------------------ Create and load a configuration file. The USER_CONFIG_FILE is set to sample.rc in the ChangeMe.F90 file. The sample.rc file is also included in the directory with the AppDriver.F90 file. ------------------------------------------------------------------------------ config = ESMF_ConfigCreate(rc=localrc) if (ESMF_LogMsgFoundError(localrc, ESMF_ERR_PASSTHRU, & ESMF_CONTEXT, rcToReturn=rc)) & call ESMF_Finalize(rc=localrc, terminationflag=ESMF_ABORT) call ESMF_ConfigLoadFile(config, USER_CONFIG_FILE, rc = localrc) if (ESMF_LogMsgFoundError(localrc, ESMF_ERR_PASSTHRU, & ESMF_CONTEXT, rcToReturn=rc)) & call ESMF_Finalize(rc=localrc, terminationflag=ESMF_ABORT) ------------------------------------------------------------------------------ Get configuration information. A configuration file like sample.rc might include: - size and coordinate information needed to create the default Grid. - the default start time, stop time, and running intervals for the main time loop. ------------------------------------------------------------------------------ call ESMF_ConfigGetAttribute(config, i_max, 'I Counts:', default=10, rc=localrc) if (ESMF_LogMsgFoundError(localrc, ESMF_ERR_PASSTHRU, & ESMF_CONTEXT, rcToReturn=rc)) & call ESMF_Finalize(rc=localrc, terminationflag=ESMF_ABORT) call ESMF_ConfigGetAttribute(config, j_max, 'J Counts:', default=40, rc=localrc) if (ESMF_LogMsgFoundError(localrc, ESMF_ERR_PASSTHRU, & ESMF_CONTEXT, rcToReturn=rc)) & call ESMF_Finalize(rc=localrc, terminationflag=ESMF_ABORT) ------------------------------------------------------------------------------ Create the top Gridded Component. ------------------------------------------------------------------------------ compGridded = ESMF_GridCompCreate(name="ESMF Gridded Component", rc=localrc) if (ESMF_LogMsgFoundError(localrc, ESMF_ERR_PASSTHRU, & ESMF_CONTEXT, rcToReturn=rc)) & call ESMF_Finalize(rc=localrc, terminationflag=ESMF_ABORT) call ESMF_LogWrite("Component Create finished", ESMF_LOG_INFO) ------------------------------------------------------------------------------ Register the set services method for the top Gridded Component. ------------------------------------------------------------------------------ call ESMF_GridCompSetServices(compGridded, SetServices, rc) if (ESMF_LogMsgFoundError(rc, "Registration failed", rc)) & call ESMF_Finalize(rc=localrc, terminationflag=ESMF_ABORT) ------------------------------------------------------------------------------ Create and initialize a Clock. ------------------------------------------------------------------------------ call ESMF_TimeIntervalSet(timeStep, s=2, rc=localrc) if (ESMF_LogMsgFoundError(localrc, ESMF_ERR_PASSTHRU, & ESMF_CONTEXT, rcToReturn=rc)) & call ESMF_Finalize(rc=localrc, terminationflag=ESMF_ABORT) call ESMF_TimeSet(startTime, yy=2004, mm=9, dd=25, rc=localrc) if (ESMF_LogMsgFoundError(localrc, ESMF_ERR_PASSTHRU, & ESMF_CONTEXT, rcToReturn=rc)) & call ESMF_Finalize(rc=localrc, terminationflag=ESMF_ABORT) call ESMF_TimeSet(stopTime, yy=2004, mm=9, dd=26, rc=localrc) if (ESMF_LogMsgFoundError(localrc, ESMF_ERR_PASSTHRU, & ESMF_CONTEXT, rcToReturn=rc)) & call ESMF_Finalize(rc=localrc, terminationflag=ESMF_ABORT) clock = ESMF_ClockCreate("Application Clock", timeStep, startTime, & stopTime, rc=localrc) if (ESMF_LogMsgFoundError(localrc, ESMF_ERR_PASSTHRU, & ESMF_CONTEXT, rcToReturn=rc)) & call ESMF_Finalize(rc=localrc, terminationflag=ESMF_ABORT) ------------------------------------------------------------------------------ Create and initialize a Grid. The default lower indices for the Grid are (/1,1/). The upper indices for the Grid are read in from the sample.rc file, where they are set to (/10,40/). This means a Grid will be created with 10 grid cells in the x direction and 40 grid cells in the y direction. The Grid section in the Reference Manual shows how to set coordinates. ------------------------------------------------------------------------------ grid = ESMF_GridCreateShapeTile(maxIndex=(/i_max, j_max/), & name="source grid", rc=localrc) if (ESMF_LogMsgFoundError(localrc, ESMF_ERR_PASSTHRU, & ESMF_CONTEXT, rcToReturn=rc)) & call ESMF_Finalize(rc=localrc, terminationflag=ESMF_ABORT) ! Attach the grid to the Component call ESMF_GridCompSet(compGridded, grid=grid, rc=localrc) if (ESMF_LogMsgFoundError(localrc, ESMF_ERR_PASSTHRU, & ESMF_CONTEXT, rcToReturn=rc)) & call ESMF_Finalize(rc=localrc, terminationflag=ESMF_ABORT) ------------------------------------------------------------------------------ Create and initialize a State to use for both import and export. In a real code, separate import and export States would normally be created. ------------------------------------------------------------------------------ defaultstate = ESMF_StateCreate("Default State", rc=localrc) if (ESMF_LogMsgFoundError(localrc, ESMF_ERR_PASSTHRU, & ESMF_CONTEXT, rcToReturn=rc)) & call ESMF_Finalize(rc=localrc, terminationflag=ESMF_ABORT) ------------------------------------------------------------------------------ Call the initialize, run, and finalize methods of the top component. When the initialize method of the top component is called, it will in turn call the initialize methods of all its child components, they will initialize their children, and so on. The same is true of the run and finalize methods. ------------------------------------------------------------------------------ call ESMF_GridCompInitialize(compGridded, defaultstate, defaultstate, & clock, rc=localrc) if (ESMF_LogMsgFoundError(rc, "Initialize failed", rc)) & call ESMF_Finalize(rc=localrc, terminationflag=ESMF_ABORT) call ESMF_GridCompRun(compGridded, defaultstate, defaultstate, & clock, rc=localrc) if (ESMF_LogMsgFoundError(rc, "Run failed", rc)) & call ESMF_Finalize(rc=localrc, terminationflag=ESMF_ABORT) call ESMF_GridCompFinalize(compGridded, defaultstate, defaultstate, & clock, rc=localrc) if (ESMF_LogMsgFoundError(rc, "Finalize failed", rc)) & call ESMF_Finalize(rc=localrc, terminationflag=ESMF_ABORT) ------------------------------------------------------------------------------ Destroy objects. ------------------------------------------------------------------------------ call ESMF_ClockDestroy(clock, rc=localrc) if (ESMF_LogMsgFoundError(localrc, ESMF_ERR_PASSTHRU, & ESMF_CONTEXT, rcToReturn=rc)) & call ESMF_Finalize(rc=localrc, terminationflag=ESMF_ABORT) call ESMF_StateDestroy(defaultstate, rc=localrc) if (ESMF_LogMsgFoundError(localrc, ESMF_ERR_PASSTHRU, & ESMF_CONTEXT, rcToReturn=rc)) & call ESMF_Finalize(rc=localrc, terminationflag=ESMF_ABORT) call ESMF_GridCompDestroy(compGridded, rc=localrc) if (ESMF_LogMsgFoundError(localrc, ESMF_ERR_PASSTHRU, & ESMF_CONTEXT, rcToReturn=rc)) & call ESMF_Finalize(rc=localrc, terminationflag=ESMF_ABORT) ------------------------------------------------------------------------------ Finalize and clean up. ------------------------------------------------------------------------------ call ESMF_Finalize() end program ESMF_AppDriver
There are a few methods that every ESMF application must contain. First, ESMF_Initialize() and ESMF_Finalize() are in complete analogy to MPI_Init() and MPI_Finalize() known from MPI. All ESMF programs, serial or parallel, must initialize the ESMF system at the beginning, and finalize it at the end of execution. The behavior of calling any ESMF method before ESMF_Initialize(), or after ESMF_Finalize() is undefined.
Second, every ESMF Component that is accessed by an ESMF application requires that its set services routine is called through ESMF_<Grid/Cpl>CompSetServices(). The Component must implement one public entry point, its set services routine, that can be called through the ESMF_<Grid/Cpl>CompSetServices() library routine. The Component set services routine is responsible for setting entry points for the standard ESMF Component methods Initialize, Run, and Finalize.
Finally, the Component library call ESMF_<Grid/Cpl>CompSetVM() can optionally be issues before calling ESMF_<Grid/Cpl>CompSetServices(). Similar to ESMF_<Grid/Cpl>CompSetServices(), the ESMF_<Grid/Cpl>CompSetVM() call requires a public entry point into the Component. It allows the Component to adjust certain aspects of its execution environment, i.e. its own VM, before it is started up.
The following sections discuss the above mentioned aspects in more detail.
INTERFACE:
subroutine ESMF_Initialize(defaultConfigFileName, defaultCalendar, & defaultLogFileName, defaultLogType, mpiCommunicator, & IOUnitLower, IOUnitUpper, vm, rc)ARGUMENTS:
character(len=*), intent(in), optional :: defaultConfigFileName type(ESMF_CalendarType), intent(in), optional :: defaultCalendar character(len=*), intent(in), optional :: defaultLogFileName type(ESMF_LogType), intent(in), optional :: defaultLogType integer, intent(in), optional :: mpiCommunicator integer, intent(in), optional :: IOUnitLower integer, intent(in), optional :: IOUnitUpper type(ESMF_VM), intent(out), optional :: vm integer, intent(out), optional :: rcDESCRIPTION:
This method must be called once on each PET before any other ESMF methods are used. The method contains a barrier before returning, ensuring that all processes made it successfully through initialization.
Typically ESMF_Initialize() will call MPI_Init() internally unless MPI has been initialized by the user code before initializing the framework. If the MPI initialization is left to ESMF_Initialize() it inherits all of the MPI implementation dependent limitations of what may or may not be done before MPI_Init(). For instance, it is unsafe for some MPI implementations, such as MPICH, to do IO before the MPI environment is initialized. Please consult the documentation of your MPI implementation for details.
Note that when using MPICH as the MPI library, ESMF needs to use the application command line arguments for MPI_Init(). However, ESMF acquires these arguments internally and the user does not need to worry about providing them. Also, note that ESMF does not alter the command line arguments, so that if the user obtains them they will be as specified on the command line (including those which MPICH would normally strip out).
By default, ESMF_Initialize() will open multiple error log files, one per processor. This is very useful for debugging purpose. However, when running the application on a large number of processors, opening a large number of log files and writing log messages from all the processors could become a performance bottleneck. Therefore, it is recommended to turn the Error Log feature off in these situations by setting defaultLogType to ESMF_LOG_NONE.
When integrating ESMF with applications where Fortran unit number conflicts exist, the optional IOUnitLower and IOUnitUpper arguments may be used to specify an alternate unit number range. See section 46.2.1 for more information on how ESMF uses Fortran unit numbers.
Before exiting the application the user must call ESMF_Finalize() to release resources and clean up ESMF gracefully.
The arguments are:
INTERFACE:
subroutine ESMF_Finalize(terminationflag, rc)ARGUMENTS:
type(ESMF_TerminationFlag), intent(in), optional :: terminationflag integer, intent(out), optional :: rcDESCRIPTION:
This must be called once on each PET before the application exits to allow ESMF to flush buffers, close open connections, and release internal resources cleanly. The optional argument terminationflag may be used to indicate the mode of termination. Note that this call must be issued only once per PET with terminationflag=ESMF_FINAL, and that this call may not be followed by ESMF_Initialize(). This last restriction means that it is not possible to restart ESMF within the same execution.
The arguments are:
Many programs call some library routines. The library documentation must explain what the routine name is, what arguments are required and what are optional, and what the code does.
In contrast, all ESMF components must be written to be called by another part of the program; in effect, an ESMF component takes the place of a library. The interface is prescribed by the framework, and the component writer must provide specific subroutines which have standard argument lists and perform specific operations. For technical reasons none of the arguments in user-provided subroutines must be declared as optional.
The only required public interface of a Component is its set services method. This subroutine must have an externally accessible name (be a public symbol), take a component as the first argument, and an integer return code as the second. Both arguments are required and must not be declared as optional. If an intent is specified in the interface it must be intent(inout) for the first and intent(out) for the second argument. The subroutine name is not predefined, it is set by the component writer, but must be provided as part of the component documentation.
The required function that the set services subroutine must provide is to specify the user-code entry points for the standard ESMF Component methods. To this end the user-written set services routine calls the ESMF_<Grid/Cpl>CompSetEntryPoint() method to set each Component entry point.
See sections 15.3.1 and 16.2.1 for examples of how to write a user-code set services routine.
Note that a component does not call its own set services routine; the AppDriver or parent component code, which is creating a component, will first call ESMF_<Grid/Cpl>CompCreate() to create an "empty" component, and then must call into ESMF_<Grid/Cpl>CompSetServices(), supplying the user-code set services routine as an argument. The framework calls into the user-code set services, after the Component's VM has been started up.
The required standard ESMF Component methods, for which user-code entry points must be set, are Initialize, Run, and Finalize. Currently optional, a Component may also set entry points for the WriteRestart and ReadRestart methods.
Sections 15.3.1 and 16.2.1 provide examples of how the entry points for Initialize, Run, and Finalize are set during the user-code set services routine, using the ESMF_<Grid/Cpl>CompSetEntryPoint() library call.
All standard user-code methods must abide exactly to the prescribed interfaces. None of the arguments must be declared as optional.
The names of the Initialize, Run, and Finalize user-code subroutines do not need to be public; in fact it is far better for them to be private to lower the chances of public symbol clashes between different components.
See sections 15.3.2, 15.3.3, 15.3.4, and 16.2.2, 16.2.3, 16.2.4 for examples of how to write entry points for the standard ESMF Component methods.
When the AppDriver or parent component code calls ESMF_<Grid/Cpl>CompCreate() it has the option to specify a petList argument. All of the parent PETs contained in this list become resources of the child component. By default all of the parent PETs are provided to the child component.
Unless the optional contextflag argument is used during ESMF_<Grid/Cpl>CompCreate(), to indicate that the child component will execute in the same VM as the parent, the child component has the option to set certain aspects of how the provided resources are to be used when executing child component methods. The resources provided via the parent PETs are the associated processing elements (PEs) and virtual address spaces (VASs).
The optional user-written set vm routine is called from the parent through the ESMF_<Grid/Cpl>CompSetVM() library code, and is the only place where the child component can set aspects of its own VM before it is started up. The child component's VM must be running before its set services routine can be called, and thus the optional ESMF_<Grid/Cpl>CompSetVM() call must be placed before ESMF_<Grid/Cpl>CompSetServices().
If called by the parent, the user-code set vm routine has the option to specify how the PETs of the child component share the provided parent PEs. Further, PETs on the same single system image can be set to run multi-threaded, within a reduced number of VAS, allowing a component to leverage shared memory concepts.
Sections 15.3.5 and 16.2.5 provide examples for simple user-written set vm routines.
DESCRIPTION:
The ESMF_GridCompType flag identifies what sort of physical domain
or computational function a particular ESMF_GridComp represents.
The flag values are purely informational; they are not used anywhere
within the framework. Use of this flag is optional.
Valid values are:
A Gridded Component is a computational entity which consumes and produces data. It uses a State object to exchange data between itself and other Components. It uses a Clock object to manage time, and a VM to describe its own and its child components' computational resources.
This section shows how to create Gridded Components. For demonstrations of the use of Gridded Components, see the system tests that are bundled with the ESMF software distribution. These can be found in the directory esmf/src/system_tests.
Every ESMF_GridComp is required to provide and document a public set services routine. It can have any name, but must follow the declaration below: a subroutine which takes an ESMF_GridComp as the first argument, and an integer return code as the second. Both arguments are required and must not be declared as optional. If an intent is specified in the interface it must be intent(inout) for the first and intent(out) for the second argument.
The set services routine must call the ESMF method ESMF_GridCompSetEntryPoint() to register with the framework what user-code subroutines should be called to initialize, run, and finalize the component. There are additional routines which can be registered as well, for checkpoint and restart functions.
Note that the actual subroutines being registered do not have to be public to this module; only the set services routine itself must be available to be used by other code.
! Example Gridded Component module ESMF_GriddedCompEx ! ESMF Framework module use ESMF_Mod implicit none public GComp_SetServices public GComp_SetVM contains subroutine GComp_SetServices(comp, rc) type(ESMF_GridComp) :: comp ! must not be optional integer, intent(out) :: rc ! must not be optional ! Set the entry points for standard ESMF Component methods call ESMF_GridCompSetEntryPoint(comp, ESMF_SETINIT, userRoutine=GComp_Init, rc=rc) call ESMF_GridCompSetEntryPoint(comp, ESMF_SETRUN, userRoutine=GComp_Run, rc=rc) call ESMF_GridCompSetEntryPoint(comp, ESMF_SETFINAL, userRoutine=GComp_Final, rc=rc) rc = ESMF_SUCCESS end subroutine
When a higher level component is ready to begin using an ESMF_GridComp, it will call its initialize routine.
The component writer must supply a subroutine with the exact interface shown below. Arguments must not be declared as optional, and the types and order must match.
At initialization time the component can allocate data space, open data files, set up initial conditions; anything it needs to do to prepare to run.
The rc return code should be set if an error occurs, otherwise the value ESMF_SUCCESS should be returned.
subroutine GComp_Init(comp, importState, exportState, clock, rc) type(ESMF_GridComp) :: comp ! must not be optional type(ESMF_State) :: importState, exportState ! must not be optional type(ESMF_Clock) :: clock ! must not be optional integer, intent(out) :: rc ! must not be optional print *, "Gridded Comp Init starting" ! This is where the model specific setup code goes. ! If the initial Export state needs to be filled, do it here. !call ESMF_StateAdd(exportState, field, rc) !call ESMF_StateAdd(exportState, bundle, rc) print *, "Gridded Comp Init returning" rc = ESMF_SUCCESS end subroutine GComp_Init
During the execution loop, the run routine may be called many times. Each time it should read data from the importState, use the clock to determine what the current time is in the calling component, compute new values or process the data, and produce any output and place it in the exportState.
When a higher level component is ready to use the ESMF_GridComp it will call its run routine.
The component writer must supply a subroutine with the exact interface shown below. Arguments must not be declared as optional, and the types and order must match.
It is expected that this is where the bulk of the model computation or data analysis will occur.
The rc return code should be set if an error occurs, otherwise the value ESMF_SUCCESS should be returned.
subroutine GComp_Run(comp, importState, exportState, clock, rc) type(ESMF_GridComp) :: comp ! must not be optional type(ESMF_State) :: importState, exportState ! must not be optional type(ESMF_Clock) :: clock ! must not be optional integer, intent(out) :: rc ! must not be optional print *, "Gridded Comp Run starting" ! call ESMF_StateGet(), etc to get fields, bundles, arrays ! from import state. ! This is where the model specific computation goes. ! Fill export state here using ESMF_StateAdd(), etc print *, "Gridded Comp Run returning" rc = ESMF_SUCCESS end subroutine GComp_Run
At the end of application execution, each ESMF_GridComp should deallocate data space, close open files, and flush final results. These functions should be placed in a finalize routine.
The component writer must supply a subroutine with the exact interface shown below. Arguments must not be declared as optional, and the types and order must match.
The rc return code should be set if an error occurs, otherwise the value ESMF_SUCCESS should be returned.
subroutine GComp_Final(comp, importState, exportState, clock, rc) type(ESMF_GridComp) :: comp ! must not be optional type(ESMF_State) :: importState, exportState ! must not be optional type(ESMF_Clock) :: clock ! must not be optional integer, intent(out) :: rc ! must not be optional print *, "Gridded Comp Final starting" ! Add whatever code here needed print *, "Gridded Comp Final returning" rc = ESMF_SUCCESS end subroutine GComp_Final
Every ESMF_GridComp can optionally provide and document a public set vm routine. It can have any name, but must follow the declaration below: a subroutine which takes an ESMF_GridComp as the first argument, and an integer return code as the second. Both arguments are required and must not be declared as optional. If an intent is specified in the interface it must be intent(inout) for the first and intent(out) for the second argument.
The set vm routine is the only place where the child component can use the ESMF_GridCompSetVMMaxPEs(), or ESMF_GridCompSetVMMaxThreads(), or ESMF_GridCompSetVMMinThreads() call to modify aspects of its own VM.
A component's VM is started up right before its set services routine is entered. ESMF_GridCompSetVM() is executing in the parent VM, and must be called before ESMF_GridCompSetServices().
subroutine GComp_SetVM(comp, rc) type(ESMF_GridComp) :: comp ! must not be optional integer, intent(out) :: rc ! must not be optional type(ESMF_VM) :: vm logical :: pthreadsEnabled ! Test for Pthread support, all SetVM calls require it call ESMF_VMGetGlobal(vm, rc=rc) call ESMF_VMGet(vm, pthreadsEnabledFlag=pthreadsEnabled, rc=rc) if (pthreadsEnabled) then ! run PETs single-threaded call ESMF_GridCompSetVMMinThreads(comp, rc=rc) endif rc = ESMF_SUCCESS end subroutine end module ESMF_GriddedCompEx
ESMF provides the concept of an Internal State that is associated with a Component. Through the Internal State API a user can attach a private data block to a Component, and later retrieve a pointer to this memory allocation. Setting and getting of Internal State information are supported from anywhere in the Component's SetServices, Initialize, Run, or Finalize code.
The code below demonstrates the basic Internal State API of ESMF_<Grid|Cpl>SetInternalState() and ESMF_<Grid|Cpl>GetInternalState(). Notice that an extra level of indirection to the user data is necessary!
! ESMF Framework module use ESMF_Mod implicit none type(ESMF_GridComp) :: comp integer :: rc, finalrc ! Internal State Variables type testData sequence integer :: testValue real :: testScaling end type type dataWrapper sequence type(testData), pointer :: p end type type(dataWrapper) :: wrap1, wrap2 type(testData), target :: data type(testData), pointer :: datap ! extra level of indirection finalrc = ESMF_SUCCESS !------------------------------------------------------------------------- call ESMF_Initialize(defaultlogfilename="InternalStateEx.Log", & defaultlogtype=ESMF_LOG_MULTI, rc=rc) if (rc .ne. ESMF_SUCCESS) finalrc = ESMF_FAILURE !------------------------------------------------------------------------- ! Creation of a Component comp = ESMF_GridCompCreate(name="test", rc=rc) if (rc .ne. ESMF_SUCCESS) finalrc = ESMF_FAILURE !------------------------------------------------------------------------- ! This could be called, for example, during a Component's initialize phase. ! Initialize private data block data%testValue = 4567 data%testScaling = 0.5 ! Set Internal State wrap1%p => data call ESMF_GridCompSetInternalState(comp, wrap1, rc) if (rc .ne. ESMF_SUCCESS) finalrc = ESMF_FAILURE !------------------------------------------------------------------------- ! This could be called, for example, during a Component's run phase. ! Get Internal State call ESMF_GridCompGetInternalState(comp, wrap2, rc) if (rc .ne. ESMF_SUCCESS) finalrc = ESMF_FAILURE ! Access private data block and verify data datap => wrap2%p if ((datap%testValue .ne. 4567) .or. (datap%testScaling .ne. 0.5)) then print *, "did not get same values back" finalrc = ESMF_FAILURE else print *, "got same values back from GetInternalState as original" endif
When working with ESMF Internal States it is important to consider the applying scoping rules. The user must ensure that the private data block, that is being referenced, persists for the entire access period. This is not an issue in the previous example, where the private data block was defined on the scope of the main program. However, the Internal State construct is often useful inside of Component modules to hold Component specific data between calls. One option to ensure persisting private data blocks is to use the Fortran SAVE attribute either on local or module variables. A second option, illustrated in the following example, is to use Fortran pointers and user controlled memory management via allocate() and deallocate() calls.
One situation where the Internal State is useful is in the creation of ensembles of the same Component. In this case it can be tricky to distinguish which data, held in saved module variables, belongs to which ensemble member - especially if the ensemble members are executing on the same set of PETs. The Internal State solves this problem by providing a handle to instance specific data allocations.
module user_mod use ESMF_Mod implicit none ! module variables private ! Internal State Variables type testData sequence integer :: testValue ! scalar data real :: testScaling ! scalar data real, pointer :: testArray(:) ! array data end type type dataWrapper sequence type(testData), pointer :: p end type
contains !--------------------------------------------------------------------
subroutine mygcomp_init(gcomp, istate, estate, clock, rc) type(ESMF_GridComp):: gcomp type(ESMF_State):: istate, estate type(ESMF_Clock):: clock integer, intent(out):: rc ! Local variables type(dataWrapper) :: wrap type(testData), pointer :: data integer :: i rc = ESMF_SUCCESS ! Allocate private data block allocate(data) ! Initialize private data block data%testValue = 4567 ! initialize scalar data data%testScaling = 0.5 ! initialize scalar data allocate(data%testArray(10)) ! allocate array data do i=1, 10 data%testArray(i) = real(i) ! initialize array data enddo ! In a real ensemble application the initial data would be set to something ! unique for this ensemble member. This could be accomplished for example ! by reading a member specific config file that was specified by the ! driver code. Alternatively, Attributes, set by the driver, could be used ! to label the Component instances as specific ensemble members. ! Set Internal State wrap%p => data call ESMF_GridCompSetInternalState(gcomp, wrap, rc) end subroutine !-------------------------------------------------------------- subroutine mygcomp_run(gcomp, istate, estate, clock, rc) type(ESMF_GridComp):: gcomp type(ESMF_State):: istate, estate type(ESMF_Clock):: clock integer, intent(out):: rc ! Local variables type(dataWrapper) :: wrap type(testData), pointer :: data logical :: match = .true. integer :: i rc = ESMF_SUCCESS ! Get Internal State call ESMF_GridCompGetInternalState(gcomp, wrap, rc) if (rc/=ESMF_SUCCESS) return ! Access private data block and verify data data => wrap%p if (data%testValue .ne. 4567) match = .false. ! test scalar data if (data%testScaling .ne. 0.5) match = .false. ! test scalar data do i=1, 10 if (data%testArray(i) .ne. real(i)) match = .false. ! test array data enddo if (match) then print *, "got same values back from GetInternalState as original" else print *, "did not get same values back" rc = ESMF_FAILURE endif end subroutine !-------------------------------------------------------------- subroutine mygcomp_final(gcomp, istate, estate, clock, rc) type(ESMF_GridComp):: gcomp type(ESMF_State):: istate, estate type(ESMF_Clock):: clock integer, intent(out):: rc ! Local variables type(dataWrapper) :: wrap type(testData), pointer :: data rc = ESMF_SUCCESS ! Get Internal State call ESMF_GridCompGetInternalState(gcomp, wrap, rc) if (rc/=ESMF_SUCCESS) return ! Deallocate private data block data => wrap%p deallocate(data%testArray) ! deallocate array data deallocate(data) end subroutine !-------------------------------------------------------------- end module
INTERFACE:
recursive function ESMF_GridCompCreate(name, gridcomptype, grid, config, & configFile, clock, petList, contextflag, rc)RETURN VALUE:
type(ESMF_GridComp) :: ESMF_GridCompCreateARGUMENTS:
character(len=*), intent(in), optional :: name type(ESMF_GridCompType), intent(in), optional :: gridcomptype type(ESMF_Grid), intent(inout), optional :: grid type(ESMF_Config), intent(inout), optional :: config character(len=*), intent(in), optional :: configFile type(ESMF_Clock), intent(inout), optional :: clock integer, intent(in), optional :: petList(:) type(ESMF_ContextFlag), intent(in), optional :: contextflag integer, intent(out), optional :: rcDESCRIPTION:
This interface creates an ESMF_GridComp object. By default, a separate VM context will be created for each component. This implies creating a new MPI communicator and allocating additional memory to manage the VM resources. When running on a large number of processors, creating a separate VM for each component could be both time and memory inefficient. If the application is sequential, i.e., each component is running on all the PETs of the global VM, it will be more efficient to use the global VM instead of creating a new one. This can be done by setting contextflag to ESMF_CHILD_IN_PARENT_VM.
The return value is the new ESMF_GridComp.
The arguments are:
INTERFACE:
subroutine ESMF_GridCompDestroy(gridcomp, rc)ARGUMENTS:
type(ESMF_GridComp) :: gridcomp integer, intent(out), optional :: rcDESCRIPTION:
Releases all resources associated with this ESMF_GridComp.
The arguments are:
INTERFACE:
recursive subroutine ESMF_GridCompFinalize(gridcomp, importState, & exportState, clock, blockingflag, phase, userRc, rc)ARGUMENTS:
type(ESMF_GridComp) :: gridcomp type(ESMF_State), intent(inout), optional :: importState type(ESMF_State), intent(inout), optional :: exportState type(ESMF_Clock), intent(inout), optional :: clock type(ESMF_BlockingFlag), intent(in), optional :: blockingflag integer, intent(in), optional :: phase integer, intent(out), optional :: userRc integer, intent(out), optional :: rcDESCRIPTION:
Call the associated user-supplied finalization code for an ESMF_GridComp.
The arguments are:
INTERFACE:
subroutine ESMF_GridCompGet(gridcomp, name, gridcomptype, grid, config, & configFile, clock, localPet, petCount, contextflag, currentMethod, & currentPhase, comptype, vm, rc)ARGUMENTS:
type(ESMF_GridComp), intent(inout) :: gridcomp character(len=*), intent(out), optional :: name type(ESMF_GridCompType), intent(out), optional :: gridcomptype type(ESMF_Grid), intent(out), optional :: grid type(ESMF_Config), intent(out), optional :: config character(len=*), intent(out), optional :: configFile type(ESMF_Clock), intent(out), optional :: clock integer, intent(out), optional :: localPet integer, intent(out), optional :: petCount type(ESMF_ContextFlag), intent(out), optional :: contextflag type(ESMF_Method), intent(out), optional :: currentMethod integer, intent(out), optional :: currentPhase type(ESMF_CompType), intent(out), optional :: comptype type(ESMF_VM), intent(out), optional :: vm integer, intent(out), optional :: rcDESCRIPTION:
Returns information about an ESMF_GridComp. For queries where the caller only wants a single value, specify the argument by name. All the arguments after the gridcomp argument are optional to facilitate this.
The arguments are:
INTERFACE:
subroutine ESMF_GridCompGetInternalState(gridcomp, dataPointer, rc)ARGUMENTS:
type(ESMF_GridComp), intent(inout) :: gridcomp type(any), pointer :: dataPointer integer, intent(out) :: rcDESCRIPTION:
Available to be called by an ESMF_GridComp at any time after ESMF_GridCompSetInternalState has been called. Since init, run, and finalize must be separate subroutines, data that they need to share in common can either be module global data, or can be allocated in a private data block and the address of that block can be registered with the framework and retrieved by this call. When running multiple instantiations of an ESMF_GridComp, for example during ensemble runs, it may be simpler to maintain private data specific to each run with private data blocks. A corresponding ESMF_GridCompSetInternalState call sets the data pointer to this block, and this call retrieves the data pointer. Note that the dataPointer argument needs to be a derived type which contains only a pointer of the type of the data block defined by the user. When making this call the pointer needs to be unassociated. When the call returns, the pointer will now reference the original data block which was set during the previous call to ESMF_GridCompSetInternalState.
Only the last data block set via ESMF_GridCompSetInternalState will be accessible.
The arguments are:
INTERFACE:
recursive subroutine ESMF_GridCompInitialize(gridcomp, importState, & exportState, clock, blockingflag, phase, userRc, rc)ARGUMENTS:
type(ESMF_GridComp) :: gridcomp type(ESMF_State), intent(inout), optional :: importState type(ESMF_State), intent(inout), optional :: exportState type(ESMF_Clock), intent(inout), optional :: clock type(ESMF_BlockingFlag), intent(in), optional :: blockingflag integer, intent(in), optional :: phase integer, intent(out), optional :: userRc integer, intent(out), optional :: rcDESCRIPTION:
Call the associated user initialization code for a GridComp.
The arguments are:
INTERFACE:
recursive function ESMF_GridCompIsPetLocal(gridcomp, rc)RETURN VALUE:
logical :: ESMF_GridCompIsPetLocalARGUMENTS:
type(ESMF_GridComp), intent(inout) :: gridcomp integer, intent(out), optional :: rcDESCRIPTION:
Inquire if this ESMF_GridComp object is to execute on the calling PET.
The return value is .true. if the component is to execute on the calling PET, .false. otherwise.
The arguments are:
INTERFACE:
subroutine ESMF_GridCompPrint(gridcomp, options, rc)ARGUMENTS:
type(ESMF_GridComp) :: gridcomp character(len = *), intent(in), optional :: options integer, intent(out), optional :: rcDESCRIPTION:
Prints information about an ESMF_GridComp to stdout.
Note: Many ESMF_<class>Print methods are implemented in C++.
On some platforms/compilers there is a potential issue with interleaving
Fortran and C++ output to stdout such that it doesn't appear in
the expected order. If this occurs, the ESMF_IOUnitFlush() method
may be used on unit 6 to get coherent output.
The arguments are:
INTERFACE:
recursive subroutine ESMF_GridCompReadRestart(gridcomp, importState, & exportState, clock, blockingflag, phase, userRc, rc)ARGUMENTS:
type(ESMF_GridComp) :: gridcomp type(ESMF_State), intent(inout), optional :: importState type(ESMF_State), intent(inout), optional :: exportState type(ESMF_Clock), intent(inout), optional :: clock type(ESMF_BlockingFlag), intent(in), optional :: blockingflag integer, intent(in), optional :: phase integer, intent(out), optional :: userRc integer, intent(out), optional :: rcDESCRIPTION:
Call the associated user read restart code for an ESMF_GridComp.
The arguments are:
INTERFACE:
recursive subroutine ESMF_GridCompRun(gridcomp, importState, exportState,& clock, blockingflag, phase, userRc, rc)ARGUMENTS:
type(ESMF_GridComp) :: gridcomp type(ESMF_State), intent(inout), optional :: importState type(ESMF_State), intent(inout), optional :: exportState type(ESMF_Clock), intent(inout), optional :: clock type(ESMF_BlockingFlag), intent(in), optional :: blockingflag integer, intent(in), optional :: phase integer, intent(out), optional :: userRc integer, intent(out), optional :: rcDESCRIPTION:
Call the associated user run code for an ESMF_GridComp.
The arguments are:
INTERFACE:
subroutine ESMF_GridCompSet(gridcomp, name, gridcomptype, grid, config, & configFile, clock, rc)ARGUMENTS:
type(ESMF_GridComp), intent(inout) :: gridcomp character(len=*), intent(in), optional :: name type(ESMF_GridCompType), intent(in), optional :: gridcomptype type(ESMF_Grid), intent(inout), optional :: grid type(ESMF_Config), intent(inout), optional :: config character(len=*), intent(in), optional :: configFile type(ESMF_Clock), intent(inout), optional :: clock integer, intent(out), optional :: rcDESCRIPTION:
Sets or resets information about an ESMF_GridComp. The caller can set individual values by specifying the arguments by name. All the arguments except gridcomp are optional to facilitate this.
The arguments are:
INTERFACE:
subroutine ESMF_GridCompSetEntryPoint(gridcomp, method, userRoutine, phase, rc)ARGUMENTS:
type(ESMF_GridComp), intent(in) :: gridcomp type(ESMF_Method), intent(in) :: method 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(in), optional :: phase integer, intent(out), optional :: rcDESCRIPTION:
Registers a user-supplied userRoutine as the entry point for one of the predefined Component methods. After this call the userRoutine becomes accessible via the standard Component method API.
The arguments are:
The Component writer must supply a subroutine with the exact interface shown above for the userRoutine argument. Arguments in userRoutine must not be declared as optional, and the types, intent and order must match.
INTERFACE:
subroutine ESMF_GridCompSetInternalState(gridcomp, dataPointer, rc)ARGUMENTS:
type(ESMF_GridComp), intent(inout) :: gridcomp type(any), pointer :: dataPointer integer, intent(out) :: rcDESCRIPTION:
Available to be called by an ESMF_GridComp at any time, but expected to be most useful when called during the registration process, or initialization. Since init, run, and finalize must be separate subroutines, data that they need to share in common can either be module global data, or can be allocated in a private data block and the address of that block can be registered with the framework and retrieved by subsequent calls. When running multiple instantiations of an ESMF_GridComp, for example during ensemble runs, it may be simpler to maintain private data specific to each run with private data blocks. A corresponding ESMF_GridCompGetInternalState call retrieves the data pointer.
Only the last data block set via ESMF_GridCompSetInternalState will be accessible.
The arguments are:
INTERFACE:
recursive subroutine ESMF_GridCompSetServices(gridcomp, userRoutine, userRc, rc)ARGUMENTS:
type(ESMF_GridComp) :: gridcomp interface subroutine userRoutine(gridcomp, rc) use ESMF_CompMod 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 :: userRc integer, intent(out), optional :: rcDESCRIPTION:
Call into user provided userRoutine which is responsible for setting Component's Initialize(), Run() and Finalize() services.
The arguments are:
The Component writer must supply a subroutine with the exact interface shown above for the userRoutine argument. Arguments in userRoutine must not be declared as optional, and the types, intent and order must match.
The userRoutine, when called by the framework, must make successive calls to ESMF_GridCompSetEntryPoint() to preset callback routines for standard Component Initialize(), Run() and Finalize() methods.
INTERFACE:
! Private name; call using ESMF_GridCompSetServices() recursive subroutine ESMF_GridCompSetServicesShObj(gridcomp, userRoutine, & sharedObj, userRc, rc)ARGUMENTS:
type(ESMF_GridComp), intent(inout) :: gridcomp character(len=*), intent(in) :: userRoutine character(len=*), intent(in), optional :: sharedObj integer, intent(out), optional :: userRc integer, intent(out), optional :: rcDESCRIPTION:
Call into user provided routine which is responsible for setting Component's Initialize(), Run() and Finalize() services. The named userRoutine must exist in the shared object file specified in the sharedObj argument. All of the platform specific details about dynamic linking and loading apply.
The arguments are:
The Component writer must supply a subroutine with the exact interface shown for userRoutine below. Arguments must not be declared as optional, and the types, intent and order must match.
INTERFACE:
interface subroutine userRoutine(gridcomp, rc) type(ESMF_GridComp) :: gridcomp ! must not be optional integer, intent(out) :: rc ! must not be optional end subroutine end interfaceDESCRIPTION:
The userRoutine, when called by the framework, must make successive calls to ESMF_GridCompSetEntryPoint() to preset callback routines for standard Component Initialize(), Run() and Finalize() methods.
INTERFACE:
recursive subroutine ESMF_GridCompSetVM(gridcomp, userRoutine, userRc, rc)ARGUMENTS:
type(ESMF_GridComp) :: gridcomp interface subroutine userRoutine(gridcomp, rc) use ESMF_CompMod 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 :: userRc integer, intent(out), optional :: rcDESCRIPTION:
Optionally call into user provided userRoutine which is responsible for for setting Component's VM properties.
The arguments are:
The Component writer must supply a subroutine with the exact interface shown above for the userRoutine argument. Arguments in userRoutine must not be declared as optional, and the types, intent and order must match.
The subroutine, when called by the framework, is expected to use any of the ESMF_GridCompSetVMxxx() methods to set the properties of the VM associated with the Gridded Component.
INTERFACE:
! Private name; call using ESMF_GridCompSetVM() recursive subroutine ESMF_GridCompSetVMShObj(gridcomp, userRoutine, sharedObj, & userRc, rc)ARGUMENTS:
type(ESMF_GridComp), intent(inout) :: gridcomp character(len=*), intent(in) :: userRoutine character(len=*), intent(in), optional :: sharedObj integer, intent(out), optional :: userRc integer, intent(out), optional :: rcDESCRIPTION:
Optionally call into user provided userRoutine which is responsible for for setting Component's VM properties. The named userRoutine must exist in the shared object file specified in the sharedObj argument. All of the platform specific details about dynamic linking and loading apply.
The arguments are:
The Component writer must supply a subroutine with the exact interface shown for userRoutine below. Arguments must not be declared as optional, and the types, intent and order must match.
INTERFACE:
interface subroutine userRoutine(gridcomp, rc) type(ESMF_GridComp) :: gridcomp ! must not be optional integer, intent(out) :: rc ! must not be optional end subroutine end interfaceDESCRIPTION:
The subroutine, when called by the framework, is expected to use any of the ESMF_GridCompSetVMxxx() methods to set the properties of the VM associated with the Gridded Component.
INTERFACE:
subroutine ESMF_GridCompSetVMMaxPEs(gridcomp, max, pref_intra_process, & pref_intra_ssi, pref_inter_ssi, rc)ARGUMENTS:
type(ESMF_GridComp), intent(inout) :: gridcomp integer, intent(in), optional :: max integer, intent(in), optional :: pref_intra_process integer, intent(in), optional :: pref_intra_ssi integer, intent(in), optional :: pref_inter_ssi integer, intent(out), optional :: rcDESCRIPTION:
Set characteristics of the ESMF_VM for this ESMF_GridComp. Attempts to associate max PEs with each PET. Only PEs that are located on the same single system image can be associated with the same PET. Within this constraint the call tries to get as close as possible to the number specified by max.
The typical use of ESMF_GridCompSetVMMaxPEs() is to allocate multiple PEs per PET in a Component for user-level threading, e.g. OpenMP.
The arguments are:
INTERFACE:
subroutine ESMF_GridCompSetVMMaxThreads(gridcomp, max, pref_intra_process, & pref_intra_ssi, pref_inter_ssi, rc)ARGUMENTS:
type(ESMF_GridComp), intent(inout) :: gridcomp integer, intent(in), optional :: max integer, intent(in), optional :: pref_intra_process integer, intent(in), optional :: pref_intra_ssi integer, intent(in), optional :: pref_inter_ssi integer, intent(out), optional :: rcDESCRIPTION:
Set characteristics of the ESMF_VM for this ESMF_GridComp. Attempts to provide max threaded PETs in each VAS. Only as many threaded PETs as there are PEs located on the same single system image can be associated with the same VAS. Within this constraint the call tries to get as close as possible to the number specified by max.
The typical use of ESMF_GridCompSetVMMaxThreads() is to run a Component multi-threaded with a groups of PETs that execute within the same virtual address space.
The arguments are:
INTERFACE:
subroutine ESMF_GridCompSetVMMinThreads(gridcomp, max, pref_intra_process, & pref_intra_ssi, pref_inter_ssi, rc)ARGUMENTS:
type(ESMF_GridComp), intent(inout) :: gridcomp integer, intent(in), optional :: max integer, intent(in), optional :: pref_intra_process integer, intent(in), optional :: pref_intra_ssi integer, intent(in), optional :: pref_inter_ssi integer, intent(out), optional :: rcDESCRIPTION:
Set characteristics of the ESMF_VM for this ESMF_GridComp. Reduces the number of threaded PETs in each VAS. The max argument may be specified to limit the maximum number of PEs that a single PET may be associated with.
The typical use of ESMF_GridCompSetVMMinThreads() is to run a Component across a set of single-threaded PETs.
The arguments are:
INTERFACE:
subroutine ESMF_GridCompValidate(gridcomp, options, rc)ARGUMENTS:
type(ESMF_GridComp) :: gridcomp character(len = *), intent(in), optional :: options integer, intent(out), optional :: rcDESCRIPTION:
Currently all this method does is to check that the gridcomp exists.
The arguments are:
INTERFACE:
subroutine ESMF_GridCompWait(gridcomp, blockingflag, userRc, rc)ARGUMENTS:
type(ESMF_GridComp), intent(inout) :: gridcomp type(ESMF_BlockingFlag), intent(in), optional :: blockingflag integer, intent(out), optional :: userRc integer, intent(out), optional :: rcDESCRIPTION:
When executing asychronously, wait for an ESMF_GridComp to return.
The arguments are:
INTERFACE:
recursive subroutine ESMF_GridCompWriteRestart(gridcomp, importState, & exportState, clock, blockingflag, phase, userRc, rc)ARGUMENTS:
type(ESMF_GridComp) :: gridcomp type(ESMF_State), intent(inout), optional :: importState type(ESMF_State), intent(inout), optional :: exportState type(ESMF_Clock), intent(inout), optional :: clock type(ESMF_BlockingFlag), intent(in), optional :: blockingflag integer, intent(in), optional :: phase integer, intent(out), optional :: userRc integer, intent(out), optional :: rcDESCRIPTION:
Call the associated user write restart code for an ESMF_GridComp.
The arguments are:
In a large, multi-component application such as a weather forecasting or climate prediction system running within ESMF, physical domains and major system functions are represented as Gridded Components (see Section 15.1). A Coupler Component, or ESMF_CplComp, arranges and executes the data transformations between the Gridded Components. Ideally, Coupler Components should contain all the information about inter-component communication for an application. This enables the Gridded Components in the application to be used in multiple contexts; that is, used in different coupled configurations without changes to their source code. For example, the same atmosphere might in one case be coupled to an ocean in a hurricane prediction model, and to a data assimilation system for numerical weather prediction in another. A single Coupler Component can couple two or more Gridded Components.
Like Gridded Components, Coupler Components have two parts, one that is provided by the user and another that is part of the framework. The user-written portion of the software is the coupling code necessary for a particular exchange between Gridded Components. This portion of the Coupler Component code must be divided into separately callable initialize, run, and finalize methods. The interfaces for these methods are prescribed by ESMF.
The term ``user-written'' is somewhat misleading here, since within a Coupler Component the user can leverage ESMF infrastructure software for regridding, redistribution, lower-level communications, calendar management, and other functions. However, ESMF is unlikely to offer all the software necessary to customize a data transfer between Gridded Components. For instance, ESMF does not currently offer tools for unit tranformations or time averaging operations, so users must manage those operations themselves.
The second part of a Coupler Component is the ESMF_CplComp derived type within ESMF. The user must create one of these types to represent a specific coupling function, such as the regular transfer of data between a data assimilation system and an atmospheric model. 2
The user-written part of a Coupler Component is associated with an ESMF_CplComp derived type through a routine called ESMF_SetServices(). This is a routine that the user must write and declare public. Inside the ESMF_SetServices() routine the user must call ESMF_SetEntryPoint() methods that associate a standard ESMF operation with the name of the corresponding Fortran subroutine in their user code. For example, a user routine called ``couplerInit'' might be associated with the standard initialize routine in a Coupler Component.
A Coupler Component manages the transformation of data between Components. It contains a list of State objects and the operations needed to make them compatible, including such things as regridding and unit conversion. Coupler Components are user-written, following prescribed ESMF interfaces and, wherever desired, using ESMF infrastructure tools.
Every ESMF_CplComp is required to provide and document a public set services routine. It can have any name, but must follow the declaration below: a subroutine which takes an ESMF_CplComp as the first argument, and an integer return code as the second. Both arguments are required and must not be declared as optional. If an intent is specified in the interface it must be intent(inout) for the first and intent(out) for the second argument.
The set services routine must call the ESMF method ESMF_CplCompSetEntryPoint() to register with the framework what user-code subroutines should be called to initialize, run, and finalize the component. There are additional routines which can be registered as well, for checkpoint and restart functions.
Note that the actual subroutines being registered do not have to be public to this module; only the set services routine itself must be available to be used by other code.
! Example Coupler Component module ESMF_CouplerEx ! ESMF Framework module use ESMF_Mod implicit none public CPL_SetServices contains subroutine CPL_SetServices(comp, rc) type(ESMF_CplComp) :: comp ! must not be optional integer, intent(out) :: rc ! must not be optional ! Set the entry points for standard ESMF Component methods call ESMF_CplCompSetEntryPoint(comp, ESMF_SETINIT, userRoutine=CPL_Init, rc=rc) call ESMF_CplCompSetEntryPoint(comp, ESMF_SETRUN, userRoutine=CPL_Run, rc=rc) call ESMF_CplCompSetEntryPoint(comp, ESMF_SETFINAL, userRoutine=CPL_Final, rc=rc) rc = ESMF_SUCCESS end subroutine
When a higher level component is ready to begin using an ESMF_CplComp, it will call its initialize routine.
The component writer must supply a subroutine with the exact interface shown below. Arguments must not be declared as optional, and the types and order must match.
At initialization time the component can allocate data space, open data files, set up initial conditions; anything it needs to do to prepare to run.
The rc return code should be set if an error occurs, otherwise the value ESMF_SUCCESS should be returned.
subroutine CPL_Init(comp, importState, exportState, clock, rc) type(ESMF_CplComp) :: comp ! must not be optional type(ESMF_State) :: importState, exportState ! must not be optional type(ESMF_Clock) :: clock ! must not be optional integer, intent(out) :: rc ! must not be optional print *, "Coupler Init starting" ! Add whatever code here needed ! Precompute any needed values, fill in any inital values ! needed in Import States rc = ESMF_SUCCESS print *, "Coupler Init returning" end subroutine CPL_Init
During the execution loop, the run routine may be called many times. Each time it should read data from the importState, use the clock to determine what the current time is in the calling component, compute new values or process the data, and produce any output and place it in the exportState.
When a higher level component is ready to use the ESMF_CplComp it will call its run routine.
The component writer must supply a subroutine with the exact interface shown below. Arguments must not be declared as optional, and the types and order must match.
It is expected that this is where the bulk of the model computation or data analysis will occur.
The rc return code should be set if an error occurs, otherwise the value ESMF_SUCCESS should be returned.
subroutine CPL_Run(comp, importState, exportState, clock, rc) type(ESMF_CplComp) :: comp ! must not be optional type(ESMF_State) :: importState, exportState ! must not be optional type(ESMF_Clock) :: clock ! must not be optional integer, intent(out) :: rc ! must not be optional print *, "Coupler Run starting" ! Add whatever code needed here to transform Export state data ! into Import states for the next timestep. rc = ESMF_SUCCESS print *, "Coupler Run returning" end subroutine CPL_Run
At the end of application execution, each ESMF_CplComp should deallocate data space, close open files, and flush final results. These functions should be placed in a finalize routine.
The component writer must supply a subroutine with the exact interface shown below. Arguments must not be declared as optional, and the types and order must match.
The rc return code should be set if an error occurs, otherwise the value ESMF_SUCCESS should be returned.
subroutine CPL_Final(comp, importState, exportState, clock, rc) type(ESMF_CplComp) :: comp ! must not be optional type(ESMF_State) :: importState, exportState ! must not be optional type(ESMF_Clock) :: clock ! must not be optional integer, intent(out) :: rc ! must not be optional print *, "Coupler Final starting" ! Add whatever code needed here to compute final values and ! finish the computation. rc = ESMF_SUCCESS print *, "Coupler Final returning" end subroutine CPL_Final
Every ESMF_CplComp can optionally provide and document a public set vm routine. It can have any name, but must follow the declaration below: a subroutine which takes an ESMF_CplComp as the first argument, and an integer return code as the second. Both arguments are required and must not be declared as optional. If an intent is specified in the interface it must be intent(inout) for the first and intent(out) for the second argument.
The set vm routine is the only place where the child component can use the ESMF_CplCompSetVMMaxPEs(), or ESMF_CplCompSetVMMaxThreads(), or ESMF_CplCompSetVMMinThreads() call to modify aspects of its own VM.
A component's VM is started up right before its set services routine is entered. ESMF_CplCompSetVM() is executing in the parent VM, and must be called before ESMF_CplCompSetServices().
subroutine GComp_SetVM(comp, rc) type(ESMF_CplComp) :: comp ! must not be optional integer, intent(out) :: rc ! must not be optional type(ESMF_VM) :: vm logical :: pthreadsEnabled ! Test for Pthread support, all SetVM calls require it call ESMF_VMGetGlobal(vm, rc=rc) call ESMF_VMGet(vm, pthreadsEnabledFlag=pthreadsEnabled, rc=rc) if (pthreadsEnabled) then ! run PETs single-threaded call ESMF_CplCompSetVMMinThreads(comp, rc=rc) endif rc = ESMF_SUCCESS end subroutine end module ESMF_CouplerEx
INTERFACE:
recursive function ESMF_CplCompCreate(name, config, configFile, clock, & petList, contextflag, rc)RETURN VALUE:
type(ESMF_CplComp) :: ESMF_CplCompCreateARGUMENTS:
character(len=*), intent(in), optional :: name type(ESMF_Config), intent(inout), optional :: config character(len=*), intent(in), optional :: configFile type(ESMF_Clock), intent(inout), optional :: clock integer, intent(in), optional :: petList(:) type(ESMF_ContextFlag), intent(in), optional :: contextflag integer, intent(out), optional :: rcDESCRIPTION:
This interface creates an ESMF_CplComp object. By default, a separate VM context will be created for each component. This implies creating a new MPI communicator and allocating additional memory to manage the VM resources. When running on a large number of processors, creating a separate VM for each component could be both time and memory inefficient. If the application is sequential, i.e., each component is running on all the PETs of the global VM, it will be more efficient to use the global VM instead of creating a new one. This can be done by setting contextflag to ESMF_CHILD_IN_PARENT_VM.
The return value is the new ESMF_CplComp.
The arguments are:
INTERFACE:
subroutine ESMF_CplCompDestroy(cplcomp, rc)ARGUMENTS:
type(ESMF_CplComp) :: cplcomp integer, intent(out), optional :: rcDESCRIPTION:
Releases all resources associated with this ESMF_CplComp.
The arguments are:
INTERFACE:
recursive subroutine ESMF_CplCompFinalize(cplcomp, importState, exportState, & clock, blockingflag, phase, userRc, rc)ARGUMENTS:
type(ESMF_CplComp) :: cplcomp type(ESMF_State), intent(inout), optional :: importState type(ESMF_State), intent(inout), optional :: exportState type(ESMF_Clock), intent(inout), optional :: clock type(ESMF_BlockingFlag), intent(in), optional :: blockingflag integer, intent(in), optional :: phase integer, intent(out), optional :: userRc integer, intent(out), optional :: rcDESCRIPTION:
Call the associated user-supplied finalization routine for an ESMF_CplComp.
The arguments are:
INTERFACE:
subroutine ESMF_CplCompGet(cplcomp, name, config, configFile, clock, & localPet, petCount, contextflag, currentMethod, currentPhase, vm, rc)ARGUMENTS:
type(ESMF_CplComp), intent(inout) :: cplcomp character(len=*), intent(out), optional :: name type(ESMF_Config), intent(out), optional :: config character(len=*), intent(out), optional :: configFile type(ESMF_Clock), intent(out), optional :: clock integer, intent(out), optional :: localPet integer, intent(out), optional :: petCount type(ESMF_ContextFlag), intent(out), optional :: contextflag type(ESMF_Method), intent(out), optional :: currentMethod integer, intent(out), optional :: currentPhase type(ESMF_VM), intent(out), optional :: vm integer, intent(out), optional :: rcDESCRIPTION:
Returns information about an ESMF_CplComp. For queries where the caller only wants a single value, specify the argument by name. All the arguments after cplcomp argument are optional to facilitate this.
The arguments are:
INTERFACE:
subroutine ESMF_CplCompGetInternalState(cplcomp, dataPointer, rc)ARGUMENTS:
type(ESMF_CplComp), intent(inout) :: cplcomp type(any), pointer :: dataPointer integer, intent(out) :: rcDESCRIPTION:
Available to be called by an ESMF_CplComp at any time after ESMF_CplCompSetInternalState has been called. Since init, run, and finalize must be separate subroutines, data that they need to share in common can either be module global data, or can be allocated in a private data block and the address of that block can be registered with the framework and retrieved by this call. When running multiple instantiations of an ESMF_CplComp, for example during ensemble runs, it may be simpler to maintain private data specific to each run with private data blocks. A corresponding ESMF_CplCompSetInternalState call sets the data pointer to this block, and this call retrieves the data pointer. Note that the dataPointer argument needs to be a derived type which contains only a pointer of the type of the data block defined by the user. When making this call the pointer needs to be unassociated. When the call returns, the pointer will now reference the original data block which was set during the previous call to ESMF_CplCompSetInternalState.
Only the last data block set via ESMF_CplCompSetInternalState will be accessible.
The arguments are:
INTERFACE:
recursive subroutine ESMF_CplCompInitialize(cplcomp, importState, & exportState, clock, blockingflag, phase, userRc, rc)ARGUMENTS:
type(ESMF_CplComp) :: cplcomp type(ESMF_State), intent(inout), optional :: importState type(ESMF_State), intent(inout), optional :: exportState type(ESMF_Clock), intent(inout), optional :: clock type(ESMF_BlockingFlag), intent(in), optional :: blockingflag integer, intent(in), optional :: phase integer, intent(out), optional :: userRc integer, intent(out), optional :: rcDESCRIPTION:
Call the associated user initialization code for a CplComp.
The arguments are:
INTERFACE:
recursive function ESMF_CplCompIsPetLocal(cplcomp, rc)RETURN VALUE:
logical :: ESMF_CplCompIsPetLocalARGUMENTS:
type(ESMF_CplComp), intent(inout) :: cplcomp integer, intent(out), optional :: rcDESCRIPTION:
Inquire if this ESMF_CplComp object is to execute on the calling PET.
The return value is .true. if the component is to execute on the calling PET, .false. otherwise.
The arguments are:
INTERFACE:
subroutine ESMF_CplCompPrint(cplcomp, options, rc)ARGUMENTS:
type(ESMF_CplComp) :: cplcomp character(len = *), intent(in), optional :: options integer, intent(out), optional :: rcDESCRIPTION:
Prints information about an ESMF_CplComp to stdout.
Note: Many ESMF_<class>Print methods are implemented in C++.
On some platforms/compilers there is a potential issue with interleaving
Fortran and C++ output to stdout such that it doesn't appear in
the expected order. If this occurs, the ESMF_IOUnitFlush() method
may be used on unit 6 to get coherent output.
The arguments are:
INTERFACE:
recursive subroutine ESMF_CplCompReadRestart(cplcomp, importState, & exportState, clock, blockingflag, phase, userRc, rc)ARGUMENTS:
type(ESMF_CplComp) :: cplcomp type(ESMF_State), intent(inout), optional :: importState type(ESMF_State), intent(inout), optional :: exportState type(ESMF_Clock), intent(inout), optional :: clock type(ESMF_BlockingFlag), intent(in), optional :: blockingflag integer, intent(in), optional :: phase integer, intent(out), optional :: userRc integer, intent(out), optional :: rcDESCRIPTION:
Call the associated user read restart code for an ESMF_CplComp.
The arguments are:
INTERFACE:
recursive subroutine ESMF_CplCompRun(cplcomp, importState, exportState, & clock, blockingflag, phase, userRc, rc)ARGUMENTS:
type(ESMF_CplComp) :: cplcomp type(ESMF_State), intent(inout), optional :: importState type(ESMF_State), intent(inout), optional :: exportState type(ESMF_Clock), intent(inout), optional :: clock type(ESMF_BlockingFlag), intent(in), optional :: blockingflag integer, intent(in), optional :: phase integer, intent(out), optional :: userRc integer, intent(out), optional :: rcDESCRIPTION:
Call the associated user run code for an ESMF_CplComp.
The arguments are:
INTERFACE:
subroutine ESMF_CplCompSet(cplcomp, name, config, configFile, clock, rc)ARGUMENTS:
type(ESMF_CplComp), intent(inout) :: cplcomp character(len=*), intent(in), optional :: name type(ESMF_Config), intent(inout), optional :: config character(len=*), intent(in), optional :: configFile type(ESMF_Clock), intent(inout), optional :: clock integer, intent(out), optional :: rcDESCRIPTION:
Sets or resets information about an ESMF_CplComp. The caller can set individual values by specifying the arguments by name. All the arguments except cplcomp are optional to facilitate this.
The arguments are:
INTERFACE:
subroutine ESMF_CplCompSetEntryPoint(cplcomp, method, userRoutine, phase, rc)ARGUMENTS:
type(ESMF_CplComp), intent (in) :: cplcomp type(ESMF_Method), intent(in) :: method 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(in), optional :: phase integer, intent(out), optional :: rcDESCRIPTION:
Registers a user-supplied userRoutine as the entry point for one of the predefined Component methods. After this call the userRoutine becomes accessible via the standard Component method API.
The arguments are:
The Component writer must supply a subroutine with the exact interface shown above for the userRoutine argument. Arguments in userRoutine must not be declared as optional, and the types, intent and order must match.
INTERFACE:
subroutine ESMF_CplCompSetInternalState(cplcomp, dataPointer, rc)ARGUMENTS:
type(ESMF_CplComp), intent(inout) :: cplcomp type(any), pointer :: dataPointer integer, intent(out) :: rcDESCRIPTION:
Available to be called by an ESMF_CplComp at any time, but expected to be most useful when called during the registration process, or initialization. Since init, run, and finalize must be separate subroutines data that they need to share in common can either be module global data, or can be allocated in a private data block and the address of that block can be registered with the framework and retrieved by subsequent calls. When running multiple instantiations of an ESMF_CplComp, for example during ensemble runs, it may be simpler to maintain private data specific to each run with private data blocks. A corresponding ESMF_CplCompGetInternalState call retrieves the data pointer.
Only the last data block set via ESMF_CplCompSetInternalState will be accessible.
The arguments are:
INTERFACE:
recursive subroutine ESMF_CplCompSetServices(cplcomp, userRoutine, userRc, rc)ARGUMENTS:
type(ESMF_CplComp) :: cplcomp interface subroutine userRoutine(cplcomp, rc) use ESMF_CompMod 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 :: userRc integer, intent(out), optional :: rcDESCRIPTION:
Call into user provided userRoutine which is responsible for for setting Component's Initialize(), Run() and Finalize() services.
The arguments are:
The Component writer must supply a subroutine with the exact interface shown above for the userRoutine argument. Arguments in userRoutine must not be declared as optional, and the types, intent and order must match.
The userRoutine, when called by the framework, must make successive calls to ESMF_CplCompSetEntryPoint() to preset callback routines for standard Component Initialize(), Run() and Finalize() methods.
INTERFACE:
! Private name; call using ESMF_CplCompSetServices() recursive subroutine ESMF_CplCompSetServicesShObj(cplcomp, userRoutine, & sharedObj, userRc, rc)ARGUMENTS:
type(ESMF_CplComp), intent(inout) :: cplcomp character(len=*), intent(in) :: userRoutine character(len=*), intent(in), optional :: sharedObj integer, intent(out), optional :: userRc integer, intent(out), optional :: rcDESCRIPTION:
Call into user provided routine which is responsible for setting Component's Initialize(), Run() and Finalize() services. The named userRoutine must exist in the shared object file specified in the sharedObj argument. All of the platform specific details about dynamic linking and loading apply.
The arguments are:
The Component writer must supply a subroutine with the exact interface shown for userRoutine below. Arguments must not be declared as optional, and the types, intent and order must match.
INTERFACE:
interface subroutine userRoutine(cplcomp, rc) type(ESMF_CplComp) :: cplcomp ! must not be optional integer, intent(out) :: rc ! must not be optional end subroutine end interfaceDESCRIPTION:
The userRoutine, when called by the framework, must make successive calls to ESMF_CplCompSetEntryPoint() to preset callback routines for standard Component Initialize(), Run() and Finalize() methods.
INTERFACE:
recursive subroutine ESMF_CplCompSetVM(cplcomp, userRoutine, userRc, rc)ARGUMENTS:
type(ESMF_CplComp) :: cplcomp interface subroutine userRoutine(cplcomp, rc) use ESMF_CompMod 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 :: userRc integer, intent(out), optional :: rcDESCRIPTION:
Optionally call into user provided userRoutine which is responsible for for setting Component's VM properties.
The arguments are:
The Component writer must supply a subroutine with the exact interface shown above for the userRoutine argument. Arguments in userRoutine must not be declared as optional, and the types, intent and order must match.
The subroutine, when called by the framework, is expected to use any of the ESMF_CplCompSetVMxxx() methods to set the properties of the VM associated with the Coupler Component.
INTERFACE:
! Private name; call using ESMF_CplCompSetVM() recursive subroutine ESMF_CplCompSetVMShObj(cplcomp, userRoutine, sharedObj, & userRc, rc)ARGUMENTS:
type(ESMF_CplComp), intent(inout) :: cplcomp character(len=*), intent(in) :: userRoutine character(len=*), intent(in), optional :: sharedObj integer, intent(out), optional :: userRc integer, intent(out), optional :: rcDESCRIPTION:
Optionally call into user provided userRoutine which is responsible for for setting Component's VM properties. The named userRoutine must exist in the shared object file specified in the sharedObj argument. All of the platform specific details about dynamic linking and loading apply.
The arguments are:
The Component writer must supply a subroutine with the exact interface shown for userRoutine below. Arguments must not be declared as optional, and the types, intent and order must match.
INTERFACE:
interface subroutine userRoutine(cplcomp, rc) type(ESMF_CplComp) :: cplcomp ! must not be optional integer, intent(out) :: rc ! must not be optional end subroutine end interfaceDESCRIPTION:
The subroutine, when called by the framework, is expected to use any of the ESMF_CplCompSetVMxxx() methods to set the properties of the VM associated with the Coupler Component.
INTERFACE:
subroutine ESMF_CplCompSetVMMaxPEs(cplcomp, max, pref_intra_process, & pref_intra_ssi, pref_inter_ssi, rc)ARGUMENTS:
type(ESMF_CplComp), intent(inout) :: cplcomp integer, intent(in), optional :: max integer, intent(in), optional :: pref_intra_process integer, intent(in), optional :: pref_intra_ssi integer, intent(in), optional :: pref_inter_ssi integer, intent(out), optional :: rcDESCRIPTION:
Set characteristics of the ESMF_VM for this ESMF_CplComp. Attempts to associate max PEs with each PET. Only PEs that are located on the same single system image can be associated with the same PET. Within this constraint the call tries to get as close as possible to the number specified by max.
The typical use of ESMF_CplCompSetVMMaxPEs() is to allocate multiple PEs per PET in a Component for user-level threading, e.g. OpenMP.
The arguments are:
INTERFACE:
subroutine ESMF_CplCompSetVMMaxThreads(cplcomp, max, pref_intra_process, & pref_intra_ssi, pref_inter_ssi, rc)ARGUMENTS:
type(ESMF_CplComp), intent(inout) :: cplcomp integer, intent(in), optional :: max integer, intent(in), optional :: pref_intra_process integer, intent(in), optional :: pref_intra_ssi integer, intent(in), optional :: pref_inter_ssi integer, intent(out), optional :: rcDESCRIPTION:
Set characteristics of the ESMF_VM for this ESMF_CplComp. Attempts to provide max threaded PETs in each VAS. Only as many threaded PETs as there are PEs located on the same single system image can be associated with the same VAS. Within this constraint the call tries to get as close as possible to the number specified by max.
The typical use of ESMF_CplCompSetVMMaxThreads() is to run a Component multi-threaded with a groups of PETs that execute within the same virtual address space.
The arguments are:
INTERFACE:
subroutine ESMF_CplCompSetVMMinThreads(cplcomp, max, pref_intra_process, & pref_intra_ssi, pref_inter_ssi, rc)ARGUMENTS:
type(ESMF_CplComp), intent(inout) :: cplcomp integer, intent(in), optional :: max integer, intent(in), optional :: pref_intra_process integer, intent(in), optional :: pref_intra_ssi integer, intent(in), optional :: pref_inter_ssi integer, intent(out), optional :: rcDESCRIPTION:
Set characteristics of the ESMF_VM for this ESMF_CplComp. Reduces the number of threaded PETs in each VAS. The max argument may be specified to limit the maximum number of PEs that a single PET may be associated with.
The typical use of ESMF_CplCompSetVMMinThreads() is to run a Component across a set of single-threaded PETs.
The arguments are:
INTERFACE:
subroutine ESMF_CplCompValidate(cplcomp, options, rc)ARGUMENTS:
type(ESMF_CplComp) :: cplcomp character(len = *), intent(in), optional :: options integer, intent(out), optional :: rcDESCRIPTION:
Currently all this method does is to check that the cplcomp exists.
The arguments are:
INTERFACE:
subroutine ESMF_CplCompWait(cplcomp, blockingflag, userRc, rc)ARGUMENTS:
type(ESMF_CplComp), intent(inout) :: cplcomp type(ESMF_BlockingFlag), intent(in), optional :: blockingflag integer, intent(out), optional :: userRc integer, intent(out), optional :: rcDESCRIPTION:
When executing asychronously, wait for an ESMF_CplComp to return.
The arguments are:
INTERFACE:
recursive subroutine ESMF_CplCompWriteRestart(cplcomp, importState, & exportState, clock, blockingflag, phase, userRc, rc)ARGUMENTS:
type(ESMF_CplComp), intent(inout) :: cplcomp type(ESMF_State), intent(inout), optional :: importState type(ESMF_State), intent(inout), optional :: exportState type(ESMF_Clock), intent(inout), optional :: clock type(ESMF_BlockingFlag), intent(in), optional :: blockingflag integer, intent(in), optional :: phase integer, intent(out), optional :: userRc integer, intent(out), optional :: rcDESCRIPTION:
Call the associated user write restart code for an ESMF_CplComp.
The arguments are:
A State contains the data and metadata to be transferred between ESMF Components. It is an important class, because it defines a standard for how data is represented in data transfers between Earth science components. The State construct is a rational compromise between a fully prescribed interface - one that would dictate what specific fields should be transferred between components - and an interface in which data structures are completely ad hoc.
There are two types of States, import and export. An import State contains data that is necessary for a Gridded Component or Coupler Component to execute, and an export State contains the data that a Gridded Component or Coupler Component can make available.
States can contain Arrays, ArrayBundles, Fields, FieldBundles, and other States. They cannot directly contain native language arrays (i.e. Fortran or C style arrays). Objects in a State must span the VM on which they are running. For sequentially executing components which run on the same set of PETs this happens by calling the object create methods on each PET, creating the object in unison. For concurrently executing components which are running on subsets of PETs, an additional method, called ESMF_StateReconcile(), is provided by ESMF to broadcast information about objects which were created in sub-components.
State methods include creation and deletion, adding and retrieving data items, adding and retrieving attributes, and performing queries.
DESCRIPTION:
Specifies the type of object being added to or retrieved from an
ESMF_State.
Valid values are:
Valid values are:
A Gridded Component generally has one associated import State and one export State. Generally the States associated with a Gridded Component will be created by the Gridded Component's parent component. In many cases, the States will be created containing no data. Both the empty States and the newly created Gridded Component are passed by the parent component into the Gridded Component's initialize method. This is where the States get prepared for use and the import State is first filled with data.
States can be created in a number of ways without the Fields, Arrays, FieldBundles, ArrayBundles, and other States they will eventually contain. They can be created with names as placeholders where these data items will eventually be. When the States are passed into the Gridded Component's initialize method, Field, FieldBundle, Array, and ArrayBundle create calls can be made in that method to replace the name placeholders with real data objects.
States can also be filled with data items that do not yet have data allocated. Fields, FieldBundles, Arrays, and ArrayBundles each have methods that support their creation without actual data allocation - the Grid and Attributes are set up but no Fortran array of data values is allocated. In this approach, when a State is passed into its associated Gridded Component's initialize method, the incomplete Arrays, Fields, FieldBundles, and ArrayBundles within the State can allocate or reference data inside the initialize method.
States are passed through the interfaces of the Gridded and Coupler Components' run methods in order to carry data between the components. While we expect a Gridded Component's import State to be filled with data during initialization, its export State will typically be filled over the course of its run method. At the end of a Gridded Component's run method, the filled export State is passed out through the argument list into a Coupler Component's run method. We recommend the convention that it enters the Coupler Component as the Coupler Component's import State. Here is it transformed into a form that another Gridded Component requires, and passed out of the Coupler Component as its export State. It can then be passed into the run method of a recipient Gridded Component as that component's import State.
While the above sounds complicated, the rule is simple: a State going into a component is an import State, and a State leaving a component is an export State.
Data items within a State can be marked needed or not needed, depending on whether they are required for a particular application configuration. If the item is marked not needed, the user can make the Gridded Component's initialize method clever enough to not allocate the data for that item at all and not compute it within the Gridded Component code. For example, some diagnostics may not be desired for all runs.
Other flags will eventually be available for data items within a State, such as data ready for reading or writing, data valid or invalid, and data required for restart or not. These are not yet fully implemented, so only the default value for each value can be set at this time.
Objects inside States are normally created in unison where each PET executing a component makes the same object create call. If the object contains data, like a Field, each PET may have a different local chunk of the entire dataset but each Field has the same name and is logically one part of a single distributed object. As States are passed between components, if any object in a State was not created in unison on all the current PETs then some PETs have no object to pass into a communication method (e.g. regrid or data redistribution). The ESMF_StateReconcile() method must be called to broadcast information about these objects to all PETs in a component; after which all PETs have a single uniform view of all objects and metadata.
If components are running in sequential mode on all available PETs and States are being passed between them there is no need to call ESMF_StateReconcile since all PETs have a uniform view of the objects. However, if components are running on a subset of the PETs, as is usually the case when running in concurrent mode, then when States are passed into components which contain a superset of those PETs, for example, a Coupler Component, all PETs must call ESMF_StateReconcile on the States before using them in any ESMF communication methods. The reconciliation process broadcasts information about objects which exist only on a subset of the PETs. On PETs missing those objects it creates a proxy object which contains any qualities of the original object plus enough information for it to be a data source or destination for a regrid or data redistribution operation. There is an option to turn off metadata reconciliation in the ESMF_StateReconcile call.
States can be created and destroyed at any time during application execution. The ESMF_StateCreate() routine can take many different combinations of optional arguments. Refer to the API description for all possible methods of creating a State. An empty State can be created by providing only a name and type for the intended State:
state = ESMF_StateCreate(statename, statetype=ESMF_STATE_IMPORT, rc=rc)
When finished with an ESMF_State, the ESMF_StateDestroy method removes it. However, the objects inside the ESMF_State created externally should be destroyed separately, since objects can be added to more than one ESMF_State.
Creation of an empty ESMF_State, and adding an ESMF_FieldBundle to it. Note that the ESMF_FieldBundle does not get destroyed when the ESMF_State is destroyed; the ESMF_State only contains a reference to the objects it contains. It also does not make a copy; the original objects can be updated and code accessing them by using the ESMF_State will see the updated version.
statename = "Ocean" state2 = ESMF_StateCreate(statename, statetype=ESMF_STATE_EXPORT, rc=rc)
bundlename = "Temperature" bundle1 = ESMF_FieldBundleCreate(name=bundlename, rc=rc) print *, "FieldBundle Create returned", rc
call ESMF_StateAdd(state2, bundle1, rc) print *, "StateAdd returned", rc
call ESMF_StateDestroy(state2, rc)
call ESMF_FieldBundleDestroy(bundle1, rc)
If a component could potentially produce a large number of optional items, one strategy is to add the names only of those objects to the ESMF_State. Other components can call framework routines to set the ESMF_NEEDED flag to indicate they require that data. The original component can query this flag and then produce only the data that is required by another component.
statename = "Ocean" state3 = ESMF_StateCreate(statename, statetype=ESMF_STATE_EXPORT, rc=rc)
dataname = "Downward wind" call ESMF_StateAdd(state3, dataname, rc)
dataname = "Humidity" call ESMF_StateAdd(state3, dataname, rc)
How to set the NEEDED state of an item.
dataname = "Downward wind" call ESMF_StateSetNeeded(state3, dataname, ESMF_NEEDED, rc)
Query an item for the NEEDED status, and creating an item on demand. Similar flags exist for "Ready", "Valid", and "Required for Restart", to mark each data item as ready, having been validated, or needed if the application is to be checkpointed and restarted. The flags are supported to help coordinate the data exchange between components.
dataname = "Downward wind" if (ESMF_StateIsNeeded(state3, dataname, rc)) then
bundlename = dataname bundle2 = ESMF_FieldBundleCreate(name=bundlename, rc=rc)
call ESMF_StateAdd(state3, bundle2, rc)
else print *, "Data not marked as needed", trim(dataname) endif
The set services routines are used to tell ESMF which routine hold the user code for the initialize, run, and finalize blocks of user level Components. These are the separate subroutines called by the code below.
! Initialize routine which creates "field1" on PETs 0 and 1 subroutine comp1_init(gcomp, istate, ostate, clock, rc) type(ESMF_GridComp) :: gcomp type(ESMF_State) :: istate, ostate type(ESMF_Clock) :: clock integer, intent(out) :: rc type(ESMF_Field) :: field1 integer :: localrc print *, "i am comp1_init" field1 = ESMF_FieldCreateEmpty(name="Comp1 Field", rc=localrc) call ESMF_StateAdd(istate, field1, rc=localrc) rc = localrc end subroutine comp1_init ! Initialize routine which creates "field2" on PETs 2 and 3 subroutine comp2_init(gcomp, istate, ostate, clock, rc) type(ESMF_GridComp) :: gcomp type(ESMF_State) :: istate, ostate type(ESMF_Clock) :: clock integer, intent(out) :: rc type(ESMF_Field) :: field2 integer :: localrc print *, "i am comp2_init" field2 = ESMF_FieldCreateEmpty(name="Comp2 Field", rc=localrc) call ESMF_StateAdd(istate, field2, rc=localrc) rc = localrc end subroutine comp2_init subroutine comp_dummy(gcomp, rc) type(ESMF_GridComp) :: gcomp integer, intent(out) :: rc rc = ESMF_SUCCESS end subroutine comp_dummy
! !PROGRAM: ESMF_StateReconcileEx - State reconciliation ! ! !DESCRIPTION: ! ! This program shows examples of using the State Reconcile function !----------------------------------------------------------------------------- ! ESMF Framework module use ESMF_Mod use ESMF_StateReconcileEx_Mod implicit none ! Local variables integer :: rc, petCount type(ESMF_State) :: state1 type(ESMF_GridComp) :: comp1, comp2 type(ESMF_VM) :: vm character(len=ESMF_MAXSTR) :: comp1name, comp2name, statename
A Component can be created which will run only on a subset of the current PET list.
! Get the global VM for this job. call ESMF_VMGetGlobal(vm=vm, rc=rc) comp1name = "Atmosphere" comp1 = ESMF_GridCompCreate(name=comp1name, petList=(/ 0, 1 /), rc=rc) print *, "GridComp Create returned, name = ", trim(comp1name) comp2name = "Ocean" comp2 = ESMF_GridCompCreate(name=comp2name, petList=(/ 2, 3 /), rc=rc) print *, "GridComp Create returned, name = ", trim(comp2name) statename = "Ocn2Atm" state1 = ESMF_StateCreate(statename, rc=rc)
Here we register the subroutines which should be called for initialization. Then we call ESMF_GridCompInitialize() on all PETs, but the code runs only on the PETs given in the petList when the Component was created.
Because this example is so short, we call the entry point code directly instead of the normal procedure of nesting it in a separate SetServices() subroutine.
! This is where the VM for each component is initialized. ! Normally you would call SetEntryPoint inside set services, ! but to make this example very short, they are called inline below. ! This is o.k. because the SetServices routine must execute from within ! the parent component VM. call ESMF_GridCompSetVM(comp1, comp_dummy, rc) call ESMF_GridCompSetVM(comp2, comp_dummy, rc) call ESMF_GridCompSetServices(comp1, comp_dummy, rc) call ESMF_GridCompSetServices(comp2, comp_dummy, rc) print *, "ready to set entry point 1" call ESMF_GridCompSetEntryPoint(comp1, ESMF_SETINIT, comp1_init, rc=rc) print *, "ready to set entry point 2" call ESMF_GridCompSetEntryPoint(comp2, ESMF_SETINIT, comp2_init, rc=rc) print *, "ready to call init for comp 1" call ESMF_GridCompInitialize(comp1, state1, rc=rc) print *, "ready to call init for comp 2" call ESMF_GridCompInitialize(comp2, state1, rc=rc)
Now we have state1 containing field1 on PETs 0 and 1, and state1 containing field2 on PETs 2 and 3. For the code to have a rational view of the data, we call ESMF_StateReconcile which determines which objects are missing from any PET, and communicates information about the object. There is the option of turning metadata reconciliation on or off with the optional parameter shown in the call below. The default behavior is for metadata reconciliation to be off. After the call to reconcile, all ESMF_State objects now have a consistent view of the data.
print *, "State before calling StateReconcile()" call ESMF_StatePrint(state1, rc=rc) call ESMF_StateReconcile(state1, vm, ESMF_ATTRECONCILE_OFF, rc=rc) print *, "State after calling StateReconcile()" call ESMF_StatePrint(state1, rc=rc)
end program ESMF_StateReconcileEx
! ESMF Framework module use ESMF_Mod implicit none ! Local variables type(ESMF_State) :: state type(ESMF_Array) :: latArray, lonArray, timeArray, humidArray, & tempArray, pArray, rhArray type(ESMF_VM) :: vm integer :: localPet, rc
------------------------------------- The following line of code will read all Array data contained in a NetCDF file, place them in ESMF_Arrays and add them to an ESMF_State. Only PET 0 reads the file; the States in the other PETs remain empty. Currently, the data is not decomposed or distributed; each PET has only 1 DE and only PET 0 contains data after reading the file. Future versions of ESMF will support data decomposition and distribution upon reading a file.
Note that the third party NetCDF library must be installed. For more details, see the "ESMF Users Guide", "Building and Installing the ESMF, Third Party Libraries, NetCDF" and the website http://www.unidata.ucar.edu/software/netcdf.
! Read the NetCDF data file into Array objects in the State on PET 0 call ESMF_StateRead(state, "io_netcdf_testdata.nc", rc=rc) ! If the NetCDF library is not present (on PET 0), cleanup and exit if (rc == ESMF_RC_LIB_NOT_PRESENT) then call ESMF_StateDestroy(state, rc=rc) goto 10 endif
Only reading data into ESMF_Arrays is supported at this time; ESMF_ArrayBundles, ESMF_Fields, and ESMF_FieldBundles will be supported in future releases of ESMF.
To see that the State now contains the same data as in the file, the following shows how to print out what Arrays are contained within the State and to print the data contained within each Array. The NetCDF utility "ncdump" can be used to view the contents of the NetCDF file. In this example, only PET 0 will contain data.
if (localPet == 0) then ! Print the names and attributes of Array objects contained in the State call ESMF_StatePrint(state, rc=rc) ! Get each Array by name from the State call ESMF_StateGet(state, "lat", latArray, rc=rc) call ESMF_StateGet(state, "lon", lonArray, rc=rc) call ESMF_StateGet(state, "time", timeArray, rc=rc) call ESMF_StateGet(state, "Q", humidArray, rc=rc) call ESMF_StateGet(state, "TEMP", tempArray, rc=rc) call ESMF_StateGet(state, "p", pArray, rc=rc) call ESMF_StateGet(state, "rh", rhArray, rc=rc) ! Print out the Array data call ESMF_ArrayPrint(latArray, rc=rc) call ESMF_ArrayPrint(lonArray, rc=rc) call ESMF_ArrayPrint(timeArray, rc=rc) call ESMF_ArrayPrint(humidArray, rc=rc) call ESMF_ArrayPrint(tempArray, rc=rc) call ESMF_ArrayPrint(pArray, rc=rc) call ESMF_ArrayPrint(rhArray, rc=rc) endif
Note that the Arrays "lat", "lon", and "time" hold spatial and temporal coordinate data for the dimensions latitude, longitude and time, respectively. These will be used in future releases of ESMF to create ESMF_Grids.
All the Array data within the State on PET 0 can be written out to a NetCDF file as follows:
! Write Arrays within the State on PET 0 to a NetCDF file call ESMF_StateWrite(state, "io_netcdf_testdata_out.nc", rc=rc)
Currently writing is limited to PET 0; future versions of ESMF will allow parallel writing, as well as parallel reading.
One important request by the user community during the ESMF object design was that there be no communication overhead or synchronization when creating distributed ESMF objects. As a consequence it is required to create these objects in unison across all PETs in order to keep the ESMF object identifiaction in sync.
When running components on subsets of the original VM all the PETs can create consistent objects but then when they are put into a State and passed to a component with a different VM and a different set of PETs, a communication call (reconcile) must be made to communicate the missing information to the PETs which were not involved in the original object creation. The reconcile call broadcasts object lists; those PETs which are missing any objects in the total list can receive enough information to reconstruct a proxy object which contains all necessary information about that object, with no local data, on that PET. These proxy objects can be queried by ESMF routines to determine the amount of data and what PETs contain data which is destined to be moved to the local PET (for receiving data) and conversely, can determine which other PETs are going to receive data and how much (for sending data).
For example, the FieldExcl system test creates 2 Gridded Components on separate subsets of PETs. They use the option of mapping particular, non-monotonic PETs to DEs. The following figures illustrate how the DEs are mapped in each of the Gridded Components in that test:
In the coupler code, all PETs must make the reconcile call before accessing data in the State. On PETs which already contain data, the objects are unchanged. On PETs which were not involved during the creation of the FieldBundles or Fields, the reconcile call adds an object to the State which contains all the same metadata associated with the object, but creates a slightly different Grid object, called a Proxy Grid. These PETs contain no local data, so the Array object is empty, and the DELayout for the Grid is like this:
The following is a simplified UML diagram showing the structure of the State class. States can contain FieldBundles, Fields, Arrays, or nested States. See Appendix A, A Brief Introduction to UML, for a translation table that lists the symbols in the diagram and their meaning.
INTERFACE:
subroutine ESMF_StateAdd(state, <item>, rc)ARGUMENTS:
type(ESMF_State), intent(inout) :: state <item>, see below for supported values integer, intent(out), optional :: rcDESCRIPTION:
Add a reference to a single <item> to an existing state. Any of the supported <item>s can be marked needed for a particular run using the ESMF_StateSetNeeded() call. The name of the <item> must be unique within the state.
One of the supported options below is to add only the name of the item to the state during a first pass. The name can be replaced with the actual <item> in a later call. When doing this, the name of the <item> provided to the state during the first pass must match the name stored in the <item> itself.
Supported values for <item> are:
The arguments are:
INTERFACE:
subroutine ESMF_StateAdd(state, <itemList>, count, rc)ARGUMENTS:
type(ESMF_State), intent(inout) :: state <itemList>, see below for supported values integer, intent(in), optional :: count integer, intent(out), optional :: rcDESCRIPTION:
Add a list of items to an ESMF_State.
Supported values for <itemList> are:
The arguments are:
INTERFACE:
function ESMF_StateCreate(stateName, statetype, & bundleList, fieldList, arrayList, nestedStateList, & nameList, itemCount, & neededflag, readyflag, validflag, reqforrestartflag, rc)RETURN VALUE:
type(ESMF_State) :: ESMF_StateCreateARGUMENTS:
character(len=*), intent(in), optional :: stateName type(ESMF_StateType), intent(in), optional :: statetype type(ESMF_FieldBundle), dimension(:), intent(inout), optional :: bundleList type(ESMF_Field), dimension(:), intent(inout), optional :: fieldList type(ESMF_Array), dimension(:), intent(in), optional :: arrayList type(ESMF_State), dimension(:), intent(in), optional :: nestedStateList character(len=*), dimension(:), intent(in), optional :: nameList integer, intent(in), optional :: itemCount type(ESMF_NeededFlag), optional :: neededflag type(ESMF_ReadyFlag), optional :: readyflag type(ESMF_ValidFlag), optional :: validflag type(ESMF_ReqForRestartFlag), optional :: reqforrestartflag integer, intent(out), optional :: rcDESCRIPTION:
Create a new ESMF_State, set default characteristics for objects added to it, and optionally add initial objects to it.
The arguments are:
INTERFACE:
recursive subroutine ESMF_StateDestroy(state, rc)ARGUMENTS:
type(ESMF_State) :: state integer, intent(out), optional :: rcDESCRIPTION:
Releases all resources associated with this ESMF_State. Actual objects added to ESMF_States will not be destroyed, it remains the user's responsibility to destroy these objects in the correct context. However, proxy objects automatically created during ESMF_StateReconcile() are destroyed when the State is destroyed.
The arguments are:
INTERFACE:
! Private name; call using ESMF_StateGet() subroutine ESMF_StateGetInfo(state, itemSearch, nestedFlag, name, statetype, itemCount, & itemNameList, stateitemtypeList, rc)ARGUMENTS:
type(ESMF_State), intent(in) :: state character (len=*), intent(in), optional :: itemSearch logical, intent(in), optional :: nestedFlag character (len=*), intent(out), optional :: name type(ESMF_StateType), intent(out), optional :: statetype integer, intent(out), optional :: itemCount character (len=*), intent(out), optional :: itemNameList(:) type(ESMF_StateItemType), intent(out), optional :: stateitemtypeList(:) integer, intent(out), optional :: rcDESCRIPTION:
Returns the requested information about this ESMF_State. The optional itemSearch argument may specify the name of an individual item to search for. When used in conjunction with the nestedFlag, nested States will also be searched.
The arguments are:
Typically, an ESMF_StateGet() information request will be performed twice. The first time, the itemCount argument will be used to query the size of arrays that are needed. Arrays can then be allocated to the correct size for itemNameList and stateitemtypeList as needed. A second call to ESMF_StateGet() will then fill in the values.
INTERFACE:
subroutine ESMF_StateGet(state, itemName, <item>, nestedStateName, rc)ARGUMENTS:
type(ESMF_State), intent(in) :: state character (len=*), intent(in) :: itemName <item>, see below for supported values character (len=*), intent(in), optional :: nestedStateName integer, intent(out), optional :: rcDESCRIPTION:
Returns an <item> from an ESMF_State by name. If the ESMF_State contains the <item> directly, only itemName is required.
If the state contains nested ESMF_States, the itemName argument may specify a fully qualified name to access the desired item with a single call. This is performed using the ``/'' character to separate the names of the intermediate State names leading to the desired item. (E.g., itemName=``state1/state12/item''.
An alternative technique for accessing a nested item which is only one level down, is to specify the nested State with the nestedStateName argument. While States can be nested to any depth, this option only searches immediate descendents. It is an error to specify a nestedStateName if the state contains no nested ESMF_States. It is an error to specify both nestedStateName and a fully qualified, nested State itemName.
Supported values for <item> are:
The arguments are:
INTERFACE:
! Private name; call using ESMF_StateGet() subroutine ESMF_StateGetItemInfo(state, name, stateitemtype, rc)ARGUMENTS:
type(ESMF_State), intent(in) :: state character (len=*), intent(in) :: name type(ESMF_StateItemType), intent(out) :: stateitemtype integer, intent(out), optional :: rcDESCRIPTION:
Returns the type for the item named name in this ESMF_State. If no item with this name exists, the value ESMF_STATEITEM_NOTFOUND will be returned and the error code will not be set to an error. Thus this routine can be used to safely query for the existance of items by name whether or not they are expected to be there. The error code will be set in case of other errors, for example if the ESMF_State itself is invalid.
The arguments are:
INTERFACE:
subroutine ESMF_StateGetNeeded(state, itemName, neededflag, rc)ARGUMENTS:
type(ESMF_State), intent(in) :: state character (len=*), intent(in) :: itemName type(ESMF_NeededFlag), intent(out) :: neededflag integer, intent(out), optional :: rcDESCRIPTION:
Returns the status of the neededflag for the data item named by itemName in the ESMF_State.
The arguments are:
INTERFACE:
function ESMF_StateIsNeeded(state, itemName, rc)RETURN VALUE:
logical :: ESMF_StateIsNeededARGUMENTS:
type(ESMF_State), intent(in) :: state character (len=*), intent(in) :: itemName integer, intent(out), optional :: rcDESCRIPTION:
Returns true if the status of the needed flag for the data item named by itemName in the ESMF_State is ESMF_STATEITEM_NEEDED. Returns false for no item found with the specified name or item marked not needed. Also sets error code if dataname not found.
The arguments are:
INTERFACE:
subroutine ESMF_StatePrint(state, options, nestedFlag, rc)ARGUMENTS:
type(ESMF_State) :: state character (len = *), intent(in), optional :: options logical, intent(in), optional :: nestedFlag integer, intent(out), optional :: rcDESCRIPTION:
Prints information about the state to stdout.
Note: Many ESMF_<class>Print methods are implemented in C++.
On some platforms/compilers there is a potential issue with interleaving
Fortran and C++ output to stdout such that it doesn't appear in
the expected order. If this occurs, the ESMF_IOUnitFlush() method
may be used on unit ESMF_IOstdout to get coherent output.
The arguments are:
INTERFACE:
subroutine ESMF_StateRead(state, fileName, fileFormat, rc)ARGUMENTS:
type(ESMF_State) :: state character (len=*), intent(in) :: fileName type (ESMF_IOFileFormat), intent(in), optional :: fileFormat integer, intent(out), optional :: rcDESCRIPTION:
Currently limited to read in all Arrays from a netCDF file and add them to a State object. Future releases will enable more items of a State to be read from a file of various formats.
Only PET 0 reads the file; the States in other PETs remain empty. Currently, the data is not decomposed or distributed; each PET has only 1 DE and only PET 0 contains data after reading the file. Future versions of ESMF will support data decomposition and distribution upon reading a file. See Section 17.3.7 for an example.
Note that the third party NetCDF library must be installed. For more details, see the "ESMF Users Guide", "Building and Installing the ESMF, Third Party Libraries, NetCDF" and the website http://www.unidata.ucar.edu/software/netcdf.
The arguments are:
INTERFACE:
subroutine ESMF_StateWrite(state, fileName, fileFormat, rc)ARGUMENTS:
type(ESMF_State) :: state character (len=*), intent(in) :: fileName type (ESMF_IOFileFormat), intent(in), optional :: fileFormat integer, intent(out), optional :: rcDESCRIPTION:
Currently limited to write out all Arrays of a State object to a netCDF file. Future releases will enable more item types of a State to be written to files of various formats.
Writing is currently limited to PET 0; future versions of ESMF will allow parallel writing, as well as parallel reading.
See Section 17.3.7 for an example.
Note that the third party NetCDF library must be installed. For more details, see the "ESMF Users Guide", "Building and Installing the ESMF, Third Party Libraries, NetCDF" and the website http://www.unidata.ucar.edu/software/netcdf.
The arguments are:
INTERFACE:
subroutine ESMF_StateReconcile(state, vm, attreconflag, rc)ARGUMENTS:
type(ESMF_State), intent(inout) :: state type(ESMF_VM), intent(in) :: vm type(ESMF_AttReconcileFlag), intent(in), optional :: attreconflag integer, intent(out), optional :: rcDESCRIPTION:
Must be called for any ESMF_State which contains ESMF objects that have not been created on all the PETs of the currently running ESMF_Component. For example, if a coupler is operating on data which was created by another component that ran on only a subset of the coupler's PETs, the coupler must make this call first before operating on any data inside that ESMF_State. After calling ESMF_StateReconcile all PETs will have a common view of all objects contained in this ESMF_State. The option to reconcile the metadata associated with the objects contained in this ESMF_State also exists. The default behavior for this capability is to not reconcile metadata unless told otherwise.
The arguments are:
NOTE: The options for attreconflag include:
INTERFACE:
subroutine ESMF_StateSetNeeded(state, itemName, neededflag, rc)ARGUMENTS:
type(ESMF_State), intent(inout) :: state character (len=*), intent(in) :: itemName type(ESMF_NeededFlag), intent(in) :: neededflag integer, intent(out), optional :: rcDESCRIPTION:
Sets the status of the needed flag for the data item named by itemName in the ESMF_State.
The arguments are:
INTERFACE:
subroutine ESMF_StateValidate(state, options, nestedFlag, rc)ARGUMENTS:
type(ESMF_State) :: state character (len = *), intent(in), optional :: options logical, intent(in), optional :: nestedFlag integer, intent(out), optional :: rcDESCRIPTION:
Validates that the state is internally consistent. Currently this method determines if the State is uninitialized or already destroyed. The method returns an error code if problems are found.
The arguments are: