INTERFACE / PUBLIC COMPONENTS / NAMELIST / FILES / REFERENCES / ERRORS / BUGS / PLANS

MODULE model_mod (pe2lyr)

Contact: Jeffrey Whitaker
Revision: $Revision: 2801 $
Source: $URL: $
Change Date: $Date: $
Change history: try "svn log" or "svn diff"

OVERVIEW

DART standard interfaces for a two-layer isentropic primitive equation model. The 17 public interfaces are standardized for all DART compliant models. These interfaces allow DART to advance the model, get the model state and metadata describing this state, find state variables that are close to a given location, and do spatial interpolation for model state variables.

This model is a 2-layer, isentropic, primitive equation model on a sphere. TODO: add more detail here, including equations, etc. Contact: Jeffrey.S.Whitaker@noaa.gov




OTHER MODULES USED

types_mod
time_manager_mod
utilities_mod
random_seq_mod
threed_sphere/location_mod



PUBLIC INTERFACE

use ensemble_manager_mod, only : get_model_size
 adv_1step
 get_state_meta_data
 model_interpolate
 get_model_time_step
 end_model
 static_init_model
 init_time
 init_conditions
 nc_write_model_atts
 nc_write_model_vars
 pert_model_state
 get_close_maxdist_init
 get_close_obs_init
 get_close_obs
 ens_mean_for_model

NOTES

Optional namelist interface &model_nml may be read from file input.nml .


PUBLIC COMPONENTS

The following routines are public.


var = get_model_size()
 integer :: get_model_size
 

Description

Returns the size of the model as an integer. For this model the default grid size is 96 (lon) by 48 (lat) by 2 levels, and 3 variables (U, V, Z) at each grid location, for a total size of 27,648. There are alternative include files which, if included at compile time instead of the default file, defines a grid at twice and 4 times this resolution. They have corresponding truncation values of T63 and T127 (the default grid uses T31).

var    The size of the model state vector.


call adv_1step(x,time)
 real(r8), dimension(:), intent(inout) :: x
 type(time_type), intent(in)           :: time
 

Description

Advances a model for a single time step. The time associated with the initial model state is also input although it is not used in this model.

x    Input state vector that is returned at the next timestep.
time    Initial time of the model state that is returned at the advanced time.


call get_state_meta_data(index_in,location,var_type)
 integer, intent(in)              :: index_in
 type(location_type), intent(out) :: location
 integer, optional, intent(out)   :: var_type
 

Description

Returns metadata about a given element, indexed by index_in, in the model state vector. The location defines where the state variable is located.

For this model, the default grid is a global lat/lon grid, 96 (lon) by 48 (lat) by 2 levels. The variable types are U, V, and Z:

Grids at twice and 4 times the resolution can be compiled in instead by using one of the alternative header files (see resolt31.h (the default), resolt63.h, and resolt127.h).

index_in    Index of state vector element about which information is requested.
location    The location of state variable element.
var_type    The type of the state variable element.


call model_interpolate(x,location,itype,obs_val,istatus)
 real(r8), dimension(:), intent(in) :: x
 type(location_type), intent(in)    :: location
 integer, intent(in)                :: itype
 real(r8), intent(out)              :: obs_val
 integer, intent(out)               :: istatus
 

Description

Given a state vector, a location, and a model state variable type, interpolates the state variable field to that location and returns the value in obs_val. The istatus variable is always returned as 0 (OK).

x    A model state vector.
location    Location to which to interpolate.
itype    Type of state field to be interpolated.
obs_val    The returned interpolated value.
istatus    Integer value returning 0 for successful, other values can be defined for various failures.


var = get_model_time_step()
 type(time_type) :: get_model_time_step
 

Description

Returns the the time step of the model; the smallest increment in time that the model is capable of advancing the state in a given implementation. For this model the default value is 20 minutes (1200 seconds), but also comes with header files with times steps of 10 and 5 minutes (for higher grid resolution and truncation constants).

var    Smallest time step of model.


call end_model()
 

Description

A stub since the pe2lyr model does no cleanup.



call static_init_model()
 

Description

Used for runtime initialization of a model, for instance calculating storage requirements, initializing model parameters, etc. This is the first call made to a model by any DART compliant assimilation routines.

In this model, it allocates space for the grid, and initializes the grid locations, data values, and various parameters, including spherical harmonic weights.



call init_time(time)
 type(time_type), intent(out) :: time
 

Description

Returns the time at which the model will start if no input initial conditions are to be used. This model sets the time to 0.

time    Initial model time.


call init_conditions(x)
 real(r8), dimension(:), intent(out) :: x
 

Description

Returns default initial conditions for model; generally used for spinning up initial model states. This model sets the default state vector based on the initialized fields in the model. (TODO: which are what?)

x    Initial conditions for state vector.


var = nc_write_model_atts(ncFileID)
 integer             :: nc_write_model_atts
 integer, intent(in) :: ncFileID
 

Description

This routine writes the model-specific attributes to a netCDF file. This includes coordinate variables and any metadata, but NOT the model state vector. This model writes out the data as U, V, and Z arrays on a lat/lon/height grid, so the attributes are organized in the same way.

The typical sequence for adding new dimensions, variables, attributes:

 NF90_OPEN             ! open existing netCDF dataset               
    NF90_redef         ! put into define mode                       
    NF90_def_dim       ! define additional dimensions (if any)     
    NF90_def_var       ! define variables: from name, type, and dims
    NF90_put_att       ! assign attribute values                    
 NF90_ENDDEF           ! end definitions: leave define mode         
    NF90_put_var       ! provide values for variable                
 NF90_CLOSE            ! close: save updated netCDF dataset        

var    Returns a 0 for successful completion.
ncFileID    netCDF file identifier.


var = nc_write_model_vars(ncFileID, statevec, copyindex, timeindex)
 integer                            :: nc_write_model_vars
 integer, intent(in)                :: ncFileID
 real(r8), dimension(:), intent(in) :: statevec
 integer, intent(in)                :: copyindex
 integer, intent(in)                :: timeindex
 

Description

This routine writes the model-specific state vector (data) to a netCDF file. This model writes out the data as U, V, and Z arrays on a lat/lon/height grid.

Typical sequence for adding new dimensions,variables,attributes:

 NF90_OPEN             ! open existing netCDF dataset               
    NF90_redef         ! put into define mode                       
    NF90_def_dim       ! define additional dimensions (if any)      
    NF90_def_var       ! define variables: from name, type, and dims
    NF90_put_att       ! assign attribute values                    
 NF90_ENDDEF           ! end definitions: leave define mode         
    NF90_put_var       ! provide values for variable                
 NF90_CLOSE            ! close: save updated netCDF dataset         

var    Returns 0 for normal completion.
ncFileID    NetCDF file identifier.
statevec    A model state vector.
copyindex    Which copy of state in output file is being written
timeindex    What time level of output is this.


call pert_model_state(state, pert_state, interf_provided)
 real(r8), dimension(:), intent(in)  :: state
 real(r8), dimension(:), intent(out) :: pert_state
 logical, intent(out)                :: interf_provided
 

Description

Given a model state vector, perturbs this vector. Used to generate initial conditions for spinning up ensembles. This model has no code to generate these values, so it returns interf_provided as .false. and the default algorithms in filter are then used by the calling code.

state    State vector to be perturbed.
pert_state    Perturbed state vector
interf_provided    Returned .false.


call ens_mean_for_model(ens_mean)
 real(r8), dimension(:), intent(in)  :: ens_mean
 

Description

Not needed by this model.

ens_mean    Ensemble mean state vector


call get_close_maxdist_init(gc, maxdist)
 type(get_close_type), intent(inout) :: gc
 real(r8),             intent(in)    :: maxdist
 

Description

In distance computations any two locations closer than the given maxdist will be considered close by the get_close_obs() routine. Pass-through to the 3-D sphere locations module. See get_close_maxdist_init() for the documentation of this subroutine.

gc    The get_close_type which stores precomputed information about the locations to speed up searching
maxdist    Anything closer than this will be considered close.


call get_close_obs_init(gc, num, obs)
 type(get_close_type), intent(inout) :: gc
 integer,              intent(in)    :: num
 type(location_type),  intent(in)    :: obs(num)
 

Description

Pass-through to the 3-D sphere locations module. See get_close_obs_init() for the documentation of this subroutine.

gc    The get_close_type which stores precomputed information about the locations to speed up searching
num    The number of items in the third argument
obs    A list of locations which will be part of the subsequent distance computations


call get_close_obs(gc, base_obs_loc, base_obs_kind, obs, obs_kind, num_close, close_ind [, dist])
 type(get_close_type), intent(in)  :: gc
 type(location_type),  intent(in)  :: base_obs_loc
 integer,              intent(in)  :: base_obs_kind
 type(location_type),  intent(in)  :: obs(:)
 integer,              intent(in)  :: obs_kind(:)
 integer,              intent(out) :: num_close
 integer,              intent(out) :: close_ind(:)
 real(r8), optional,   intent(out) :: dist(:)
 

Description

Given a location and kind, compute the distances to all other locations in the obs list. The return values are the number of items which are within maxdist of the base, the index numbers in the original obs list, and optionally the distances. The gc contains precomputed information to speed the computations.

Pass-through to the 3-D sphere locations module. See get_close_obs() for the documentation of this subroutine.

gc    The get_close_type which stores precomputed information about the locations to speed up searching
base_obs_loc    Reference location. The distances will be computed between this location and every other location in the obs list
base_obs_kind  The kind of base_obs_loc
obs    Compute the distance between the base_obs_loc and each of the locations in this list
obs_kind    The corresponding kind of each item in the obs list
num_close    The number of items from the obs list which are within maxdist of the base location
close_ind    The list of index numbers from the obs list which are within maxdist of the base location
dist    If present, return the distance between each entry in the close_ind list and the base location. If not present, all items in the obs list which are closer than maxdist will be added to the list but the overhead of computing the exact distances will be skipped.


NAMELIST

We adhere to the F90 standard of starting a namelist with an ampersand '&338;' and terminating with a slash '/'.

 namelist / model_nml / 
                                                                                
                                                               
 

Discussion

This model currently has no values settable by namelist.

The namelist would be read in a file called input.nml





FILES




REFERENCES

Zou, X., Barcilon, A., Navon, I.M., Whitaker, J., Cacuci, D.G.. 1993: An Adjoint Sensitivity Study of Blocking in a Two-Layer Isentropic Model. Monthly Weather Review: Vol. 121, No. 10, pp. 2833-2857.


ERROR CODES and CONDITIONS





KNOWN BUGS




FUTURE PLANS