MODULE model_mod (WRF)

DART project logo

Jump to DART Documentation Main Index
version information for this file:
$Id: model_mod.html 7668 2015-03-05 23:08:55Z nancy $

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

Overview

DART interface module for the WRF model. This page documents the details of the module compiled into DART that interfaces with the WRF data in the state vector. There is additional overview and tutorial documentation for running a WRF/DART assimilation on this web page:

http://www.image.ucar.edu/wrfdart/tutorial

The WRF model_mod code reads the state data from a filter initial condition/restart file, typically called filter_ics or filter_restart. It also requires a NetCDF file named wrfinput_d01 in the current directory (and d02, etc. for multiple domain cases). This file must be at the same resolution and have the same surface elevation data as the files converted to create the DART initial conditions. No data will be read from this file, but the grid information must match exactly. See the wrf_to_dart documentation for more information on how to generate DART IC/restart files given WRF NetCDF files.

The model interface code supports WRF configurations with multiple domains. Data for all domains is read into the DART state vector. During the computation of the forward operators (getting the estimated observation values from each ensemble member), the search starts in the domain with the highest number, which is generally the finest nest or one of multiple finer nests. The search stops as soon as a domain contains the observation location, working its way from largest number to smallest number domain, ending with domain 1. For example, in a 4 domain case the data in the state vector that came from wrfinput_d04 is searched first, then wrfinput_d03, wrfinput_d02, and finally wrfinput_d01. The forward operator is computed from the first domain grid that contains the lat/lon of the observation. During the assimilation phase, when the state values are adjusted based on the correlations and assimilation increments, all points in all domains that are within the localization radius are adjusted, regardless of domain. The impact of an observation on the state depends only on the distance between the observation and the state vector point, and the regression coefficient based on the correlation between the distributions of the ensemble of state vector points and the ensemble of observation forward operator values.

The fields from WRF that are copied into the DART state vector are controlled by namelist. See below for the documentation on the &model_nml entries. The state vector should include all fields needed to restart a WRF run. There may be additional fields needed depending on the microphysics scheme selected. See the ascii file wrf_state_variables_table in the models/wrf directory for a list of fields that are often included in the DART state.

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 a variety of variables required in observational operators.

[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
   default_state_variables     = .true.,
   wrf_state_variables         = 'NULL',
   wrf_state_bounds            = 'NULL',
   num_domains                 = 1,
   output_state_vector         = .false.,
   calendar_type               = 3,
   num_moist_vars              = 3,
   surf_obs                    = .true.,
   soil_data                   = .true.,
   h_diab                      = .false.,
   assimilation_period_seconds = 21600,
   adv_mod_command             = './wrf.exe',
   allow_obs_below_vol         = .false.,
   vert_localization_coord     = 3,
   center_search_half_length   = 500000.,
   center_spline_grid_scale    = 10,
   circulation_pres_level      = 80000.0,
   circulation_radius          = 108000.0,
   sfc_elev_max_diff           = -1.0,
   polar                       = .false.,
   periodic_x                  = .false.,
   periodic_y                  = .false.,
   scm                         = .false.  
/

# Notes for model_nml:
# (1) vert_localization_coord must be one of:
#     1 = model level
#     2 = pressure
#     3 = height
#     4 = scale height
# (2) see bottom of this file for explanations of polar, periodic_x, 
#     periodic_y, and scm
# (3) calendar = 3 is GREGORIAN, which is what WRF uses.
# (4) if 'default_state_variables' is .true. the model_mod.f90 code will
#     fill the state variable table with the following wrf vars: 
#        U, V, W, PH, T, MU
#     you must set it to false before you change the value 
#     of 'wrf_state_variables' and have it take effect.
# (5) the format for 'wrf_state_variables' is an array of 5 strings:
#     wrf netcdf variable name, dart KIND_xxx string, type string (must be 
#     unique, will soon be obsolete, we hope), 'UPDATE', and '999' if the 
#     array is part of all domains.  otherwise, it is a string with the domain
#     numbers (e.g. '12' for domains 1 and 2, '13' for domains 1 and 3).
#   example:
# wrf_state_variables='U','KIND_U_WIND_COMPONENT','TYPE_U','UPDATE','999',
#                     'V','KIND_V_WIND_COMPONENT','TYPE_V','UPDATE','999',
#                     'W','KIND_VERTICAL_VELOCITY','TYPE_W','UPDATE','999',
#                     'T','KIND_POTENTIAL_TEMPERATURE','TYPE_T','UPDATE','999',
#                     'PH','KIND_GEOPOTENTIAL_HEIGHT','TYPE_GZ','UPDATE','999',
#                     'MU','KIND_PRESSURE','TYPE_MU','UPDATE','999',
#                     'QVAPOR','KIND_VAPOR_MIXING_RATIO','TYPE_QV','UPDATE','999',
#                     'QCLOUD','KIND_CLOUD_LIQUID_WATER','TYPE_QC','UPDATE','999',
#                     'QRAIN','KIND_RAINWATER_MIXING_RATIO','TYPE_QR','UPDATE','999',
#                     'U10','KIND_U_WIND_COMPONENT','TYPE_U10','UPDATE','999',
#                     'V10','KIND_V_WIND_COMPONENT','TYPE_V10','UPDATE','999',
#                     'T2','KIND_TEMPERATURE','TYPE_T2','UPDATE','999',
#                     'TH2','KIND_POTENTIAL_TEMPERATURE','TYPE_TH2','UPDATE','999',
#                     'Q2','KIND_SPECIFIC_HUMIDITY','TYPE_Q2','UPDATE','999',
#                     'PSFC','KIND_PRESSURE','TYPE_PS','UPDATE','999',
# (6) the format for 'wrf_state_bounds' is an array of 4 strings:
#     wrf netcdf variable name, minimum value, maximum value, and either
#     FAIL or CLAMP.  FAIL will halt the program if an out of range value
#     is detected.  CLAMP will set out of range values to the min or max.
#     The special string 'NULL' will map to plus or minus infinity and will
#     not change the values.  arrays not listed in this table will not
#     be changed as they are read or written.
#
#
# polar and periodic_x are used in global wrf.  if polar is true, the 
# grid interpolation routines will wrap over the north and south poles.  
# if periodic_x is true, when the east and west edges of the grid are
# reached the interpolation will wrap.  note this is a separate issue
# from regional models which cross the GMT line; those grids are marked
# as having a negative offset and do not need to wrap; this flag controls
# what happens when the edges of the grid are reached.

# the scm flag is used for the 'single column model' version of WRF.
# it needs the periodic_x and periodic_y flags set to true, in which
# case the X and Y directions are periodic; no collapsing of the grid
# into a single location like the 3d-spherical polar flag implies.


Item Type Description
default_state_variables logical If .true., the dart state vector contains the fields U, V, W, PH, T, MU, in that order, and only those. Any values listed in the wrf_state_variables namelist item will be ignored.
wrf_state_variables character(:, 5) A 2D array of strings, 5 per wrf array to be added to the dart state vector. If default_state_variables is .true., this is ignored. When .false., this list of array names controls which arrays and the order that they are added to the state vector. The 5 strings are:
  1. WRF field name - must match netcdf name exactly
  2. DART KIND name - must match a valid DART KIND_xxx exactly
  3. TYPE_NN - will hopefully be obsolete, but for now NN should match the field name.
  4. the string UPDATE. at some future point, non-updatable fields may become part of the state vector.
  5. A numeric string listing the domain numbers this array is part of. The specical string 999 means all domains. For example, '12' means domains 1 and 2, '13' means 1 and 3.
wrf_state_bounds character(:, 4) A 2D array of strings, 4 per wrf array. During the copy of data to and from the wrf netcdf file, variables listed here will have minimum and maximum values enforced. The 4 strings are:
  1. WRF field name - must match netcdf name exactly
  2. Minimum -- specified as a string but must be a numeric value (e.g. '0.1') Can be 'NULL' to allow any minimum value.
  3. Maximum -- specified as a string but must be a numeric value (e.g. '0.1') Can be 'NULL' to allow any maximum value.
  4. Action -- valid strings are 'CLAMP', 'FAIL'. 'FAIL' means if a value is found outside the range, the code fails with an error. 'CLAMP' simply sets the out of range values to the given minimum or maximum without error.
output_state_vector logical Controls the output to netCDF files. If .true., output the raw dart state vector. If .false., output the prognostic flavor (gridded data) for easier plotting (recommended).
num_domains integer Total number of WRF domains, including nested domains.
calendar_type integer Calendar type. Should be 3 (GREGORIAN) for WRF.
assimilation_period_seconds integer The time (in seconds) between assimilations. This is modified if necessary to be an integer multiple of the underlying model timestep.
periodic_x logical If .true., the grid is periodic in longitude, and points above the last grid cell and points below the first grid cell are wrapped. Note this is not the same as a grid which crosses the prime meridian. WRF handles that with an offset in longitude and points beyond the last grid index are outside the domain.
periodic_y logical Used for the Single Column Model to make the grid wrap in Y (see scm below). This is NOT the same as wrapping in latitude (see polar below).
polar logical If .true., points at the poles are wrapped across the grid. It is not clear this is a good idea since the grid is degnerate here.
scm logical If .true. the Single Column Model is assumed. The grid is a single vertical column, and there are 9 cells arranged in a 3x3 grid. See the WRF documentation for more information on this configuration. periodic_x and periodic_y should also be .true. in this case.
sfc_elev_max_diff real(r8) If > 0, the maximum difference, in meters, between an observation marked as a 'surface obs' as the vertical type (with the surface elevation, in meters, as the numerical vertical location), and the surface elevation as defined by the model. Observations further away from the surface than this threshold are rejected and not assimilated. If the value is negative, this test is skipped.
allow_obs_below_vol logical If .false. then if an observation with a vertical coordinate of pressure or height (i.e. not a surface observation) is below the lowest 3d sigma level, it is outside the field volume and the interpolation routine rejects it. If this is set to .true. and the observation is above the surface elevation but below the lowest field volume level, the code will extrapolate downward from data values at levels 1 and 2.
center_search_half_length real(r8) The model_mod now contains two schemes for searching for a vortex center location. If the old scheme is compiled in, then this and the center_spline_grid_scale namelist items are used. (Search code for 'use_old_vortex'.) Half length (in meters) of a square box for searching the vortex center.
center_spline_grid_scale integer The model_mod now contains two schemes for searching for a vortex center location. If the old scheme is compiled in, then this and the center_search_half_length namelist items are used. (Search code for 'use_old_vortex'.) Ratio of refining grid for spline-interpolation in determining the vortex center.
circulation_pres_level real(r8) The model_mod now contains two schemes for searching for a vortex center location. If the new scheme is compiled in, then this and the circulation_radius namelist items are used. (Search code for 'use_old_vortex'.) Pressure, in pascals, of the level at which the circulation is computed when searching for the vortex center.
circulation_radius real(r8) The model_mod now contains two schemes for searching for a vortex center location. If the new scheme is compiled in, then this and the circulation_pres_level namelist items are used. (Search code for 'use_old_vortex'.) Radius, in meters, of the circle over which the circulation calculation is done when searching for the vortex center.
vert_localization_coord integer Vertical coordinate for vertical localization.
  • 1 = model level
  • 2 = pressure (in pascals)
  • 3 = height (in meters)
  • 4 = scale height (unitless)
surf_obs logical DEPRECATED -- has no effect on the code and will be removed.
soil_data logical DEPRECATED -- has no effect on the code and will be removed.
h_diab logical DEPRECATED -- has no effect on the code and will be removed.
num_moist_vars integer DEPRECATED -- has no effect on the code and will be removed.
adv_mod_command character(len=32) DEPRECATED -- A variable with the same name has moved to the &dart_to_wrf namelist and does the same function from there. This one has no effect on the code and will be removed soon.


[top]

OTHER MODULES USED

types_mod
time_manager_mod
threed_sphere/location_mod
utilities_mod
obs_kind_mod
map_utils
netcdf
typesizes

PUBLIC INTERFACES

use model_mod, only : get_model_size
 get_state_meta_data
 model_interpolate
 get_model_time_step
 static_init_model
 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
 adv_1step
 end_model
 init_time
 init_conditions

The last 4 interfaces are only required for low-order models where advancing the model can be done by a call to a subroutine. The WRF model only advances by executing the program wrf.exe. Thus the last 4 interfaces only appear as stubs in the wrf module.

The interface pert_model_state is presently not provided for WRF. The initial ensemble has to be generated off-line. If coherent structures are not required, the filter can generate an ensemble with uncorrelated random Gaussian noise of 0.002. This is of course not appropriate for a model like WRF which has variables expressed in a wide range of scales. It is thus recommended to generate the initial ensemble off-line, perhaps with the tools provided in models/wrf/PERTURB/3DVAR-COVAR.

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. This includes all nested domains.

model_size The length of the model state vector.


call get_state_meta_data (index_in, location [, var_type_out, id_out] )
integer,             intent(in)  :: index_in
type(location_type), intent(out) :: location
integer, optional,   intent(out) :: var_type_out
integer, optional,   intent(out) :: id_out

Returns metadata about a given element, indexed by index_in, in the model state vector. The location defines where the state variable is located while the type of the variable (for instance temperature, or u wind component) is returned by var_type. The integer values used to indicate different variable types in var_type are themselves defined as public interfaces to model_mod if required. The last optional argument is the wrf domain identification number - obviously this is unique to the WRF version of this required routine.

index_in    Index of state vector element about which information is requested.
location Returns location of indexed state variable. The location should use a location_mod that is appropriate for the model domain. For realistic atmospheric models, for instance, a three-dimensional spherical location module that can represent height in a variety of ways is provided.
var_type_out Returns the type of the indexed state variable as an optional argument.
id_out Returns the wrf domain identification number of the indexed state variable as an optional argument.


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

Given model state, returns the value of observation type interpolated to a given location by a method of the model's choosing. All observation kinds defined in obs_kind_mod are supported. In the case where the observational operator is not defined at the given location (e.g. the observation is below the model surface or outside the domain), obs_val is returned as -888888.0 and istatus = 1. Otherwise, istatus = 0. The interpolation is performed in the domain with the highest resolution containing the observation.

x A model state vector.
location    Location to which to interpolate.
obs_kind Integer indexing which type of observation is to be interpolated.
obs_val The interpolated value from the model.
istatus Integer flag indicating the result of the interpolation.


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

Returns the model base time step as a time_type. For the model wrf, it returns the time step used for domain 1 (usually the largest time step among all domains because domain 1 is the coarser grid). The time step is read from the wrfinput_d01 file (used to be namelist.input). In the long run, a more general extended interface may be required that specifies the models range of time stepping possibilities.

var    Smallest time step of model.


call static_init_model()

Used for runtime initialization of the model. This is the first call made to the model by any DART compliant assimilation routine. It reads the model namelist parameters, set the calendar type (the GREGORIAN calendar is used with the WRF model), and determine the dart vector length. This subroutine requires that wrfinput_d01, wrfinput_d02, ... (one file for each domain) be present in the working directory to retrieve model information (grid dimensions and spacing, time step, pressure at the top of the model, map projection parameters, etc).



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.


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. Models should override this if some structure is required for perturbations or if the magnitude of perturbations in DART is too large. It is thus recommended to generate the initial ensemble off-line, perhaps with the tools provided in models/wrf/PERTURB/3DVAR-COVAR.

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


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

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)

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

Calls the 3-D sphere locations module to get a list of potentially close state vector points, and again for unassimilated observations. See get_close_obs() for the documentation of the locations module version of this code. Then, if vertical localization is enabled, this code converts all vertical locations to the selected vertical type (&model_nml::vert_localization_coord). It then computes a real 3D distance and returns it to the calling code.



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

A local copy is kept and used during other computations in the model_mod code.

ens_mean    Ensemble mean state vector


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

This operation is not defined for the WRF model. If called it will throw a fatal error.

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


call end_model( )

Called when use of a model is completed to clean up storage, etc. A stub is provided for the WRF model.



call init_time(i_time)
type(time_type),        intent(in)  ::  i_time 

Not supported for the WRF model.



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

Not supported for the WRF model. Will throw a fatal error if called.

x Model state vector.

[top]

FILES

[top]

REFERENCES

http://www.mmm.ucar.edu/wrf/users/docs/user_guide/contents.html
[top]

ERROR CODES and CONDITIONS

RoutineMessageComment
static_init_model
nc_write_model_vars
num_moist_vars is too large The maximum number of moist variable in WRFV2.1 is 7
static_init_model 'Please put wrfinput_d0'//idom//' in the work directory.' One of the wrfinput_d0# is missing in the work directory
static_init_model Map projection no supported Try PROJ_LATLON(0), PROJ_LC(1), PROJ_PS(2), PROJ_MERC(3)
static_init_model
nc_write_model_atts
nc_write_model_vars
Various NetCDF-f90 interface error messages From one of the NetCDF calls in the named routine
get_state_meta_data dart index out of range Unlikely. Would indicate a serious bug in the code
model_interpolate
get_dist_wrf
wrong option for which_vert See which_vert description in location/threed_sphere/location_mod.html
model_interpolate 'do not recognize obs kind ',obs_kind See list in 'use obs_kind_mod' statement in model_mod.f90
get_wrf_index 'Indices ',i,j,k,' exceed grid dimensions: ',#1,#2,#3 One of the grid indices exceeds the corresponding dimension for the var_type input. Unlikely to happen but would indicate a serious bug in the code
get_dist_wrf Unable to define vloc The vertical location is below the model surface or above the model lid
nc_write_model_atts Time dimension ID # must match Unlimited Dimension ID # NetCDF file writing error
read_dt_from_wrf_nml Please put namelist.input in the work directory The file namelist.input is missing in the work directory
read_dt_from_wrf_nml 'max_dom in namelist.input = ',max_dom'num_domains in input.nml = ',num_domains'Make them consistent.' The number of WRF domains in namelist.input and in input.nml do not match

KNOWN BUGS

Only the Lambert projection (MAP_PROJ = 1) is working. Other map projections (0=none, 2=polar, 3=Mercator) are likely to break the code.

[top]

FUTURE PLANS

none at this time

[top]

PRIVATE COMPONENTS


ind = get_wrf_index(i,j,k,var_type,id)
integer                               :: get_wrf_index 
integer,                  intent(in)  :: i,j,k,var_type,id

Given grid indices, variable type, and domain identification number, returns the index in the model state vector as an integer.



dist = get_dist_wrf( i, j, k, var_type, o_loc, id, x )
real(r8)                              :: get_dist_wrf 
integer,                  intent(in)  :: i,j,k,var_type,id

type(location_type),      intent(in)  :: o_loc 
real(r8), dimension(:),   intent(in)  :: x 

Given grid indices, variable type, and domain identification number, computes the distance to the observation at location o_loc.



call get_wrf_horizontal_location( i, j, var_type, id, long, lat )
integer,                  intent(in)  :: i,j,var_type,id 
real(r8),                 intent(out) :: long, lat 

Given horizontal grid indices, variable type, and domain identification number, returns the longitude and latitude in degrees.



call toGrid (x, j, dx, dxm)
real(r8),                 intent(in)  :: x 
real(r8),                 intent(out) :: dx, dxm 
integer ,                 intent(out) :: j 

Given position x, find nearest grid point j (but smaller than or equal to x) and calculate its distance to grid j and j+1.



call pres_to_zk(pres, mdl_v, n3, zk)
real(r8),                  intent(in)  :: pres 
integer                    intent(in)  :: n3 
real(r8), dimension(0:n3), intent(in)  :: mdl_v 
real(r8),                  intent(out) :: zk 

Calculate the position zk on half (mass) levels in the profile mdl_v corresponding to pressure location pres.



call height_to_zk(obs_v, mdl_v, n3, zk)
real(r8),                  intent(in)  :: obs_v 
integer                    intent(in)  :: n3 
real(r8), dimension(0:n3), intent(in)  :: mdl_v 
real(r8),                  intent(out) :: zk 

Calculate the position zk on half (mass) levels in the profile mdl_v corresponding to height location obs_v.



call get_model_pressure_profile(i,j,dx,dy,dxm,dym,n,x,id,v_p)
integer                   intent(in)  :: i,j,n,id 
real(r8),                 intent(in)  :: dx,dy,dxm,dym 
real(r8), dimension(:),   intent(in)  :: x 
real(r8), dimension(0:n), intent(out) :: v_p 

Extract the pressure profile at position (i+dx, j+dy) on the non-staggered vertical grid (half (mass) levels).



vloc = model_pressure(i,j,k,id,var_type,x)
real(r8),                             :: model_pressure 
integer                   intent(in)  :: i,j,k,id,var_type

real(r8), dimension(:),   intent(in)  :: x 

Calculate the pressure at grid point (i,j,k), domain id. The grid is defined according to var_type.



pres = model_pressure_t(i,j,k,id,x)
real(r8)                              :: model_pressure_t 
integer                   intent(in)  :: i,j,k,id 
real(r8), dimension(:),   intent(in)  :: x 

Calculate the pressure at grid point (i,j,k), domain id, on mass point (half (mass) levels, T-point).



rho = model_rho_t(i,j,k,id,x)
real(r8)                              :: model_rho_t 
integer                   intent(in)  :: i,j,k,id 
real(r8), dimension(:),   intent(in)  :: x 

Calculate the total density at grid point (i,j,k), domain id, on mass point (half (mass) levels, T-point).



call get_model_height_profile(i,j,dx,dy,dxm,dym,n,x,id,fld)
integer                   intent(in)  :: i,j,n,id 
real(r8),                 intent(in)  :: dx,dy,dxm,dym 
real(r8), dimension(:),   intent(in)  :: x 
real(r8), dimension(n),   intent(out) :: fld 

Extract the height profile at position (i+dx, j+dy) on the non-staggered vertical grid.



vloc = model_height(i,j,k,id,var_type,x)
real(r8)                              :: model_height
integer                   intent(in)  :: i,j,k,id,var_type
real(r8), dimension(:),   intent(in)  :: x

Calculate the height at grid point (i,j,k), domain id. The grid is defined according to var_type.



call read_dt_from_wrf_nml( )

Read the wrf model time step from namelist.input for all domains.


[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: Hui Liu, Glen Romine, Chris Synder, David Dowell
Revision: $Revision: 7668 $
Source: $URL: https://svn-dares-dart.cgd.ucar.edu/DART/releases/Lanai/models/wrf/model_mod.html $
Change Date: $Date: 2015-03-05 16:08:55 -0700 (Thu, 05 Mar 2015) $
Change history:  try "svn log" or "svn diff"