MODULE model_mod (Lorenz_04)

DART project logo

Jump to DART Documentation Main Index
version information for this file:
$Id: model_mod.html 6382 2013-08-07 20:45:16Z nancy $

NAMELIST / INTERFACES / FILES / REFERENCES / ERRORS / PLANS / TERMS OF USE

Overview

DART interface module for the Lorenz-05 models II and III. 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.

The reference for these models is Lorenz, E.N., 2005: Designing chaotic models. J. Atmos. Sci., 62, 1574-1587.

Model II is a single-scale model, similar to Lorenz 96, but with spatial continuity in the waves. Model III is a two-scale model. It is fudamentally different from the Lorenz 96 two-scale model because of the spatial continuity and the fact that both scales are projected onto a single variable of integration. The scale separation is achived by a spatial filter and is therefore not perfect (i.e. there is leakage). The slow scale in model III is model II, and thus model II is a deficient form of model III. The basic equations are documented in Lorenz (2005) and also in the model_mod.f90 code. The user is free to choose model II or III with a Namelist variable.



[top]

NAMELIST

This namelist is read from the file input.nml. Namelists start with an ampersand '&' and terminate with a slash '/'. Character strings that contain a '/' must be enclosed in quotes to prevent them from prematurely terminating the namelist.

&model_nml
   model_size        = 960,
   forcing           = 15.00,
   delta_t           = 0.001,
   space_time_scale  = 10.00,
   coupling          = 3.00,
   K                 = 32,
   smooth_steps      = 12,
   time_step_days    = 0,
   time_step_seconds = 3600,
   model_number      = 3  
/

Contents Type Description
model_size integer Number of variables in model
forcing real(r8) Forcing, F, for model
delta_t real(r8) Non-dimensional timestep
space_time_scale real(r8) Determines temporal and spatial relationship between fast and slow variables (model III)
coupling real(r8) Linear coupling between fast and slow variables (model III)
K integer Determines the wavenumber of the slow variables (K=1, smooth_steps=0 reduces model II to Lorenz 96)
smooth_steps integer Determines filter length to separate fast and slow scales
time_step_days integer Arbitrary real time step days
time_step_seconds integer Arbitrary real time step seconds (could choose this for proper scaling)
model_number integer 2 = single-scale, 3 = 2-scale. (This follows the notation in the paper.)


[top]

OTHER MODULES USED

types_mod
time_manager_mod
oned/location_mod
utilities_mod


[top]

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
 nc_read_model_vars
 pert_model_state
 get_close_maxdist_init
 get_close_obs_init
 get_close_obs
 ens_mean_for_model

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

A note about documentation style. Optional arguments are enclosed in brackets [like this].


model_size = get_model_size( )
integer :: get_model_size

Returns the length of the model state vector as an integer. Settable by namelist, the default is 960.

model_size 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 a model for a single time step using a four-step Runga-Kutta. The time argument is passed in but not used by this model.

x State vector of length model_size.
time    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 

Returns metadata about a given element, indexed by index_in, in the model state vector. The location defines where the state variable is located (at present, a variety of simple location models for support of gridpoint models are provided) while the type of the variable (for instance temperature, or u wind component) is returned by var_type. In this model there are no types so it always returns 1.

index_in    Index of state vector element about which information is requested.
location The location of state variable element.
var_type Returns the type (always 1) of the indexed state variable as an optional argument.


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 model state, returns the value of variable itype interpolated to a given location.

x A model state vector.
location    Location to which to interpolate.
itype Integer indexing which type of state variable is to be interpolated. Is ignored for this model since there is only a single type of variable.
obs_val The interpolated value from the model.
istatus Quality control information about the observation of the model state. 0 == successful.


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

Returns the time step (forecast length) of the model; The time step defaults to 1 hour but is settable by namelist.

var    Smallest time step of model.


call static_init_model()

Reads the namelist, allocates space for the state vector, initializes the locations, generates some initial parameters.



call end_model()

A stub routine in this model.



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 is used to spin-up the model from rest. Sets time to 0 in this model.

time    Initial model time.


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

Sets the entire state variable to the value of the forcing and then slightly perturbs the first element.

x    Initial conditions for state vector.


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

Function to write model specific attributes to a netCDF file. At present, DART is using the NetCDF format to output diagnostic information. This is not a requirement, and models could choose to provide output in other formats. This function writes the metadata associated with the model to a NetCDF file opened to a file identified by ncFileID.

ncFileID    Integer file descriptor to previously-opened netCDF file.
ierr 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

Writes a copy of the state variables to a netCDF file. Multiple copies of the state for a given time are supported, allowing, for instance, a single file to include multiple ensemble estimates of the state.

ncFileID file descriptor to previously-opened netCDF file.
statevec A model state vector.
copyindex    Integer index of copy to be written.
timeindex The timestep counter for the given state.
ierr Returns 0 for normal completion.


ierr = nc_read_model_vars(ncFileID, statevec, copyindex, timeindex)
integer                               :: ierr
integer,                  intent(in)  :: ncFileID
real(r8), dimension(:),   intent(out) :: statevec
integer,                  intent(in)  :: copyindex
integer,                  intent(in)  :: timeindex

A general interface routine to read the state vector at a particular time for a particular ensemble member.

ierr Returned error code. 0 == normal termination.
ncFileID Integer file descriptor connected to NetCDF file.
statevec State vector.
copyindex    Integer index to designate copy to be read.
timeindex Integer index of desired time to be read.

This interface is essentially untested in the Lorenz 04/05 model. It began life as a stub, but there is so little to it ... TJH 05 Aug 2005.



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, produces a perturbed model state. This is used to generate initial ensemble conditions perturbed around some control trajectory state when one is preparing to spin-up ensembles. A DART compliant model can choose not to provide an implementation of this algorithm and use the default mechanism in DART by simply returning .false. as a returned value for the interf_provided argument. In this case, DART perturbs the state to generate ensemble members by adding a random sample from a N(0.0, 0.002) distribution independently to each state variable. This model returns .false. for interf_provided and uses the default DART perturbation routine.

state State vector to be perturbed
pert_state Perturbed state vector is returned
interf_provided    Returns false to have DART perturb state

The Lorenz 04/05 model relies implicitly on the default DART perturbation mechanism, i.e. this routine always returns .false..



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

Pass-through to the 1-D locations module. See get_close_maxdist_init() for the documentation of this subroutine.



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 1-D 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(:)

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



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

A NULL INTERFACE in this model.

ens_mean    State vector containing the ensemble mean.


[top]

FILES

filename purpose
input.nml to read the model_mod namelist
Prior_Diag.nc the time-history of the model state before assimilation
Posterior_Diag.nc  the time-history of the model state after assimilation
dart_log.out [default name] the run-time diagnostic output
dart_log.nml [default name] the record of all the namelists actually USED - contains the default values


[top]

REFERENCES

  1. Lorenz, E.N., 2005: Designing chaotic models. J. Atmos. Sci., 62, 1574-1587.


[top]

ERROR CODES and CONDITIONS

Routine Message Comment
nc_write_model_atts
nc_write_model_vars
Various netCDF-f90 interface error messages From one of the netCDF calls in the named routine

KNOWN BUGS

none at this time



[top]

FUTURE PLANS

none at this time



[top]

Terms of Use

DART software - Copyright 2004 - 2013 UCAR.
This open source software is provided by UCAR, "as is",
without charge, subject to all terms of use at
http://www.image.ucar.edu/DAReS/DART/DART_download

Contact: DART core group
Revision: $Revision: 6382 $
Source: $URL: https://svn-dares-dart.cgd.ucar.edu/DART/releases/classic/models/lorenz_04/model_mod.html $
Change Date: $Date: 2013-08-07 14:45:16 -0600 (Wed, 07 Aug 2013) $
Change history: try "svn log" or "svn diff"