MODULE obs_def_1d_state_mod

DART project logo

Jump to DART Documentation Main Index
version information for this file:
$Id: obs_def_1d_state_mod.html 6380 2013-08-05 23:47:11Z nancy $

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

Overview

There is a flexible mechanism built into the DART framework for defining, at compile time, the list of observation types to be supported by the executables. This can be changed at any time by adding or removing items from a namelist and rebuilding the programs. The special obs_def module being described here is used by the program preprocess to insert appropriate code sections into a DEFAULT_obs_def module and a DEFAULT_obs_kind module to generate obs_def_mod.f90 and obs_kind_mod.f90 that are used by filter and other DART programs. This module is intended to provide a prototype for those developing more complicated specialized observation definition modules.

This is an extended format Fortran 90 module that provides the definition for observation types designed for use with idealized low-order models that use the 1D location module and can be thought of as having a state vector that is equally spaced on a 1D cyclic domain. The two observation types are: a straight linear interpolation to a point on a [0,1] domain called a RAW_STATE_VARIABLE; and an area-weighted 'integral' of the state variable over some part of the cyclic 1D domain. The second type can be used to do idealized studies related to remote sensing observations that are best thought of as weighted integrals of some quantity over a finite volume.

The second observation type is referred to as a RAW_STATE_1D_INTEGRAL. It has an associated half_width and localization type (see the cov_cutoff module documentation) and a number of points at which to compute the associated integral by quadrature. The location of the observation defines the center of mass of the integral. The integral is centered around the location and extends outward on each side to 2*half_width. The weight associated with the integral is defined by the weight of the localization function (for instance Gaspari Cohn) using the same localization options as defined by the cov_cutoff module. The number of points are used to equally divide the range for computing the integral by quadrature.

This module was partly motivated to provide a prototype for those developing more complicated specialized observation definition modules. Special observation modules like this contain Fortran 90 code PLUS additional specially formatted commented code that is used to guide the preprocess program in constructing the obs_def_mod.f90 and obs_kind_mod.f90. These specially formatted comments are most conveniently placed at the beginning of the module and comprise seven sections, each beginning and ending with a special F90 comment line that must be included VERBATIM.

The seven sections and their specific instances for the 1d_raw_state_mod are:

  1. A list of all observation types defined by this module and their associated generic kinds (see obs_kind module). The header line is followed by lines that have the observation type name (an all caps Fortran 90 identifier) and their associated generic kind identifier from the obs_kind module. If there is no special processing needed for an observation type and no additional data needed beyond the standard contents of an observation then a third word on the line, COMMON_CODE, will instruct the preprocess program to automatically generate all stubs and code needed for this type. For observation types needing special code or additional data, this word should not be specified and the user must supply the code manually.
    ! BEGIN DART PREPROCESS KIND LIST 
    ! RAW_STATE_VARIABLE,    KIND_RAW_STATE_VARIABLE,   COMMON_CODE
    ! RAW_STATE_1D_INTEGRAL, KIND_1D_INTEGRAL 
    ! END DART PREPROCESS KIND LIST 
    


  2. A list of all the use statements that the completed obs_def_mod.f90 must have in order to use the public interfaces provided by this special obs_def module. This section is optional if there are no external interfaces.
    ! BEGIN DART PREPROCESS USE OF SPECIAL OBS_DEF MODULE 
    !   use obs_def_1d_state_mod, only : write_1d_integral, read_1d_integral,  &
    !                                    interactive_1d_integral, get_expected_1d_integral, &
    !                                    set_1d_integral
    ! END DART PREPROCESS USE OF SPECIAL OBS_DEF MODULE 
    


  3. Case statement entries for each observation type defined by this special obs_def module stating how to compute the forward observation operator. Not permitted if COMMON_CODE was specified (because the case statement will be automatically generated), otherwise there must be a case statement entry for each type of observation.
    ! BEGIN DART PREPROCESS GET_EXPECTED_OBS_FROM_DEF
    !         case(RAW_STATE_1D_INTEGRAL) 
    !            call get_expected_1d_integral(state, location, obs_def%key, obs_val, istatus) 
    ! END DART PREPROCESS GET_EXPECTED_OBS_FROM_DEF 
    


  4. Case statement entries for each observation type defined by this special obs_def module stating how to read any extra required information from an obs sequence file. Not permitted if COMMON_CODE was specified (because the case statement will be automatically generated), otherwise there must be a case statement entry for each type of observation defined even if no special action is required. (In that case put a 'continue' statement as the body of the case instead of a subroutine call.)
    ! BEGIN DART PREPROCESS READ_OBS_DEF 
    !      case(RAW_STATE_1D_INTEGRAL) 
    !         call read_1d_integral(obs_def%key, ifile, fform) 
    ! END DART PREPROCESS READ_OBS_DEF 
    


  5. Case statement entries for each observation type defined by this special obs_def module stating how to write any extra required information to an obs sequence file. Not permitted if COMMON_CODE was specified (because the case statement will be automatically generated), otherwise there must be a case statement entry for each type of observation defined even if no special action is required. (In that case put a 'continue' statement as the body of the case instead of a subroutine call.)
    ! BEGIN DART PREPROCESS WRITE_OBS_DEF 
    !      case(RAW_STATE_1D_INTEGRAL) 
    !         call write_1d_integral(obs_def%key, ifile, fform) 
    ! END DART PREPROCESS WRITE_OBS_DEF 
    


  6. Case statement entries for each observation type defined by this special obs_def module stating how to interactively create any extra required information. Not permitted if COMMON_CODE was specified (because the case statement will be automatically generated), otherwise there must be a case statement entry for each type of observation defined even if no special action is required. (In that case put a 'continue' statement as the body of the case instead of a subroutine call.)
    ! BEGIN DART PREPROCESS INTERACTIVE_OBS_DEF 
    !      case(RAW_STATE_1D_INTEGRAL) 
    !         call interactive_1d_integral(obs_def%key) 
    ! END DART PREPROCESS INTERACTIVE_OBS_DEF 
    


  7. Any executable F90 module code must be tagged with the following comments. All lines between these markers will be copied, verbatim, to obs_def_mod.f90. This section is not required if there are no observation-specific subroutines.
    ! BEGIN DART PREPROCESS MODULE CODE
    module obs_def_1d_state_mod
    
    ... (module executable code)
    
    end module obs_def_1d_state_mod
    ! END DART PREPROCESS MODULE CODE
    



[top]

OTHER MODULES USED

types_mod
utilities_mod
location_mod (1d_location_mod_only)
time_manager_mod
assim_model_mod
cov_cutoff_mod
[top]

PUBLIC INTERFACES

use obs_def_mod, only : write_1d_integral
 read_1d_integral
 interactive_1d_integral
 get_expected_1d_integral
 set_1d_integral

call write_1d_integral(igrkey, ifile, fform)
integer,          intent(in) :: igrkey
integer,          intent(in) :: ifile
character(len=*), intent(in) :: fform

Writes out the extra information for observation with unique identifier key for a 1d_integral observation kind. This includes the half-width, localization type and number of quadrature points for this observation.

igrkey   Unique integer key associated with the 1d integral observation being processed. This is not the same as the key that all types of observations have and uniquely distinguishes all observations from each other; this is a key that is only set and retrieved by this code for 1d integral observations. It is stored in the obs_def derived type, not in the main obs_type definition.
ifile   Unit number on which observation sequence file is open
fform   String noting whether file is opened for 'formatted' or 'unformatted' IO.


call read_1d_integral(igrkey, ifile, fform)
integer,          intent(out) :: igrkey
integer,          intent(in)  :: ifile
character(len=*), intent(in)  :: fform

Reads the extra information for observation with unique identifier key for a 1d_integral observation kind. This information includes the half-width, localization type and number of quadrature points for this observation. The key that is returned is uniquely associated with the definition that has been created and is used by this module to keep track of the associated parameters for this observation.

igrkey   Unique integer key associated with the observation being processed.
ifile   Unit number on which observation sequence file is open
fform   String noting whether file is opened for 'formatted' or 'unformatted' IO.


call interactive_1d_integral(igrkey)
integer, intent(out) :: igrkey

Uses input from standard in to define the characteristics of a 1D integral observation. The key that is returned is uniquely associated with the definition that has been created and can be used by this module to keep track of the associated parameters (half_width, localization kind, number of quadrature points) for this key.

igrkey   Unique identifier associated with the created observation definition in the obs sequence.


call get_expected_1d_integral(state, location, igrkey, val, istatus)
real(r8), intent(in)            :: state
type(location_type), intent(in) :: location
integer, intent(in)             :: igrkey
real(r8), intent(out)           :: val
integer, intent(out)            :: istatus

Computes the forward observation operator for a 1d integral observation. Calls the interpolate() routine multiple times to invoke the forward operator code in whatever model this has been compiled with.

state   Model state vector (or extended state vector).
location   Location of this observation.
igrkey   Unique integer key associated with this observation.
val   Returned value of forward observation operator.
istatus   Returns 0 if forward operator was successfully computed, else returns a positive value. (Negative values are reserved for system use.)


call set_1d_integral(integral_half_width, num_eval_pts, localize_type, igrkey, istatus)
real(r8), intent(in)  :: integral_half_width
integer,  intent(in)  :: num_eval_pts
integer,  intent(in)  :: localize_type
integer,  intent(out) :: igrkey
integer,  intent(out) :: istatus

Available for use by programs that create observations to set the additional metadata for these observation types. This information includes the integral half-width, localization type and number of quadrature points for this observation. The key that is returned is uniquely associated with the definition that has been created and should be set in the obs_def structure by calling set_obs_def_key(). This key is different from the main observation key which all observation types have. This key is unique to this observation type and is used when reading in the observation sequence to match the corresponding metadata with each observation of this type.

integral_half_width   Real value setting the half-width of the integral.
num_eval_pts   Integer, number of evaluation points. 5-20 recommended.
localize_type   Integer localization type: 1=Gaspari-Cohn; 2=Boxcar; 3=Ramped Boxcar
igrkey   Unique integer key associated with the observation being processed.
istatus   Return code. 0 means success, any other value is an error

[top]

NAMELIST

This module has no namelist.

[top]

FILES

[top]

REFERENCES

  1. none
[top]

ERROR CODES and CONDITIONS

RoutineMessageComment
interactive_1d_integral Out of space, max_1d_integral_obs limit NNNN (currently 1000). There is only room for a fixed number of 1d integral observations. The max number is defined by max_1d_integral_obs. Set this to a larger value if more are needed.

KNOWN BUGS

none at this time

[top]

FUTURE PLANS

none at this time

[top]

Terms of Use

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

Contact: DART core group
Revision: $Revision: 6380 $
Source: $URL: https://svn-dares-dart.cgd.ucar.edu/DART/releases/Lanai/obs_def/obs_def_1d_state_mod.html $
Change Date: $Date: 2013-08-05 17:47:11 -0600 (Mon, 05 Aug 2013) $
Change history:  try "svn log" or "svn diff"