MODULE model_mod (CAM)

DART project logo

Jump to DART Documentation Main Index
version information for this file:
$Id: model_mod.html 11427 2017-03-31 21:24:30Z nancy@ucar.edu $

NAMELIST / SETUP VARIATIONS / INTERFACES / PRIVATE COMPONENTS / DISCUSSION / FILES / REFERENCES / CONTRIBUTORS / ERRORS / BUGS / PLANS / TERMS OF USE

Overview

The DART system supports data assimilation into the Community Atmosphere Model, CAM, which is the atmospheric component of the Community Earth System Model (CESM). This DART interface is being used by graduate students, post-graduates, and scientists at universities and research labs to conduct data assimilation reseearch. Others are using the products of data assimilation (analyses), which were produced here at NCAR using CESM+DART, to conduct related research. The variety of research can be sampled on the DART Publications page.

"CAM" refers to a family of related atmospheric components, which can be built with 2 independent main characteristics. CESM labels these as:

   'resolution' = horizontal (not vertical) grid AND dynamical core (fluid dynamics equations on that grid)
      3 supported dycores; eulerian, FV, and SE 
      SE refined grids will have some support, but by their nature invite the use 
        of user defined grid and map files.
   'compset' = vertical grid AND parameterizations (aka physics):
      Parameterization is the equations describing a physical process like
         convection, radiation, chemistry, ...
      Vertical grid is determined by the needs of the chosen parameterizations.
         Spacing and height of the top (ptop) vary.
      The combinations of parameterizations and vertical grids are named:
         CAM3.5, CAM5, CAM#, ...
         WACCM, WACCM#, WACCM-X, 
         CAM-Chem,
         ...
There are minor characteristics choices within each of these, but only chemistry choices in WACCM and CAM-Chem have an impact on DART. As of April, 2015, all of these variants are handled by the same model_mod.f90, namelist, and build scripts, with differences in the assimilation set up described here.

This DART+CAM interface has the following features.

In addition to the standard DART package there are ensembles of initial condition files at the large file website http://www.image.ucar.edu/pub/DART/CAM/ that are helpful for interfacing CAM with DART. In the current (2015) mode, CESM+DART can easily be started from a single model state, which is perturbed to create an ensemble of the desired size. A spin-up period is then required to allow the ensemble members to diverge.

Sample sets of observations, which can be used with CESM+DART assimilations, can be found at http://www.image.ucar.edu/pub/DART/Obs_sets/ of which the NCEP BUFR observations are the most widely used.

Experience on a variety of machines has shown that it is a very good idea to make sure your run-time environment has the following:

limit stacksize unlimited
limit datasize unlimited

This page contains the documentation for the DART interface module for the CAM and WACCM models, using the dynamical cores listed above. This implementation uses the CAM initial files (not restarts) for transferring the model state to/from the filter. This may change in future versions, but probably only for CAM-SE. The reasons for this include:

  1. The contents of the restart files vary depending on both the model release version and the physics packages selected.
  2. There is no metadata describing the variables in the restart files. Some information can be tracked down in the atm.log file, but not all of it.
  3. The restart files (for non-chemistry model versions) are much larger than the initial files (and we need to deal with an ensemble of them).
  4. The temperature on the restart files is virtual equivalent potential temperature, which requires (at least) surface pressure, specific humidity, and sensible temperature to calculate.
  5. CAM does not call the initialization routines when restart files are used, so fields which are not modified by DART may be inconsistent with fields which are.
  6. If DART modifies the contents of the .r. restart file, it might also need to modify the contents of the .rs. restart file, which has similar characteristics (1-3 above) to the .r. file.

The DART interfaces to CAM and many of the other CESM components have been integrated with the CESM set-up and run scripts. Unlike previous versions of DART-CAM, CESM runs using its normal scripts, then stops and calls a DART script, which runs a single assimilation step, then returns to the CESM run script to continue the model advances. See the CESM interface documentation for more information on running DART with CESM. Due to the complexity of the CESM software environment, the versions of CESM which can be used for assimilation are more restricted than previously. Each supported CESM version has similar, but unique, sets of set-up scripts and CESM SourceMods. Those generally do not affect the cam/model_mod.f90 interface. Current (April, 2015) set-up scripts are:

The DART state vector should include all prognostic variables in the CAM initial files which cannot be calculated directly from other prognostic variables. In practice the state vector sometimes contains derived quantities to enable DART to compute forward operators (expected observation values) efficiently. The derived quantities are often overwritten when the model runs the next timestep, so the work DART does to update them is wasted work.

Expected observation values on pressure, scale height, height or model levels can be requested from model_interpolate. Surface observations can not yet be interpolated, due to the difference between the model surface and the earth's surface where the observations are made. Model_interpolate can be queried for any (non-surface) variable in the state vector (which are variables native to CAM) plus pressure on height levels. The default state vector is PS, T, U, V, Q, CLDLIQ, CLDICE and any tracers or chemicals needed for a given study. Variables which are not in the initial file can be added, but minor modifications to model_mod.f90 and CAM may be necessary.

The 16 public interfaces in model_mod are standardized for all DART compliant models. These interfaces allow DART to 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 by 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. The values shown here are the default values.

&model_nml
   output_state_vector = .false.,
   model_config_file   = 'caminput.nc',
   cam_phis            = 'cam_phis.nc',
   state_num_0d        = 0,
   state_num_1d        = 0,
   state_num_2d        = 1,
   state_num_3d        = 6,
   state_names_0d      = '',
   state_names_1d      = '',
   state_names_2d      = 'PS'
   state_names_3d      = 'T', 'U', 'V', 'Q', 'CLDLIQ','CLDICE',
   which_vert_1d       = -2 ,
   which_vert_2d       = -1 , 
   which_vert_3d       = 6*1  ,
   pert_names          = '',
   pert_sd             = -888888.0,
   pert_base_vals      = -888888.0,
   vert_coord          = 'pressure'
   highest_obs_pressure_Pa   = 1000.0,
   highest_state_pressure_Pa = 9400.0,
   max_obs_lat_degree  = 90.0,
   Time_step_seconds   = 21600,
   Time_step_days      = 0,
   print_details       = .false.
   model_version       = '4.0',
   impact_only_same_kind = '',
   /


The names of the fields to put into the state vector come from the CAM initial file field names. The specification of lists of names and numbers for the various dimensions enables the very flexible definition of the state vector. It can be done via the namelist, instead of recompiling DART for each different set. In the CESM+DART framework, filter initial condition files are created based on the state vector defined in the namelist.

The dimension of these lists is currently hardwired to size 100. If more fields need to be assimilated (e.g. many chemical species), look for the integer parameter MAX_STATE_NAMES in the source code and change it to a large enough value and recompile DART.

The values for which_vert_#d are described in location_mod: location_mod.html.

Item Type Description
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).
model_config_file character(len=128) CAM initial file used to provide configuration information, like the resolution of the grid, number of vertical levels, whether fields are staggered or not, etc.
cam_phis character(len=128) CAM topography file for the Eulerian and Finite Volume versions. Not used in the Spectral Element version.
cs_grid_file character(len=128) CAM Spectral Element grid file. Not used in the Eulerian or Finite Volume versions.
homme_map_file character(len=128) CAM Spectral Element mapping file written by CAM-SE. Not used in the Eulerian or Finite Volume versions.
state_num_#d,
#=0,1,2,3
integer Numbers of fields of various dimensions to put into the state vector. Note that CAM-SE fields have only 1 horizontal dimension on the initial files, due to the cubed-sphere grid being not logically rectangular. "Dimensions" in cam/model_mod.f90 generally refers to dimensions on the initial file, not spatial dimensions of the fields.
state_names_#d,
#=0,1,2,3
character(len=8), dimension(100) Names of fields of various dimensions to put into the state vector.
which_vert_#d,
#=1,2,3
integer, dimension(100) Vertical location types of fields in state_names_#d. See the 3D sphere location documentation for the mapping of integer values to vertical location types.
pert_names character(len=8), dimension(100) To make filter generate an ensemble from a single model state by randomly perturbing it, list the field(s) to be perturbed here (and see the DART namelist settings for a perturbed initial ensemble).
pert_sd real(r8), dimension(100) If positive, it's the standard deviation of the perturbation for each field in the pert_names list (filter). Unused unless pert_names is set.
pert_base_vals real(r8), dimension(100) If pert_sd is positive, this the list of values to which the field(s) listed in pert_names will be reset if filter is told to create an ensemble from a single state vector. Otherwise, it's is the list of values to use for each ensemble member when perturbing the single field named in pert_names. Unused unless pert_names is set and pert_base_vals is not the DART missing value.
max_obs_lat_degree real(r8) Observations closer to the poles than this latitude will be ignored.
vert_coord character(len=8) The vertical coordinate to which all vertical locations are converted in model_mod. "log_invP" ("scale height" especially for WACCM) or "pressure".
highest_obs_pressure_Pa real(r8) Observations higher than this pressure are ignored. NOTE that this has a non-backwards incompatible change from previous versions. It is now specified in Pascals, not millibars. Divide by 100 to convert the units.
highest_state_pressure_Pa real(r8) Influence of all observations on model points higher than this is reduced. NOTE that this has a non-backwards incompatible change from previous versions. It is now specified in Pascals, not millibars. Divide by 100 to convert the units to mb. Details of calculating the minimum value recommended for a given vertical grid and set of DART namelist variables can be found in model_top_issues.pdf,
Time_step_seconds real(r8) Minimum forecast duration (the part < 1 day). CESM1_2_1_setup_{hybrid,advanced} assume that this is 6 hours (21600 s).
Time_step_days real(r8) Minimum forecast duration (the part > 24*3600 sec)
print_details logical If true, print out detailed information about the sizes, shapes, offsets, etc. of items in the CAM state vector. If false, just print out the names of the fields as they are read into the state vector.
model_version character(len=128) Not used in CESM+DART; deprecated. The number of the CAM version being used, i.e. '3.0.7'.
impact_only_same_kind character(len=32) Name of one observation kind which can only affect state variables of the same kind.


Setup Variations

The variants of CAM require slight changes to the setup scripts (in $DART/models/cam/shell_scripts) and in the namelists (in $DART/models/cam/work/input.nml). From the DART side, assimilations can be started from a pre-existing ensemble, or an ensemble can be created from a single initial file before the first assimilation. In addition, there are setup differences between 'perfect model' runs, which are used to generate synthetic observations, and assimilation runs. Those differences are extensive enough that they've been coded into separate setup scripts:

Since the CESM compset and resolution, and the initial ensemble source are essentially independent of each other, changes for each of those may need to be combined to perform the desired setup.

:

The default values in work/input.nml and shell_scripts/CESM1_2_1_setup_{pmo,hybrid} are set up for a CAM-FV, single assimilation cycle using the default values as found in model_mod.f90 and starting from a single model state, which must be perturbed into an ensemble. The following are suggestions for setting it up for other assimilations. Namelist variables listed here might be in any namelist within input.nml.

CAM-FV

If built with the FV dy-core, the number of model top levels with extra diffusion in CAM is controlled by div24del2flag. The recommended minium values of highest_state_pressure_Pa come from that variable, and cutoff*vert_normalization_X:

 
   2    ("div2") -> 2 levels  -> highest_state_pressure_Pa =  9400. Pa
   4,24 ("del2") -> 3 levels  -> highest_state_pressure_Pa = 10500. Pa
   vert_coord          = 'pressure'
   state_num_1d        = 0,
   state_num_2d        = 1,
   state_num_3d        = 6,
   state_names_1d      = ''
   state_names_2d      = 'PS'
   state_names_3d      = 'T', 'US', 'VS', 'Q', 'CLDLIQ', 'CLDICE'
   which_vert_1d       = 0,
   which_vert_2d       = -1,
   which_vert_3d       = 6*1,
   highest_state_pressure_Pa = 9400. or 10500. 

CAM-SE

There's an existing ensemble, so see Continuing to start from it instead of a single state. To set up a "1-degree" CAM-SE assimilation CESM1_2_1_setup_hybrid:

   setenv resolution  ne30_g16  
   setenv refcase     SE30_Og16
   setenv refyear     2005
   setenv refmon      08
   setenv refday      01
input.nml:
   approximate_distance = .FALSE.
   vert_coord          = 'pressure'
   state_num_1d        = 1,
   state_num_2d        = 6,
   state_num_3d        = 0,
   state_names_1d      = 'PS'
   state_names_2d      = 'T','U','V','Q','CLDLIQ','CLDICE'
   state_names_3d      = ''
   which_vert_1d       = -1,
   which_vert_2d       = 6*1,
   which_vert_3d       = 0,
   highest_obs_pressure_Pa   = 1000.,
   highest_state_pressure_Pa = 10500.,

Variable resolution CAM-SE

To set up a variable resolution CAM-SE assimilation (as of April 2015) there are many changes to both the CESM code tree and the DART setup scripts. This is for very advanced users, so please contact dart @ ucar dot edu or raeder @ ucar dot edu for scripts and guidance.

WACCM

WACCM[#][-X] has a much higher top than the CAM versions, which requires the use of scale height as the vertical coordinate, instead of pressure, during assimilation. One impact of the high top is that the number of top model levels with extra diffusion in the FV version is different than in the low-topped CAM-FV, so the div24del2flag options lead to the following minimum values for highest_state_pressure_Pa:

 
   2    ("div2") -> 3 levels  -> highest_state_pressure_Pa = 0.01 Pa
   4,24 ("del2") -> 4 levels  -> highest_state_pressure_Pa = 0.02 Pa
The best choices of vert_normalization_scale_height, cutoff, and highest_state_pressure_Pa are still being investigated (April, 2015), and may depend on the observation distribution being assimilated.

WACCM is also typically run with coarser horizontal resolution. There's an existing 2-degree ensemble, so see Continuing to start from it, instead of a single state. If you use this, ignore any existing inflation restart file and tell DART to make its own in the first cycle in input.nml:

   inf_initial_from_restart    = .false.,                 .false.,
   inf_sd_initial_from_restart = .false.,                 .false.,

In any case, make the following changes (or similar) to convert from a CAM setup to a WACCM setup. CESM1_2_1_setup_hybrid:

   setenv compset     F_2000_WACCM
   setenv resolution  f19_f19  
   setenv refcase     FV1.9x2.5_WACCM4
   setenv refyear     2008
   setenv refmon      12
   setenv refday      20
input.nml:
   vert_normalization_scale_height = 2.5
   vert_coord                = 'log_invP'
   highest_obs_pressure_Pa   = .001,
   highest_state_pressure_Pa = .01,
If built with the SE dy-core (warning; experimental), then 4 levels will have extra diffusion, and also see here.

If there are problems with instability in the WACCM foreasts, try changing some of the following parameters in either the user_nl_cam section of the setup script or input.nml.

Continuing after the first cycle

After the first forecast+assimilation cycle, using an ensemble created from a single file, it is necessary to change to the 'continuing' mode, where CAM will not perform all of its startup procedures and DART will use the most recent ensemble. This example applies to an assimiation using prior inflation (inf_... = .true.). If posterior inflation were needed, then the 2nd column of infl_... would be set to "true".

input.nml:
   start_from_restart       = .true.,
   restart_in_file_name     = "filter_ics",
   single_restart_file_in  = .false.,

   inf_initial_from_restart    = .true.,                 .false.,
   inf_sd_initial_from_restart = .true.,                 .false.,

Combining multiple cycles into one job

CESM1_2_1_setup_{hybrid,pmo} are set up in the default cycling mode, where each submitted job performs one model advance and one assimilation, then resubmits the next cycle as a new job. For long series of cycles, this can result in a lot of time waiting in the queue for short jobs to run. This can be prevented by using the 'cycles' scripts generated by CESM1_2_1_setup_advanced (instead of ..._hybrid). This mode is described in the models/cam/doc/README.

[top]

OTHER MODULES USED

mpi_utilities_mod
netcdf
obs_kind_mod
random_seq_mod
threed_cartesian/xyz_location_mod
threed_sphere/location_mod
time_manager_mod
types_mod
typeSizes
utilities_mod

[top]

PUBLIC INTERFACES

FILTER INTERFACE

Here they are listed alphabetically, except for the last 4 interfaces, which are only required for low-order models where advancing the model can be done by a call to a subroutine. The last 4 interfaces only appear as stubs in the CAM module. Following this list they are listed in the order in which they appear in model_mod.f90, which is separated into sections corresponding to the main functions of model_mod.f90: static_init_model, I/O, model_interpolate, vector<->field translations, get_close_obs, and utility routines.

use model_mod, only : ens_mean_for_model
 get_close_maxdist_init
 get_close_obs_init
 get_close_obs
 get_model_time_step
 get_model_size
 get_state_meta_data
 model_interpolate
 nc_write_model_atts
 nc_write_model_vars
 pert_model_state
 static_init_model
 adv_1step
 end_model
 init_time
 init_conditions

The interface pert_model_state is provided for CAM, and it allows each field of the state vector to be randomly perturbed using a separate standard deviation. Even this may not be flexible enough to handle variables such as specific humidity, which can vary by orders of magnitude from the surface to the top of the model.

OTHER PUBLIC INTERFACES

The list of interfaces above is part of the list of public routines in this module. The rest of the public list are used by programs other than filter; dart_to_cam, cam_to_dart, etc.

use model_mod, only : model_type
 prog_var_to_vector
 read_cam_init
 vector_to_prog_var
 write_cam_init

Namelist interface &model_nml is read from file input.nml.



model_type

A type for the CAM model. It consists of pointers to model variables having ranks (non-time dimensions) of 0...3. Note that CAM-SE fields have rank = spatial_dims - 1.

call static_init_model( )

Used for run time initialization of the model. This is the first call made to the model by any DART compliant assimilation routine. It reads the model_mod namelist parameters, sets the calendar type (the GREGORIAN calendar is used with the CAM model), and determines the dart vector length, among other things. This subroutine requires that caminput.nc, (or the names in namelist variable model_config_file) be present in the working directory to retrieve model information (grid dimensions and spacing including the vertical hybrid coordinate coefficients, time step, and Gaussian weights, etc).

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

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 DART KIND of the variable (for instance temperature, or u wind component) is returned by var_kind.

index_in Index into the long state vector.
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_kind Returns the type of the indexed state variable as an optional argument.

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

Saves a local copy of the ensemble means which can be used for computing vertical heights, for example.

ens_mean   Ensemble mean state vector

model_size = get_model_size( )
integer                               ::  get_model_size 

Returns the length of the model state vector as an integer.

ts = get_model_time_step( )
type(time_type)                       ::  get_model_time_step 

Returns the nominal useful forecast length for CAM, not CAM's internal time step (dtime). get_model_time_step is set by the input.nml &model_nml:Time_step_seconds,Time_step_days

get_model_time_step also indirectly defines the assimilation window. Every observation within +/- half the get_model_time_step of the current model time will be assimilated. The current CESM+DART setup scripts also implicitly use this to determine which set of observations to use. If the observation sequence file has observations beyond the end of the current assimilation window, DART will try to advance CAM. This is not supported under the current configuration. Similarly, if the observation sequence file has observations before the assimilation window, DART will fail because it cannot move the model backwards in time. Consequently, the observation sequence files have to be preprocessed to contain just the observations to be assimilated without advancing the model.

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 identified by ncFileID. DART uses the NetCDF format to output assimilation diagnostic information.

ncFileId    Integer file descriptor opened to NetCDF file.
ierr Returned error code.

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    Integer file descriptor to opened NetCDF file.
statevec State vector.
copyindex Integer index to which copy is to be written.
timeindex Integer index of which time in the file is being written.
ierr Returned error code. success == 0, failure == -1

call read_cam_init(file_name,var)
character(len, intent(in)     :: file_name
type(model_type), intent(out) :: var

Reads state vector fields from a CAM initial file. Fields are specified in model_mod.nml.

file_name   CAM initial file name.
var   Structure to hold the state vector fields read from file_name.

call write_cam_init(file_name,model_time, var)
character, intent(in)                     :: file_name
type(time_type), intent(in)               :: model_time 
type(model_type), allocatable, intent(in) :: var

Write fields that have been updated by assimilation to the CAM 'initial' file.

file_name   Name of CAM initial file to which CAM state 'var' will be written.
model_time   Time of state to be written.
var   Structure containing all the fields of the state vector.

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

Given a model state, returns the value of the variable encoded in obs_kind interpolated to a given location by a method of the model's choosing. Currently DART KINDs: QTY_U_WIND_COMPONENT, QTY_V_WIND_COMPONENT, QTY_SURFACE_PRESSURE, QTY_TEMPERATURE, QTY_SPECIFIC_HUMIDITY, QTY_PRESSURE are supported, but others can be added. QTY_PRESSURE does not have a corresponding field on CAM initial files, but is routinely calculated in CAM using its own subroutines. Interpolation of this field has been incorporated in order to facilitate assimilations of observations which require it, such as GPS radio occultation. CAM-SE requires a completely different interpolation algorithm from the logically rectangular grids of the eulerian and FV dycores. So model_interpolate is now simply calls either interp_lonlat or interp_cubed_sphere.

If the interpolation is valid, istatus = 0. In the case where the observational operator is not defined at the given location (e.g. the observation is below the lowest model level or above the top level), interp_val is returned as MISSING_R8 and istatus > 0. CAM is highly damped in the upper levels of the model, which may require the exclusion of otherwise valid observations above a certain level, which can be specified in the model_mod namelist variable highest_obs_pressure_Pa . Such cases return istatus = 2, and also do the interpolation and return the value, which is NOT used by filter. Eventually such quality control may be moved to another module, but for now is performed in subroutine get_val_pressure. DART reserves the use of istatus < 0.

st_vec    Model state vector.
location Location to which to interpolate.
obs_kind Integer indexing which type of observation is to be interpolated.
interp_val Value interpolated to location.
istatus Integer flag indicating the success of the interpolation.

call solve_quadratic(a, b, c, r1, r2)
real(r8), intent(in)  :: a, b, c
real(r8), intent(out) :: r1, r2

Calculate the roots (r1, r2) of the quadratic equation defined by coefficients a, b, and c.


call interp_lonlat(st_vec, location, obs_kind, interp_val, istatus)
real(r8),            intent(in) :: st_vec(:)
type(location_type), intent(in) :: location
integer,             intent(in) :: obs_kind
real(r8),           intent(out) :: interp_val
integer,            intent(out) :: istatus

Find the 8 corners of the lon-lat-lev grid cell, which encloses an observation at 'location', and interpolate the values of variable 'obs_kind' to that location. Variables defined as in model_interpolate.


call prog_var_to_vector(var,x)
type(model_type), allocatable, intent(in)        :: var
real(r8), allocatable, dimension(:), intent(out) :: x

Insert CAM fields (N-Dimensional) into DART state vector array (1-D)

var   Structure containing all the fields of the state vector.
x   DART state vector.

call vector_to_prog_var(x,var)
real(r8), allocatable, dimension(:), intent(in) :: x
type(model_type), allocatable, intent(out)      :: var

Extract CAM fields (N-Dimensional) from state vector (1-D).

x   DART state vector.
var   Structure containing all the fields of the state vector.

call get_close_obs(filt_gc, base_obs_loc, base_obs_kind, obs, obs_kind, num_close, close_ind [, dist])
type(get_close_type), intent(in)  :: filt_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(:)

First calls the location_mod:get_close_obs() and then updates the distances to account for the extra diffusion in CAM's highest levels.

CAM uses a hybrid vertical coordinate, which requires the surface pressure beneath a point in order to determine the point's vertical location. The ensemble members may each have a different surface pressure, which could cause an observation to have different distances to the surrounding grid points. This can cause the observation to be excluded (distance > 2*cutoff) in some members, but not others, which is an ambiguity we'd like to avoid. The ensemble mean state is used to a provide a consistent vertical location for the whole ensemble. All observations on height and level are handled automatically, using the ensemble mean for calculations.

Due to the extra diffusion at high levels in CAM (see model_interpolate) there is also code which reduces the influence of all observations on model points above some altitude. Currently namelist variable highest_state_pressure_Pa controls this. The default value in model_mod (for CAM), and the values for WACCM given in model_top_issues.pdf, gradually reduce the influence of observations from full strength at highest_state_pressure_Pa to 0 at the model levels where extra diffusion is applied (2 or 3 for CAM, 3 or 4 for WACCM).

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 routine of the same name. See location_mod:get_close_maxdist_init() for the documentation of this subroutine.

gc    A get_close type whose maxdist will be defined here.
maxdist    The maximum distance around a location at which a lon-lat grid box can be considered close to a location.
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 location_mod:get_close_obs_init() for the documentation of this subroutine.

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 ensemble initial conditions perturbed around some control trajectory state when one is preparing to spin-up ensembles. interf_provided is .true. for this cam/model_mod.f90.

This implementation allows each field in the state vector to be randomly perturbed with a separate standard deviation. The fields to be perturbed, and the associated standard deviations are specified in the model_mod namelist using state_names_pert and state_names_sd. The entries in state_names_pert should be in the same order as those in state_names_#d (# = 0,1,2,3 in that order). As in the default, the perturbations sd*N[0,1] are added onto the basic state field, so sd should not be a percentage of the basic state field, but an actual physical size.

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

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

This operation is not defined for the CAM model. This interface is only required if `synchronous' model state advance is supported (the model is called directly as a Fortran90 subroutine from the assimilation programs). This is generally not the preferred method for large models and a stub for this interface is provided for the CAM model.

x State vector of length model_size.
time    Gives time of the initial model state. Needed for models that have real time state requirements, for instance the computation of radiational parameters. Note that DART provides a time_manager_mod module that is used to support time computations throughout the facility.

call end_model( )

Called when use of a model is completed to clean up storage, etc.

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 frequently used to spin-up models from rest, but is not meaningfully supported for the CAM model.

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

Returns default initial conditions for model; generally used for spinning up initial model states. For the CAM model it's just a stub because initial state is always to be provided from input files.

x Model state vector.

[top]

PRIVATE COMPONENTS

Here they are listed alphabetically. Following this list they are listed in the order in which they appear in model_mod.f90, which is separated into sections corresponding to the main functions of model_mod.f90: static_init_model, I/O, model_interpolate, vector<->field translations, get_close_obs, and utility routines.

 convert_vert
 coord_val
 coord_index
 dcz2
 end_grid_1d_instance
 end_model_instance
 find_name
 get_val
 get_val_level
 get_val_height
 get_val_pressure
 gph2gmh
 gravity
 index_from_grid
 create_grid_1d_instance
 init_model_instance
 interp_lonlat
 map_kinds
 model_heights
 nc_read_global_att
 read_cam_coord
 read_cam_init_size
 read_cam_2Dreal
 order_state_fields
 plevs_cam
 set_ps_arrays
 trans_coord
 verify_namelist
 write_cam_coord

call read_cam_init_size(ncfileid)
integer, intent(out)      :: ncfileid

Gets the number of lons, lats and levels from a netcdf CAM initial file and stores in global storage.

ncfileid   NetCDF file ID for CAM the initial file.

call trans_coord(ncfileid)
integer, intent(in) :: ncfileid

Rearranges sizes of the coordinates of each variable as found on caminput.nc file (stored in variable f_dim_#d) into the order model_mod wants to see them for the state (s_dim_#d); lev, lon, lat. # = 1,2,3 for the numbers of dimensions of the variables.

ncfileid   The file ID number of the caminput.nc file.

Uses the model_version from input.nml:model_nml to figure out the order of lev, lat, and lon on the caminput.nc file. Designed to be flexible enough to handle staggered and unstaggered U and V, 2D fields which may be (lev,lat) instead of (lon,lat), and unstructured grids which may come from HOMME or the CLM.

call read_cam_2Dreal(file_name, cfield)
character, len=*,  intent(in)                :: file_name
character, len=8,  intent(in)                :: cfield

Read a 2-D field from a NetCDF file, so far just surface height (phis) from cam_phis.nc.

file_name   The name of the file file containing field "cfield".
cfield   Name of the field (on the file) to be read.

call nc_read_model_atts(att,att_vals,nflds)
character, intent(in)                :: att
character, dimension(:), intent(out) :: att_vals
integer, intent(in)                  :: nflds

Reads the value of an attribute for each of the fields in cflds.

att   The name of an attribute of the fields on a CAM initial file.
att_vals   The values which that attribute has for each of the nflds fields needed for the state vector.
nflds   The number of fields in the state vector.

call nc_read_global_att(ncFileID, att, att_val)
integer,           intent(in)  :: ncFileID
character (len=*), intent(in)  :: att
integer,           intent(out) :: att_val

Reads the value of a global attribute. Needed to get grid information from CAM-SE (cubed sphere) initial files.

ncFileID   The NetCDF identifier of the file on which to find the attribute.
att   The name of the global attribute.
att_val   The value of the attribute.

call read_cam_coord(ncfileid,cfield,var)
integer, intent(in)             :: ncfileid
character(len=8), intent(in)    :: cfield
type(grid_1d_type), intent(out) :: var

Reads a coordinate array and metadata from a CAM initial file.

ncfileid   Length of var.
cfield   Name of var, as found on the CAM initial file.
var   A coordinate array from a CAM initial file.

call create_grid_1d_instance(length, num_atts, var)
integer,            intent(in )   :: length
integer,            intent(in )   :: num_atts
type(grid_1d_type), intent(inout) :: var

Allocate space for the variable of defined type grid_1d_type. Put vector length and the number of attributes in structure. These are coordinate variables read from caminput.nc, and have lots of characteristics which are helpful to keep together in one structure. The rest of var is filled in read_cam_coord.

length   The length of var.
num_atts   The number of attributes of this coordinate, from caminput.nc.
var   The coordinate variable to read from caminput.nc.
  type grid_1d_type
     private
     character (len=8)       :: label          ! e.g. 'lat     '
     integer                 :: dim_id         ! NetCDF dimension ID from the file
     integer                 :: length         ! number of elements in the coordinate array
     real(r8)                :: resolution     ! spacing between elements, or 0. for irreg.
     real(r8), allocatable   :: vals(:)        ! coordinate values
     integer                 :: num_atts       ! number of NetCDF attributes
     character (len=NF90_MAX_NAME), allocatable :: atts_names(:)  ! names of those attributes
     character (len=NF90_MAX_NAME), allocatable :: atts_vals(:)   ! values of those attributes
  end type grid_1d_type

call end_grid_1d_instance(var)
type(grid_1d_type), intent(inout) :: var

Deallocate the array components of the variable of defined type grid_1d_type.

var   Coordinate variable and metadata.

call order_state_fields()

Fills cflds with state_names for use in I/O of caminput.nc Also assigns field TYPE_ variables for use by get_state_meta_data, and other routines.


call map_kinds()

Makes an array (dart_to_cam_types) of 'locations within the state vector' of all the available obs KINDs that come from obs_kind_mod. Also maps the local model_mod TYPEs onto the DART KINDs by the same mechanism.

The obs kind that's needed will be the index into this array, the corresponding value will be the position of that field (not individual variable) within the state vector according to state_name_Xd. This subroutine is called from static_init_model, and arrays dart_to_cam_types and cam_to_dart_kinds are global, so they will not have to be recomputed for every obs.

call verify_namelist()

Checks the values of several of the variables read from input.nml:model_nml.

call write_cam_coord_def(ncFileID, c_name, coord, dim_id, c_id)
integer,            intent(in)  :: ncFileID
character (len=8),  intent(in)  :: c_name
type(grid_1d_type), intent(in)  :: coord
integer,            intent(in)  :: dim_id
integer,            intent(out) :: c_id

Define a coordinate variable on the P[oste]rior_Diag.nc file and write some of its attributes there.

ncFileID   Integer ID number of the NetCDF file.
c_name   Coordinate name.
coord   Coordinate values and metadata.
dim_id   The dimension ID of the coordinate on the NetCDF file.
c_id   Coordinate variable ID number on the NetCDF file.

call get_val_level(st_vec,lon_index,lat_index,level_index,obs_kind,val,istatus)
real(r8), intent(in), allocatable, dimension(:)  :: st_vec
integer, intent(in)                             :: lon_index
integer, intent(in)                             :: lat_index
integer,  intent(in)                             :: level_index
integer, intent(in)                             :: obs_kind
real(r8), intent(out)                           :: val
integer, intent(out)                            :: istatus

Gets the value on the model level for variable obs_kind at the lon_index, lat_index horizontal grid point. Some observations may be excluded here by model_nml restrictions. For examples, observations below the lowest model level are excluded, including surface pressure. So are observations above the lower of highest_obs_pressure_Pa and the model's highest level.

st_vec DART state vector.
lon_index   Index of longitude of this val.
lat_index   Index of latitude of this val.
level_index Index of the level for which this obs_kind must be retrieved.
obs_kind The DART KIND of variable being interpolated, e.g. QTY_TEMPERATURE.
val The value of the obs_kind variable interpolated from state vector x.
istatus Flag containing the status of the interpolation. 0; all good, 1; can't interpolate, 2; can interpolate but don't use the value, < 0; reserved for DART library use.

call get_val_pressure(st_vec,lon_index,lat_index,pressure,obs_kind,val,istatus)
real(r8), intent(in), allocatable, dimension(:) :: st_vec
integer, intent(in)                             :: lon_index
integer, intent(in)                             :: lat_index
real(r8), intent(in)                            :: pressure
integer, intent(in)                             :: obs_kind
real(r8), intent(out)                           :: val
integer, intent(out)                           :: istatus

Gets the vertically interpolated value on pressure for variable obs_kind at the lon_index, lat_index horizontal grid point. Some quality control of observations is done here, but may be moved later. For now, observations above the lower of highest_obs_pressure_Pa and the model's highest level, and below the lowest model level are excluded.

val    The value of the obs_kind variable interpolated from state vector.
st_vec   DART state vector.
lon_index    Index of longitude of this val.
lat_index    Index of latitude of this val.
pressure    The pressure to which this obs_kind must be interpolated.
obs_kind    The kind of variable being interpolated, i.e. T.
istatus    Flag containing the status of the interpolation. 0; all good, 1; can't interpolate, 2; can interpolate but don't use the value, < 0; reserved for DART library use.

call get_val_height(st_vec,lon_index,lat_index,height,obs_kind,val,istatus)
real(r8), intent(in), allocatable, dimension(:) :: st_vec
integer, intent(in)                             :: lon_index
integer, intent(in)                             :: lat_index
real(r8), intent(in)                            :: height
integer, intent(in)                             :: obs_kind
real(r8), intent(out)                           :: val
integer, intent(out)                            :: istatus

Gets the vertically interpolated value on height for variable obs_kind at the lon_index, lat_index horizontal grid point. Some quality control of observations is done here, but may be moved later. For now, observations above the lower of highest_obs_height_m and the model's highest level, and below the lowest model level are excluded.

val   The value of the obs_kind variable interpolated from state vector st_vec.
st_vec   DART state vector.
lon_index   Index of longitude of this val.
lat_index   Index of latitude of this val.
height   The height to which this obs_kind must be interpolated.
obs_kind   The kind of variable being interpolated, i.e. T.
istatus   Flag containing the status of the interpolation. 0; all good, 1; can't interpolate, 2; can interpolate but don't use the value, < 0; reserved for DART library use.

call get_val(st_vec,lon_index,lat_index,level,obs_kind,val,istatus)
real(r8), intent(in), allocatable, dimension(:) :: st_vec
integer, intent(in)                             :: lon_index
integer, intent(in)                             :: lat_index
integer, intent(in)                             :: level
integer, intent(in)                             :: obs_kind
real(r8), intent(out)                           :: val
integer, intent(out)                            :: istatus

Extracts the value of a field at a specified location from the DART state vector.

st_vec    DART state vector.
lon_index    Index of longitude of this val.
lat_index    Index of latitude of this val.
level    The level of the obs_kind variable desired.
obs_kind    The kind of variable being interpolated, i.e. T.
val    The value of the obs_kind variable at a grid point, from state vector x.
istatus    Flag containing the status of the interpolation. 0; all good, 1; can't interpolate, 2; can interpolate but don't use the value, < 0; reserved for DART library use.

call convert_vert(old_array,old_which,old_loc,old_kind,new_array,new_which)
real(r8), dimension(3), intent(in)    :: old_array
integer,                intent(in)    :: old_which
type(location_type),    intent(in)    :: old_loc
integer,                intent(in)    :: old_kind
real(r8), dimension(3), intent(inout) :: new_array
integer,                intent(out)   :: new_which

Uses model information and subroutines to convert the vertical location of an observation (prior, model state variable, or actual ob) into the standard vertical coordinate (pressure or scale height). Accommodates staggered grid US and VS of the FV core.

old_array   The location array of the incoming ob.
old_which   The vertical coordinate type of the incoming ob.
old_loc   The DART location of the incoming ob.
old_kind   The DART KIND of the incoming ob.
new_array   The location array of the converted ob.
new_which   The vertical coordinate type of the converted ob.
call init_closest_node()

Initialize the "get_close" structure for the cubed sphere grid. The structure is used to find a small set of candidate "close" nodes to a location.

closest = find_closest_node(lat, lon)
real(r8), intent(in)  :: lat, lon

Determine the index or name of the closest node (corner) to the given point.

closest    The name of the cubed sphere node (cell corner) closest to (lon,lat).
lat, lon    The latitude and longitude (degrees) of the point of interest.

index = index_from_grid(lev_ind, lon_ind, lat_ind, ifld)
integer, intent(in) :: lev_ind
integer, intent(in) :: lon_ind
integer, intent(in) :: lat_ind
integer, intent(in) :: ifld
integer             :: index_from_grid

Function to generate the state vector index corresponding to the grid location and variable given.

lev_ind   Level of the desired variable.
lon_ind   Longitude of the desired variable.
lat_ind   Latitude of the desired variable.
ifld   TYPE_ (PS or T or ...) index of the variable whose index is needed.

Called by many routines.

name = find_name(nam, list)
character (len=*),              intent(in) :: nam
character (len=*), dimension(:),intent(in) :: list

Function to return the index of a character string as found within a list of character strings, typically variable or dimension names.

nam   The string to be found within list.
list   The list which hopefully contains nam.

call coord_val(dim_name, indx, lon_val, lat_val, lev_val)
character (len=*), intent(in)    :: dim_name
integer,           intent(in)    :: indx
real(r8),          intent(inout) :: lon_val
real(r8),          intent(inout) :: lat_val
real(r8),          intent(inout) :: lev_val

Given the name of the coordinate to be searched and the index into that array, returns the coordinate value in either lon_val, lat_val, or lev_val.

dim_name   The name of the dimension being searched for the value at an index.
indx   The index of the variable whose location is needed.
lon_val   (Possibly) the location of the grid point in the longitudinal direction.
lat_val   (Possibly) the location of the grid point in the latitudinal direction.
lev_val   (Possibly) the location of the grid point in the vertical direction.

All 3 _val arguments are present so that this routine can return the value in the coordinate that the calling routine wants it to be, and search/placement doesn't have to be done there. Used by get_state_meta_data and model_interpolate.

call coord_index(dim_name, val, indx , other_indx)
character (len=*), intent(in)  :: dim_name
real(r8),          intent(in)  :: val
integer,           intent(out) :: indx
integer, optional, intent(out) :: other_indx

Given the name of the coordinate to be searched and the coordinate value, Returns the index of the closest coordinate value. Optionally returns the next closest index too, which may be < or > the closest.

dim_name   The name of the dimension being searched for the value at an index.
val   The coordinate value whose nearest neighbors are needed.
indx   The index whose value is closest to val.
other_indx   The index whose value is 2nd closest to val.

Used by get_state_meta_data. Uses coordinate metadata to search differently, depending on whether coordinate is regularly spaced or not.

call set_ps_arrays(vec)
real(r8), intent(in) :: vec(:)

Subroutine to put the ensemble average PS into a globally defined array. Also provides PS on the FV staggered grids, if needed.

vec   The ensemble average state vector.

call plevs_cam(p_surf, num_levs, pmid)
real(r8), dimension(:), intent(in)  :: p_surf
integer, intent(in)                             :: num_levs
real(r8), dimension(:), intent(out) :: pmid

Define the pressures of the CAM data levels (layer midpoints) from the coordinate definitions and the surface pressure.

p_surf   Surface pressure at this latitude and longitude (Pa).
num_levs   Number of levels of pressures to calculate
pmid   Pressures at the CAM "midpoint" levels.

call model_heights(num_levs,vec,p_surf,base_obs_loc,model_h,istatus)
integer,                intent(in)  :: num_levs
real(r8), dimension(:), intent(in)  :: vec
real(r8),               intent(in)  :: p_surf
type(location_type),    intent(in)  :: base_obs_loc
real(r8), dimension(:), intent(out) :: model_h
integer,                intent(out) :: istatus

This routine calculates geometrical height (m) at mid-layers of the CAM model

num_levs    Number of vertical levels.
vec   DART state vector.
p_surf   Surface pressure at this location.
base_obs_loc   The (horizontal) location at which heights will be calculated.
model_h   geometrical height at midlayer (m).
istatus   Status of the calculation returned to the calling routine.

call dcz2(kmax,p_surf,phis0,tv,hprb,hybrid_As,hybrid_Bs,pmln,pterm,z2)
integer,                       intent(in)  :: kmax
real(r8),                      intent(in)  :: p_surf
real(r8),                      intent(in)  :: phis0
real(r8), dimension(kmax),     intent(in)  :: tv
real(r8),                      intent(in)  :: hprb
real(r8), dimension(2,kmax+1), intent(in)  :: hybrid_As
real(r8), dimension(2,kmax+1), intent(in)  :: hybrid_Bs
real(r8), dimension(kmax+1),   intent(out) :: pmln
real(r8), dimension(kmax),     intent(out) :: pterm
real(r8), dimension(kmax),     intent(out) :: z2
Compute geopotential height for a CESM hybrid coordinate column.
All arrays except hybrid_As, hybrid_Bs are oriented top to bottom.
hybrid_[AB]s first subscript:
  = 1 for layer interfaces
  = 2 for layer midpoints
hybrid_As coord coeffs for P0 reference pressure term in plevs_cam
hybrid_Bs coord coeffs for surface pressure term in plevs_cam (in same format as hybrid_As)
kmax   Number of vertical levels.
p_surf   Surface pressure (pascals).
phis0   Surface geoptential.
tv   Virtual temperature, top to bottom.
hprb   Hybrid base pressure.
hybrid_As   Hybrid A coefficients (multiply P0).
hybrid_Bs   Hybrid B coefficients (multiply surface pressure).
pmln   Natural logs of midpoint pressures.
pterm   Pressure profile.
z2   Geopotential height, top to bottom.

mean_sea_level_alt = gph2gmh(h, lat)
real(r8), intent(out)                           :: gph2gmh
real(r8), intent(in)                            :: h
real(r8), intent(in)                            :: lat

Convert a geopotential altitude to mean sea level altitude.

gph2gmh   MSL altitude, in km.
h   geopotential altitude (in km).
lat   latitude in degrees.

call gravity (xlat,alt,galt)
real(r8), intent(in)                            :: xlat
real(r8), intent(in)                            :: alt
real(r8), intent(out)                           :: galt

This subroutine computes the Earth's gravity at any altitude and latitude (Bill Schreiner, May/1995). The model assumes the Earth is an oblate spheriod rotating at a the Earth's spin rate. The model was taken from "Geophysical Geodesy, Kurt Lambeck, 1988".

xlat   latitude in radians.
alt   altitude above the reference ellipsiod, km.
galt   gravity at the given lat and alt, km/sec*sec.

call init_model_instance(var)
type(model_type), allocatable, intent(out) :: var

Initializes an instance of a cam model state variable; all the fields specified for the state vector.

var   Structure which contains all the fields of various dimensions which are in the DART state vector.
type model_type
    private
    real(r8), allocatable :: vars_0d(:)           ! scalars
    real(r8), allocatable :: vars_1d(:, :)        ! vectors
    real(r8), allocatable :: vars_2d(:, :, :)     ! 2-D fields
    real(r8), allocatable :: vars_3d(:, :, :, :)  ! 3-D fields
end type model_type
call end_model_instance(var)
type(model_type), allocatable, intent(inout) :: var

Ends an instance of a cam model_type state variable

var   Structure which contains all the fields included in the state vector.

[top]

DISCUSSION

Many CAM initial file variables are already handled in the model_mod. Here is a list of others, which may be used in the future. Each would need to have a DART KIND associated with it in model_mod.

  Atmos
     CLOUD:       "Cloud fraction" ;
     QCWAT:       "q associated with cloud water" ;
     TCWAT:       "T associated with cloud water" ;
     CWAT:        "Total Grid box averaged Condensate Amount (liquid + ice)" ;
    also? LCWAT

  pbl
     PBLH:        "PBL height" ;
     QPERT:       "Perturbation specific humidity (eddies in PBL)" ;
     TPERT:       "Perturbation temperature (eddies in PBL)" ;

  Surface
     LANDFRAC:    "Fraction of sfc area covered by land" ;
     LANDM:       "Land ocean transition mask: ocean (0), continent (1), transition (0-1)" ;
       also LANDM_COSLAT
     ICEFRAC:     "Fraction of sfc area covered by sea-ice" ;
     SGH:         "Standard deviation of orography" ;
     Z0FAC:       "factor relating z0 to sdv of orography" ;
     TS:          "Surface temperature (radiative)" ;
     TSOCN:       "Ocean tempertare" ;
     TSICE:       "Ice temperature" ;
     TSICERAD:    "Radiatively equivalent ice temperature" ;

  Land/under surface
     SICTHK:      "Sea ice thickness" ;
     SNOWHICE:    "Water equivalent snow depth" ;
     TS1:         "subsoil temperature" ;
     TS2:         "subsoil temperature" ;
     TS3:         "subsoil temperature" ;
     TS4:         "subsoil temperature" ;

  Other fields are not included because they look more CLM oriented.

  Other fields which users may add to the CAM initial files are not listed here.
  

[top]

FILES

[top]

REFERENCES

Ave Arellano did the first work with CAM-Chem, assimilating MOPPITT CO observations into CAM-Chem. Jerome Barre and Benjamin Gaubert took up the development work from Ave, and prompted several additions to DART, as well as cam/model_mod.

Nick Pedatella developed the first vert_coord = 'log_invP' capability to enable assimilation using WACCM and scale height vertical locations.

[top]

ERROR CODES and CONDITIONS

RoutineMessageComment
read_cam_init_size
read_cam_init
write_cam_init
nc_read_model_atts
nc_write_model_atts
nc_write_model_vars
read_cam_horiz
read_cam_coord
Various NetCDF-f90 interface error messages From one of the NetCDF calls in the named routine
prog_var_to_vector
vector_to_prog_var
get_state_meta_data
get_val
model_get_close_states
scalar and 1-D vector components of state vector are not coded into this routine Only 2D and 3D fields can be part of state vector so far.
prog_var_to_vector
vector_to_prog_var
indx # and model_size # must be equal indx was tallied during insertion of fields into state vector
model_get_close_states which_vert = # not handled in model_get_close_states See which_vert description in location/threed_sphere/location_mod.html
nc_write_model_atts Time dimension ID # must match Unlimited Dimension ID # NetCDF file writing error
order_state_fields nfld = #, nflds = # must be equal Mismatch(es) of state_names_#d and state_num_#d in model_nml

[top]

KNOWN BUGS

The nonlocal GPS observation operator requires that highest_obs_pressure_Pa be ≤ 10000 Pa (100 mb). The local GPS observation operator doesn't have this restriction.

GPS observations require that highest_obs_pressure_Pa be <= 10000 Pa (100 mb).

[top]

FUTURE PLANS

Nitty gritty: Efficiency possibilities

[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: Kevin Raeder
Revision: $Revision: 11427 $
Source: $URL: https://svn-dares-dart.cgd.ucar.edu/DART/releases/Manhattan/models/cam-fv/model_mod.html $
Change Date: $Date: 2017-03-31 15:24:30 -0600 (Fri, 31 Mar 2017) $
Change history:  try "svn log" or "svn diff"