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

MODULE model_mod (PBL 1D)

Contact: Jeff Anderson
Revision: $Revision: 2801 $
Source: $URL: http://subversion.ucar.edu/DAReS/DART/trunk/models/PBL_1d/model_mod.html $
Change Date: $Date: 2007-04-04 22:11:48 -0600 (Wed, 04 Apr 2007) $
Change history: try "svn log" or "svn diff"

OVERVIEW

DART interface module for the PBL 1D 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. The distinctive part of the model interfaces is the namelist.

The PBL model is a single column version of the WRF model.




OTHER MODULES USED

use types_mod
use time_manager_mod
use location_mod (threed-sphere)
use utilities_mod
use sort_mod
use random_seq_mod
use obs_kind_mod
use map_utils

also these WRF modules:

use module_model_constants
use module_initialize
use module_namelist
use module_wrf



PUBLIC INTERFACE

use module_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
The following are also public routines 
but are not used by DART code.
 get_pblh
 pert_params_time
 adjust_param_spread
 real_obs_period
 start_real_obs
 synchronize_mavail
 ok_to_nudge

NOTES

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




PUBLIC COMPONENTS

fill in text here


var = get_model_size()
 integer :: get_model_size
 

Description

Returns the size of the model as an integer.

get_model_size    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

Does a single timestep advance of the model. The input value of the vector x is the starting condition and x is 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 of filter or if the program integrate_model is to be used to advance the model state as a separate executable. In this model this is a call to code in wrf.F .

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

Given an integer index into the state vector structure, returns the associated location. A second intent(out) optional argument kind can be returned if the model has more than one type of field (which is true in this model.)

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 should be returned as 0 unless there is some problem in computing the interpolation in which case an alternate value should be returned. The itype variable is a model specific integer that specifies the type of field (of which there are many in this model).

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.

get_model_time_step    Time step of model.


call end_model()
 

Description

A NULL INTERFACE (stub) in this model at this point.



call static_init_model()
 

Description

Called to do one time initialization of the model. Reads the namelist, initializes random sequences, sets up many local derived types with sizes and info about the WRF data, allocates space and fills it in.



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

Description

Sets starting time based on namelist values.

time    Initial model time.


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

Description

Returns a model state vector, x, that is some sort of appropriate initial condition for starting up a long integration of the model. Most of the work in this model is done by calling the init_wrf() subroutine from the WRF module.

x    Initial conditions for state vector.


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

Description

Writes the model-specific attributes to a netCDF file.

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  

nc_write_model_atts    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

Writes the model variables to a netCDF file.

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  

nc_write_model_vars    Returns 0 for normal completion.
ncFileID    NetCDF file id.
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 calls init_conditions() to generate a perturbed state.

state    State vector to be perturbed.
pert_state    Perturbed state vector
interf_provided    Returned true from this model


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

Description

Pass-through to the 3-D sphere 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)
 

Description

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

Description

Pass-through to the 3-D sphere 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
 

Description

A NULL INTERFACE in this model.

ens_mean    Ensemble mean state vector


var = get_pblh(x)
 real(r8)              :: get_pblh
 real(r8), intent(in)  :: x(:)
 

Description

Get diagnosed PBL height.
var    Diagnosed PBL height.
x    State vector

subroutine pert_params_time(x)
   real(r8), intent(inout)  :: x(:)

Description

Adds a little noise to any stochastic parameters. This is separate in case we want a different distribution from the initial perturbations.

subroutine adjust_param_spread(ens, nens)
   integer, intent(in) :: nens
   real(r8), intent(inout) :: ens(:,:)

Description

Given the ensemble, find the parameters and adjust the spread to the initial values IF IT IS LESS. In other words, this is only a lower bound.

subroutine synchronize_mavail(obs, ens, ens_size, domain_size)
   integer,                  intent(in)    :: ens_size, domain_size
   real(r8),                 intent(in)    :: obs
   real(r8), dimension(ens_size,domain_size), intent(inout) :: ens

Description

Argument ens is dimension(n_ensemble,n_state). Want to either update them all independently (they are all being nudged) or pick the first one.

function ok_to_nudge(obs_kind,state_index)
   integer, intent(in)  :: obs_kind, state_index
   logical              :: ok_to_nudge

Description

Checks for physical quantity consistency between obs and state.

integer :: real_obs_period = 1800

Description

Public variable.

integer :: start_real_obs = 1800

Description

Public variable.




NAMELIST

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

 namelist / model_nml /   &
         num_est_params, est_param_types, pert_param_sd, &
         pert_init_sd, pert_init_beta_1, pert_init_beta_2, &
         maintain_initial_spread, dist_shape, pert_param_min, &
         pert_param_max, real_obs_period, start_real_obs
 

Discussion

Need descriptions here (in table form).

This namelist is read from the file input.nml





FILES




REFERENCES


ERROR CODES and CONDITIONS







KNOWN BUGS

1. The model cannot successfully use existing initial conditions. Make sure start_from_restart = .false. at all times.

2. If you want to try some parameter estimation experiments, you must become familiar with the gen_init external utility to create a filter_ics file with the proper state size.

3. Ideal runs (gen_init=.FALSE.) may not work for the MYJ schemes. This will be corrected shortly.




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_state will allow those who do not wish to support enhanced interfaces to add NULL interfaces by simply pasting in an interface block.