pe2lyr
Attention
pe2lyr
works with versions of DART before Manhattan (9.x.x) and has yet to be updated. If you are interested in
using pe2lyr
with more recent versions of DART, contact DAReS staff to assess the feasibility of an update.
Until that time, you should consider this documentation as out-of-date.
Overview
DART standard interfaces for a two-layer isentropic primitive equation model.
The 16 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 interfaces
use model_mod, only : |
get_model_size |
adv_1step |
|
get_state_meta_data |
|
model_interpolate |
|
get_model_time_step |
|
static_init_model |
|
end_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 |
A note about documentation style. Optional arguments are enclosed in brackets [like this].
model_size = get_model_size( )
integer :: get_model_size
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).
|
The length of the model state vector. |
call adv_1step(x, time)
real(r8), dimension(:), intent(inout) :: x
type(time_type), intent(in) :: time
Advances the model for a single time step. The time associated with the initial model state is also input although it is not used for the computation.
|
State vector of length model_size. |
|
Specifies time of the initial model state. |
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
index_in
, in the model state vector. The location
defines where the state variable is located.1 = TYPE_u
2 = TYPE_v
901 = TYPE_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 of state vector element about which information is requested. |
|
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
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).
|
A model state vector. |
|
Location to which to interpolate. |
|
Type of state field to be interpolated. |
|
The interpolated value from the model. |
|
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
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).
|
Smallest time step of model. |
call static_init_model()
call end_model()
A stub since the pe2lyr model does no cleanup.
call init_time(time)
type(time_type), intent(out) :: time
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.
|
Initial model time. |
call init_conditions(x)
real(r8), dimension(:), intent(out) :: x
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?)
|
Initial conditions for state vector. |
ierr = nc_write_model_atts(ncFileID)
integer :: nc_write_model_atts
integer, intent(in) :: ncFileID
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.
|
Integer file descriptor to previously-opened netCDF file. |
|
Returns a 0 for successful completion. |
ierr = 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
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.
|
file descriptor to previously-opened netCDF file. |
|
A model state vector. |
|
Integer index of copy to be written. |
|
The timestep counter for the given state. |
|
Returns 0 for normal completion. |
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
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 vector to be perturbed. |
|
Perturbed state vector |
|
Returned false; interface is not implemented. |
call get_close_maxdist_init(gc, maxdist)
type(get_close_type), intent(inout) :: gc
real(r8), intent(in) :: maxdist
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.
|
The get_close_type which stores precomputed information about the locations to speed up searching |
|
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)
Pass-through to the 3-D sphere locations module. See get_close_obs_init() for the documentation of this subroutine.
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(:)
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.call ens_mean_for_model(ens_mean)
real(r8), dimension(:), intent(in) :: ens_mean
Stub only. Not needed by this model.
|
State vector containing the ensemble mean. |
This model currently has no values settable by namelist.
Files
The model source is in pe2lyr_mod.f90, and the spherical harmonic code is in spharmt_mod.f90. The various resolution settings are in resolt31.h, resolt63.h, and resolt127.h.
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.
Private components
N/A