MODULE model_mod

DART project logo

Jump to DART Documentation Main Index
version information for this file:
$Id: model_mod.html 11400 2017-03-24 22:38:14Z thoar@ucar.edu $

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

Overview

DART interface module for the dynamical core of the GFDL AM2 Bgrid model. This model is subroutine callable from DART and can be run in a similar fashion to low-order models that produce diagnostic output files with multiple assimilation times per file.

The Bgrid model was originally configured as a comprehensive atmospheric model (Anderson et al. 2004). All of that code remains in the directories under the models/bgrid_solo directory, however, much of the capability has been disabled by code modification. What is left is a dry dynamical core for a model with no dirunal cycle at equinox with Held-Suarez forcing (Held and Suarez 1994). The default settings are for a model with a 60x30 horizontal grid and 5 vertical levels. This is close to the smallest version that has somewhat realistic baroclinic instability resulting in mid-latitude 'storm tracks'. The model resolution can be changed with the entries in the bgrid_cold_start_nml. It may be necessary to change the model time step to maintain stability for larger model grids. The model state variables are the gridded surface pressure, temperature, and u and v wind components.

Some examples of ways in which this model can be configured and modified to test DART assimilation capabilities are documented in Anderson et al. (2005).

Several programs that generate interesting observation sequences are available in the models/bgrid_solo directory. These programs take interactive user input and create a text file that can be piped into program create_obs_sequence to create obs_sequence files. These can serve as examples for users who are interested in designing their own custom obs_sequence files.
Program column_rand creates an obs_sequence with randomly located columns of observations (essentially synthetic radiosondes) that observe surface pressure along with temperature and wind components at all model levels.
Program id_set_def_stdin generates an obs_sequence file that observes every state variable with error variance of 10000 for surface pressure and 1.0 for temperature and wind components.
Program ps_id_stdin generates an obs_sequence that observes every surface pressure variable for the default model size (30x60) with an error variance of 100.
Program ps_rand_local generates a set of randomly located surface pressure observations with an interactively specified error variance. It also allows the observations to be confined to a rectangular subdomain.

[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 
    current_time =  0, 0, 0, 0
    override = .false.,
    dt_atmos = 3600,
    days     = 10,
    hours    = 0,
    minutes  = 0,
    seconds  = 0,
    noise_sd = 0.0,
    dt_bias  = -1,
    state_variables = 'ps', 'QTY_SURFACE_PRESSURE',
                      't',  'QTY_TEMPERATURE',
                      'u',  'QTY_U_WIND_COMPONENT',
                      'v',  'QTY_V_WIND_COMPONENT',
    template_file = 'perfect_input.nc'
  /

# only used if initial conditions file not specified in run
&bgrid_cold_start_nml
    nlon = 60,
    nlat = 30,
    nlev = 5,
    equal_vert_spacing = .true.
  /

#Values set in hs_forcing_nml are from:
#Held, I. M., and M. J. Suarez, 1994: A proposal for the intercomparison of the dynamical cores of atmospheric general circulation models.
#Bulletin of the American Meteorological Society, 75(10), 1825-1830.

&hs_forcing_nml
    delh      =  60.,
    t_zero    = 315.,
    t_strat   = 200.,
    delv      =  10.,
    eps       =   0.,
    ka        = -40.,
    ks        =  -4.,
    kf        =  -1.,
    sigma_b   =  .7,
    do_conserve_energy = .false.
  /

&bgrid_core_driver_nml
    damp_coeff_wind   = 0.10,
    damp_coeff_temp   = 0.10,
    damp_coeff_tracer = 0.10,
    advec_order_wind   = 4,
    advec_order_temp   = 2,
    advec_order_tracer = 2,
    num_sponge_levels = 1,
    sponge_coeff_wind   = 1.00,
    sponge_coeff_temp   = 1.00,
    sponge_coeff_tracer = 1.00,
    num_fill_pass = 2,
    decomp = 0,0,
    num_adjust_dt = 3,
    num_advec_dt  = 3,
    halo = 1,
    do_conserve_energy = .false.
  /

&bgrid_integrals_nml
    file_name  = 'dynam_integral.out',
    time_units = 'days',
    output_interval = 1.00
  /

Item Type Description
current_time(4) integer Specifies the initial time of the Bgrid model internal clock. The four integer values are the day, hour, minutes, and seconds. The default version of the Bgrid model has neither a diurnal or seasonal cycle, so these can all be set to 0, the default value.
override logical If true, then the initial model date is taken from namelist entry current_time, even if an atmos_model.res file is found in directory INPUT. For most DART applications, atmospheric restart values are coming from DART files and no INPUT directory is used.
dt_atmos integer Model timestep in seconds.
noise_sd real(r8) Standard deviation of random perturbations to the time tendency of temperature applied at each timestep. Each gridpoint value of the computed temperature tendency is multiplied by 1+N(0, noise_sd) before the updated values of temperature are computed.
dt_bias integer Allows a simple mechanism to simulate model error. If dt_bias is non-zero, the assimilation programs believe that each model advance changes the time by dt_bias. However, internally the bgrid model is moving things forward by dt_atmos. By running perfect_model_obs with one time step for the internal bgrid clock (for instance dt_atmos = 3600, dt_bias = 3600), and filter with another (dt_atmos = 3000, and dt_bias = 3600) model error is simulated.
state_variables(:,2) character(len=129) Strings that identify the bgrid_solo variables that should be part of the DART state vector. The first column is the netCDF variable name, the second column is the corresponding DART quantity.
template_file character(len=256) This is the name of the file that specifies the resolution of the variables DART uses to create the DART state vector. If template_file = "null" the &bgrid_cold_start_nml namelist variables are used to specify the resolution. The actual input filenames for filter and perfect_model_obs come from their respective namelists. The resolutions in the file specified in template_file must match the resolutions of the variables in the input filenames. To start an experiment with a new model resolution, set template_file to "null" and set the resolutions in bgrid_cold_start_nml.

The following values are set in the bgrid_cold_start_nml
nlon integer The number of longitudes on the model grid.
nlat integer The number of latitudes on the model grid.
nlev integer The number of model levels.
equal_vertical_spacing logical Model levels are equally spaced in pressure if true.

The Held-Suarez forcing details can be modified with the hs_forcing_nml namelist using the documentation in Held and Suarez (1994).

Model dynamics can be adjusted with the bgrid_core_driver_nml following the documentation in the references and internal documentation in the bgrid code.

[top]

OTHER MODULES USED

bgrid_core_driver_mod
bgrid_prog_var_mod
bgrid_horiz_mod
bgrid_change_grid_mod
bgrid_vert_mod
bgrid_halo_mod
hs_forcing_mod
location_mod
random_seq_mod
types_mod
utilities_mod
obs_kind_mod
mpi_utilities_mod
ensemble_manager_mod
distributed_state_mod
state_structure_mod
dart_time_io_mod
time_manager_mod
fms_mod

Know that these modules have their own dependencies

[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
 init_time
 init_conditions
 nc_write_model_atts
 nc_write_model_vars
 pert_model_copies
 get_close_maxdist_init
 get_close_obs_init
 get_close_obs
 get_close_state_init
 get_close_state
 query_vert_localization_coord
 vert_convert
 read_model_time
 write_model_time
 end_model

A namelist interface &model_nml may be defined by the module, in which case it will be read from file input.nml. The details of the namelist are always model-specific (there are no generic namelist values).

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. Required.

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

Does a single timestep advance of the model. The input value of the vector x is the starting condition and x must be updated to reflect the changed state after a timestep. The time argument is intent in and is used for models that need to know the date/time to compute a timestep, for instance for radiation computations. This interface is only called if the namelist parameter async is set to 0 in perfect_model_obs or filter or if the program integrate_model is to be used to advance the model state as a separate executable. If one of these options is not going to be used (the model will only be advanced as a separate model-specific executable), this can be a NULL INTERFACE. (The subroutine name must still exist, but it can contain no code and it will not be called.)

x State vector of length model_size.
time    Specifies time of the initial model state.


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

Given a handle to a state structure and an integer index into the state vector, returns the associated location. A second intent(out) optional argument returns the generic quantity of this item, e.g. QTY_TEMPERATURE, QTY_DENSITY, QTY_SALINITY, QTY_U_WIND_COMPONENT. This interface is required to be functional for all applications.

state_handle    The handle to the state structure containing the state vector about which information is requested.
index_in Index of state vector element about which information is requested.
location The location of state variable element.
var_type The generic kind of the state variable element.


call model_interpolate(state_handle, x, location, itype, obs_val, istatus)
type(ensemble_type),    intent(in)  :: state_handle
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 handle containing information for a state vector, a state vector, a location, and a model state variable kind interpolates the state variable field to that location and returns the value in obs_val. The istatus variable should be returned as 0 unless there is some problem in computing the interpolation in which case a positive value should be returned. The itype variable is one of the KIND parameters defined in the obs_kind_mod.f90 file and defines which generic kind of item is being interpolated. In low-order models that have no notion of kinds of variables this argument may be ignored. For applications in which only perfect model experiments with identity observations (i.e. only the value of a particular state variable is observed), this can be a NULL INTERFACE. Otherwise it is required (which is the most common case).

state_handle    The handle to the state structure containing information about the state vector about which information is requested.
x A model state vector.
location    Location to which to interpolate.
itype Quantity of state field to be interpolated.
obs_val The interpolated value from the model.
istatus Integer value returning 0 for success. Other values can be defined for various failures.


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

Returns the time step (forecast length) of the model; the smallest increment in time that the model is capable of advancing the state in a given implementation. The actual value may be set by the model_mod namelist (depends on the model). This interface is required for all applications.

var    Smallest time step of model.


call static_init_model()

Called to do one time initialization of the model. As examples, might define information about the model size or model timestep. In models that require pre-computed static data, for instance spherical harmonic weights, these would also be computed here. Can be a NULL INTERFACE for the simplest models.



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

Companion interface to init_conditions. Returns a time that is somehow appropriate for starting up a long integration of the model. At present, this is only used if the perfect_model_obs namelist parameter read_input_state_from_file = .false. If this option should not be used in perfect_model_obs, calling this routine should issue a fatal error.

time    Initial model time.


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

Returns a model state vector, x, that is some sort of appropriate initial condition for starting up a long integration of the model. At present, this is only used if the perfect_model_obs namelist parameter read_input_state_from_file = .false. If this option should not be used in perfect_model_obs, calling this routine should issue a fatal error.

x    Initial conditions for state vector.


ierr = nc_write_model_atts(ncFileID)
integer, intent(in)  :: ncFileID
logical, intent(out) :: model_mod_writes_state_variables
integer              :: nc_write_model_atts

This routine writes the model-specific attributes to netCDF files that DART creates. This includes coordinate variables and any metadata, but NOT the actual model state vector. models/template/model_mod.f90 contains code that can be used for any model as-is.

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, kind, 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        
ncFileID Integer file descriptor to previously-opened netCDF file.
model_mod_writes_state_variables    logical descriptor to flag whether or not the model_mod will write the actual state variables or the DART-intrinsic routines will write them.
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

This routine may be used to write the model-specific state vector (data) to a netCDF file. Only used if model_mod_writes_state_variables = .true.

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, kind, 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         
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.


call pert_model_copies(state_ens_handle, ens_size, pert_amp, interf_provided)
type(ensemble_type), intent(inout) :: state_ens_handle
integer,             intent(in)    :: ens_size
real(r8),            intent(in)    :: pert_amp
logical,             intent(out)   :: interf_provided

Given an ensemble handle, the ensemble size, and a perturbation amplitude; perturb the ensemble. Used to generate initial conditions for spinning up ensembles. If the model_mod does not want to do this, instead allowing the default algorithms in filter to take effect, interf_provided =&nbps;.false. and the routine can be trivial. Otherwise, interf_provided must be returned as .true.

state_ens_handle handle containing and ensemble of state vectors to be perturbed.
ens_size the number of ensemble members to perturb.
pert_amp the amplitude of the perturbations. The interpretation is based on the model-specific implementation.
interf_provided    Returns false if model_mod cannot do this, else true.


call get_close_maxdist_init(gc, maxdist [, maxdist_list])
type(get_close_type), intent(inout) :: gc
real(r8),             intent(in)    :: maxdist
real(r8), optional,   intent(in)    :: maxdist_list(:)

In distance computations any two locations closer than the given maxdist will be considered close by the get_close_obs() routine. In general this is a PASS-THROUGH ROUTINE. It is listed on the use line for the locations_mod, and in the public list for this module, but has no subroutine declaration and no other code in this module:

use location_mod, only: get_close_maxdist_init

public :: get_close_maxdist_init

The location module has code which stores maxdist in the gc derived type. However, if the model needs to alter the value or wants to supply an alternative implementation it can intercept the call like so:

use location_mod, only: &
        lm_get_close_maxdist_init => get_close_maxdist_init

public :: get_close_maxdist_init

In this case a local get_close_maxdist_init() routine must be supplied. To call the original code in the location module use:

call lm_get_close_maxdist_init(gc, mymaxdist)

This subroutine will be called before get_close_obs_init and get_close_obs.

In most cases the PASS-THROUGH ROUTINE will be used.

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.
maxdist_list(:)    An array of distances - one for each specific type. This allows different types to have different localization behavior.


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)

In general this is a PASS-THROUGH ROUTINE. The default routine in the location module precomputes information to accelerate the distance computations done by get_close_obs(). Like the other PASS-THROUGH ROUTINES it is listed on the use line for the locations_mod, and in the public list for this module, but has no subroutine declaration and no other code in this module:

use location_mod, only: get_close_obs_init

public :: get_close_obs_init

The location module code bins the list of locations and precomputes maximum possible distances between bins. However, if the model needs to alter the values or wants to supply an alternative implementation it can intercept the call like so:

use location_mod, only: &
        lm_get_close_obs_init => get_close_obs_init
        
public :: get_close_obs_init

In this case a local get_close_obs_init() routine must be supplied. To call the original code in the location module use:

call lm_get_close_obs_init(gc, num, obs)

This subroutine will be called after get_close_maxdist_init and before get_close_obs.

In most cases the PASS-THROUGH ROUTINE will be used.

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_loc, obs_kind, num_close, close_ind, dist, state_handle)
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_loc(:)
integer,              intent(in)  :: obs_kind(:)
integer,              intent(out) :: num_close
integer,              intent(out) :: close_ind(:)
real(r8),             intent(out) :: dist(:)
type(ensemble_type),  intent(in)  :: state_handle

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.

In general this is a PASS-THROUGH ROUTINE. It is listed on the use line for the locations_mod, and in the public list for this module, but has no subroutine declaration and no other code in this module:

use location_mod, only: get_close_obs

public :: get_close_obs

However, if the model needs to alter the values or wants to supply an alternative implementation it can intercept the call like so:

use location_mod, only: &
        lm_get_close_obs => get_close_obs
        
public :: get_close_obs

In this case a local get_close_obs() routine must be supplied. To call the original code in the location module use:

call lm_get_close_obs(gc, base_obs_loc, ...)

This subroutine will be called after get_close_maxdist_init and get_close_obs_init.

In most cases the PASS-THROUGH ROUTINE will be used, but some models need to alter the actual distances depending on the observation or state vector kind, or based on the observation or state vector location. It is reasonable in this case to leave get_close_maxdist_init() and get_close_obs_init() as pass-through routines and intercept only get_close_obs(). The local get_close_obs() can first call the location mod routine and let it return a list of values, and then inspect the list and alter or remove any entries as needed. See the CAM and WRF model_mod files for examples of this use.

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_loc(:) 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.
state_handle    The handle to the state structure containing information about the state vector about which information is requested.


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

In general this is a PASS-THROUGH ROUTINE. The default routine in the location module precomputes information to accelerate the distance computations done by get_close_state(). Like the other PASS-THROUGH ROUTINES it is listed on the use line for the locations_mod, and in the public list for this module, but has no subroutine declaration and no other code in this module:

use location_mod, only: get_close_state_init

public :: get_close_state_init

The location module code bins the list of locations and precomputes maximum possible distances between bins. However, if the model needs to alter the values or wants to supply an alternative implementation it can intercept the call like so:

use location_mod, only: &
        loc_get_close_state_init => get_close_state_init
        
public :: get_close_state_init

In this case a local get_close_state_init() routine must be supplied. To call the original code in the location module use:

call loc_get_close_state_init(gc, num, obs)

This subroutine will be called after get_close_maxdist_init and before get_close_state.

In most cases the PASS-THROUGH ROUTINE will be used.

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_state(gc, base_obs_loc, base_obs_kind, state_loc, state_kind, num_close, close_ind, dist, state_handle)
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)  :: state_loc(:)
integer,              intent(in)  :: state_kind(:)
integer,              intent(out) :: num_close
integer,              intent(out) :: close_ind(:)
real(r8),             intent(out) :: dist(:)
type(ensemble_type),  intent(in)  :: state_handle

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

In general this is a PASS-THROUGH ROUTINE. It is listed on the use line for the locations_mod, and in the public list for this module, but has no subroutine declaration and no other code in this module:

use location_mod, only: get_close_state

public :: get_close_state

However, if the model needs to alter the values or wants to supply an alternative implementation it can intercept the call like so:

use location_mod, only: &
        lm_get_close_state => get_close_state
        
public :: get_close_state

In this case a local get_close_state() routine must be supplied. To call the original code in the location module use:

call loc_get_close_state(gc, base_obs_loc, ...)

This subroutine will be called after get_close_maxdist_init and get_close_state_init.

In most cases the PASS-THROUGH ROUTINE will be used, but some models need to alter the actual distances depending on the observation or state vector kind, or based on the observation or state vector location. It is reasonable in this case to leave get_close_maxdist_init() and get_close_state_init() as pass-through routines and intercept only get_close_state(). The local get_close_state() can first call the location mod routine and let it return a list of values, and then inspect the list and alter or remove any entries as needed. See the CAM and WRF model_mod files for examples of this use.

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
state_loc(:) Compute the distance between the base_obs_loc and each of the locations in this list
state_kind(:) The corresponding kind of each item in the state_loc list
num_close The number of items from the state_loc list which are within maxdist of the base location
close_ind(:) The list of index numbers from the state_loc 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 state_loc list which are closer than maxdist will be added to the list but the overhead of computing the exact distances will be skipped.
state_handle    The handle to the state structure containing information about the state vector about which information is requested.


ivert = query_vert_localization_coord()

Returns the integer code for the vertical coordinate system used for localization.

ivert The integer code for the vertical coordinate system.


call vert_convert(state_handle, location, obs_kind, istatus)
type(ensemble_type), intent(in)  :: state_handle
type(location_type), intent(in)  :: location
integer,             intent(in)  :: obs_kind
integer,             intent(out) :: istatus

Converts the state to the desired vertical localization coordinate system. Some models (toy models with no 'real' observations) will not need this. Most (real) models have observations in one or more coordinate systems (pressure, height) and the model is generally represented in only one coordinate system. To be able to interpolate the model state to the observation location, or to compute the true distance between the state and the observation, it is necessary to convert everything to one coodinate system.

state_handle    The handle to that state.
location the desired location
obs_kind the quantity defining the part of state to convert.
istatus Specifies the success or failure of the vertical conversion. If istatus = 0, the conversion was a sucess. Any other value is a failure.


model_time = read_model_time(filename)
character(len=*), intent(in) :: filename
type(time_type)              :: model_time

Reads the valid time of the model state in a netCDF file. There is a default routine in assimilation_code/modules/io/dart_time_io_mod.f90 that can be used as a pass-through. That routine will read the last timestep of a 'time' variable - which is the same strategy used for reading netCDF files that have multiple timesteps in them. If your model has some other representation of time (i.e. it does not use a netCDF variable named 'time') - you will have to write this routine.

ncid handle to an open netCDF file
dart_time    Specifies the (last) time of the model state.


call write_model_time(ncid, dart_time)
integer,          intent(in) :: ncid
type(time_type),  intent(in) :: dart_time

Writes the assimilation time to a netCDF file. There is a default routine in assimilation_code/modules/io/dart_time_io_mod.f90 that can be used as a pass-through. If your model has some other representation of time (i.e. it does not use a netCDF variable named 'time') - you will have to write this routine.

ncid handle to an open netCDF file
dart_time    Specifies the time of the assimilation (the current time step).


call end_model()

Does any shutdown and clean-up needed for model. Can be a NULL INTERFACE if the model has no need to clean up storage, etc.


[top]

FILES

Files used or created by perfect_model_obs
 
filename purpose
input.nml to read the namelist specifying the run-time control
perfect_input.nc the initial model state (if not cold-start)
obs_seq.in the metadata for the desired synthetic observations
perfect_output.nc the time-history of the model state - the True State
obs_seq.out the synthetic observations created from the True State
dart_log.[out,nml] run-time output
 
Files used or created by filter
 
filename purpose
filter_input_list.txt The file specifying the filename(s) of the initial ensemble. The default filename for the initial ensemble is filter_input.nc
obs_seq.out the observations to be assimilated
preassim.nc the ensemble information after prior inflation but before assimilation
postassim.nc the ensemble information after assimilation but before posterior inflation
filter_output_list.txt    the file specifying the filename(s) of the time-evolution of the ensemble after all inflations and assimilation - the 'posterior'. The default filename for the output ensemble is filter_output.nc
obs_seq.final the original observations and (optionally) the expected observations from the ensemble of model states. This may contain both the prior and posterior expected observation values.
dart_log.[out,nml] run-time output

[top]

REFERENCES

  1. Anderson, J. L. and Coauthors, 2004: The new GFDL global atmosphere and land model AM2-LM2: Evaluation with prescribed SST simulations. J. Climate, 17, 4641-4673.
  2. Anderson, J. L., Wyman, B., Zhang, S. & Hoar, T., 2005: Assimilation of surface pressure observations using an ensemble filter in an idealized global atmospheric prediction system, J. Atmos. Sci., 62, 2925-2938.
  3. Held, I. M., and M. J. Suarez, 1994: A proposal for the intercomparison of the dynamical cores of atmospheric general circulation models. Bulletin of the American Meteorological Society, 75(10), 1825-1830.

[top]

ERROR CODES and CONDITIONS

KNOWN BUGS

none at this time

[top]

FUTURE PLANS

It is likely that a number of additional optional interfaces will be added to the model_mod structure. For instance, hints about how to divide the state vector into regions for parallel assimilation will need to be obtained from the model. It is planned that the interf_provided mechanism used in pert_model_copies will allow those who do not wish to support enhanced interfaces to add NULL interfaces by simply pasting in an interface block.

[top]

Terms of Use

DART software - Copyright 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: your_name_here
Revision: $Revision: 11400 $
Source: $URL: https://svn-dares-dart.cgd.ucar.edu/DART/releases/Manhattan/models/bgrid_solo/model_mod.html $
Change Date: $Date: 2017-03-24 16:38:14 -0600 (Fri, 24 Mar 2017) $
Change history:  try "svn log" or "svn diff"