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

MODULE location_mod (threed_sphere)

Contact: Jeff Anderson, Tim Hoar
Revision: $Revision: 2944 $
Source: $URL: http://subversion.ucar.edu/DAReS/DART/trunk/location/threed_sphere/location_mod.html $
Change Date: $Date: 2007-05-23 15:05:17 -0600 (Wed, 23 May 2007) $
Change history: try "svn log" or "svn diff"

OVERVIEW

Optional arguments are enclosed in brackets [like this].

Provides a representation of a physical location on a 3-D spherical shell. A type that abstracts the location is provided along with operators to compute the distance between two locations. This is a member of a class of similar location modules that provide the same abstraction for different represenations of physical space (for instance a one-dimensional periodic domain for low-order modules).




OTHER MODULES USED

types_mod
utilities_mod
random_seq_mod



PUBLIC INTERFACE

use location_mod, only : location_type
 get_dist
 get_location
 set_location
 write_location
 read_location
 interactive_location
 set_location_missing
 vert_is_undef
 vert_is_surface
 vert_is_pressure
 vert_is_level
 vert_is_height
 query_location
 LocationDims
 LocationName
 LocationLName
 horiz_dist_only
 get_close_obs
 get_close_type
 get_close_maxdist_init
 get_close_obs_init
 VERTISUNDEF
 VERTISSURFACE
 VERTISLEVEL
 VERTISPRESSURE
 VERTISHEIGHT
 operator(==)
 operator(/=)

NOTES

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




PUBLIC ENTITIES




type location_type
   private
   real(r8) :: lon, lat, vloc
   integer  :: which_vert
end type location_type

Detailed Description

Provides an abstract representation of physical location on a three-d spherical shell.

Component Description
lon longitude in radians
lat latitude in radians
vloc vertical location, units as selected by which_vert
which_vert type of vertical location: -2=no vert location; -1=surface; 1=level; 2=pressure; 3=height


var = get_location(loc)
 real(r8), dimension(3), intent(inout) :: get_location
 type(location_type), intent(in)       :: loc
 

Description

Extracts the longitude and latitude (converted to degrees) and the vertical location from a location type and returns in a 3 element real array.

get_location    The longitude and latitude (in degrees) and vertical location
loc    A location type


var = set_location(lon, lat, vert_loc, which_vert)
 type(location_type)     :: set_location
 real(r8), intent(in)    :: lon
 real(r8), intent(in)    :: lat
 real(r8), intent(in)    :: vert_loc
 integer, intent(in)     :: which_vert
 

Description

Returns a location type with the input longitude and latitude (input in degrees) and the vertical location of type specified by which_vert.

set_location    A location type
lon    Longitude in degrees
lat    Latitude in degrees
vert_loc    Vertical location consistent with which_vert
which_vert    The vertical location type


call write_location(file, loc [,fform])
integer,               intent(in)      ::  file 
type(location_type),   intent(in)      ::  loc 
character(len=*), optional, intent(in) ::  fform 

Description

Given an integer IO channel of an open file and a location, writes the location to this file. The fform argument controls whether write is "FORMATTED" or "UNFORMATTED" with default being formatted.

file the unit number of the open file.
loc location type to be written.
fform Format specifier (FORMATTED or UNFORMATTED).

Notes

Eventually, a more general file descriptor type should replace the use of the integer unit number.


var = get_dist( loc1, loc2, no_vert)
real(r8)                             ::  get_dist 
type(location_type),     intent(in)  ::  loc1, loc2 
logical, optional                    ::  no_vert 

Description

Returns the distance between two locations. If horiz_only is set to true just computes great circle distance on sphere. If horiz_only is false, then computes an ellipsoidal distance with the horizontal component as above and the vertical distance determined by the types of the locations and the normalization constants set by the namelist for the different vertical coordinate types. The vertical normalization gives the vertical distance that is equally weighted as a horizontal distance of 1 radian. If no_vert is present, it overrides the value in the namelist and controls whether vertical distance is included or not.

loc1     first location.
loc2 second location.
no_vert If true, no vertical component to distance. If false, vertical component is included.
var distance between loc1 and loc2.


var = read_location(ifile [,fform])
 type(location_type)                    :: read_location
 integer, intent(in)                    :: ifile
 character(len=*), optional, intent(in) :: fform
 

Description

Reads a location_type from a file open on channel ifile using format fform (default is formatted).

read_location    Returned location type read from file
ifile    Integer channel opened to a file to be read
fform    Optional format specifier (FORMATTED or UNFORMATTED)


call interactive_location(location [,set_to_default])
 type(location_type), intent(out) :: location
 logical, optional, intent(in)    :: set_to_default
 

Description

Use standard input to define a location type. With set_to_default true get one with all elements set to 0.

location    Location created from standard input
set_to_default    If true, sets all elments of location type to 0


var = set_location_missing()
 type(location_type), intent(in) :: set_location_missing
 

Description

Returns a location with all elements set to missing values defined in types module.

set_location_missing    A location with all elements set to missing values


var = vert_is_undef(loc)
 logical                         :: vert_is_undef
 type(location_type), intent(in) :: loc
 

Description

Returns true if which_vert has no vertical location (undefined), else false.

vert_is_undef    Returns true if vertical coordinate is not applicable.
loc    A location type


var = vert_is_surface(loc)
 logical                         :: vert_is_surface
 type(location_type), intent(in) :: loc
 

Description

Returns true if which_vert is for surface, else false.

vert_is_surface    Returns true if vertical coordinate type is surface
loc    A location type


var = vert_is_pressure(loc)
 logical                         :: vert_is_pressure
 type(location_type), intent(in) :: loc
 

Description

Returns true if which_vert is for pressure, else false.

vert_is_pressure    Returns true if vertical coordinate type is pressure
loc    A location type


var = vert_is_level(loc)
 logical                         :: vert_is_level
 type(location_type), intent(in) :: loc
 

Description

Returns true if which_vert is for level, else false.

vert_is_level    Returns true if vertical coordinate type is level
loc    A location type


var = vert_is_height(loc)
 logical                         :: vert_is_height
 type(location_type), intent(in) :: loc
 

Description

Returns true if which_vert is for height, else false.

vert_is_height    Returns true if vertical coordinate type is height
loc    A location type


var = query_location(loc [,attr])
 real(r8), intent(out)                  :: query_location
 type(location_type), intent(in)        :: loc
 character(len=*), optional, intent(in) :: attr
 

Description

Returns the value of which_vert, latitude, longitude, or vertical location from a location type as selected by the string argument attr. If attr is not present or if it is 'WHICH_VERT', the value of which_vert is converted to real and returned. Otherwise, attr='LON' returns longitude, attr='LAT' returns latitude and attr='VLOC' returns the vertical location.

query_location    Returns longitude, latitude, vertical location, or which_vert (converted to real)
loc    A location type
attr    Selects 'WHICH_VERT', 'LON', 'LAT' or 'VLOC'


loc1 == loc2
 type(location_type), intent(in)       :: loc1, loc2
 

Description

Returns true if the two location types have identical values, else false.


loc1 /= loc2
 type(location_type), intent(in)       :: loc1, loc2
 

Description

Returns true if the two location types do NOT have identical values, else false.



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

Description

Initialize storage for efficient identification of locations close to a given location. Allocates storage for keeping track of which 'box' each observation in the list is in.

gc    Structure that contains data to efficiently find locations close to a given location.
num    The number of locations in the list.
obs    The locations of each element in the list, not used in 1D implementation.


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), dimension(:), intent(in) :: obs
 integer, dimension(:), intent(in)             :: obs_kind
 integer, intent(out)                          :: num_close
 integer, dimension(:), intent(out)            :: close_ind
 real(r8), dimension(:), intent(out)           :: dist
 

Description

Given a single location and a list of other locations, returns the indices of all the locations close to the single one along with the number of these and the distances for the close ones.

gc    Structure to allow efficient identification of locations close to a given location.
base_obs_loc    Single given location.
base_obs_kind    Kind of the single location.
obs    List of observations from which close ones are to be found.
obs_kind    Kind associated with observations in obs list.
num_close    Number of observations close to the given location.
close_ind    Indices of those locations that are close.
dist    Distance between given location and the close ones identified in close_ind.


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

Description

Sets the threshhold distance. Anything closer that this is deemed to be close.

gc    Data for efficiently finding close locations.
maxdist    Anything closer than this distance is a close location.




type get_close_type
   private
   integer  :: num
   real(r8) :: maxdist
   integer, pointer :: lon_offset(:, :)
   integer, pointer :: obs_box(:)
   integer, pointer :: count(:, :)
   integer, pointer :: start(:, :)
end type get_close_type

Description

Provides a structure for doing efficient computation of close locations. Doesn't do anything in the 1D implementation except provide appropriate stubs.

Component Description
num Number of observations in list
maxdist Threshhold distance. Anything closer is close.
lon_offset Dimensioned nlon by nlat. For a given offset in longitude boxes and difference in latitudes, gives max distance from base box to a point in offset box.
obs_box Dimensioned num. Gives index of what box each observation is in.
count Dimensioned nlon by nlat. Number of obs in each box.
start Dimensioned nlon by nlat. Index in straight storage list where obs in each box start.


Global variable: integer :: VERTISUNDEF

Global variable: integer :: VERTISSURFACE

Global variable: integer :: VERTISLEVEL

Global variable: integer :: VERTISPRESSURE

Global variable: integer :: VERTISHEIGHT

Constant parameters used to denote the different vertical types.


Global variable: integer :: LocationDims

Contains the number of real values in a location type. Useful for output routines that must deal transparently with many different location modules.


Global variable: character(len=129) :: LocationName

A parameter set to "loc3Dsphere" used to identify this location module in netCDF output metadata.


Global variable: character(len=129) :: LocationLName

A parameter set to "threed sphere locations: lon, lat, lev or pressure" used to identify this location module in netCDF output long name metadata.


Global variable: logical :: horiz_dist_only

This is a namelist variable that determines whether distance computation is horizontal only or not.



NAMELIST

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

 namelist / location_nml /  &
   horiz_dist_only, vert_normalization_pressure, vert_normalization_height, &
   vert_normalization_level, approximate_distance, nlon, nlat, output_box_info
 

Discussion

Items in this namelist either control the way in which distances are computed and/or influence the code performance.

If horiz_distance_only is .false. then the appropriate normalization constant determines the relative impact of vertical and horizontal separation. Since only a single localization distance is specified, and the vertical scales might have very different distance characteristics, the vert_normalization_xxx values can be used to scale the vertical appropriately to control the desired influence of observations in the vertical.

For search efficiency all locations are pre-binned. The surface of the sphere is divided up into nlon by nlat boxes and the index numbers of all items (both state vector entries and observations) are stored in the appropriate box. To locate all points close to a given location, only the locations listed in the boxes within the search radius must be checked.

The default values have given good performance on many of our existing model runs, but for tuning purposes the box counts have been added to the namelist to allow adjustment. By default the code prints some summary information about how full the average box is, how many are empty, and how many items were in the box with the largest count. The namelist value output_box_info can be set to .true. to get even more information about the box statistics. The best performance will be obtained somewhere between two extremes; the worst extreme is all the points are located in just a few boxes. This degenerates into a (slow) linear search through the index list. The other extreme is a large number of empty or sparsely filled boxes. The overhead of creating, managing, and searching a long list of boxes will impact performance. The best performance lies somewhere in the middle, where each box contains a reasonable number of values, more or less evenly distributed across boxes. The absolute numbers for best performance will certainly vary from case to case.

For latitude, the nlat boxes are distributed evenly across the actual extents of the data. (Locations are in radians, so the maximum limits are the poles at -PI/2 and +PI/2). For longitude, the code automatically determines if the data is spread around more than half the sphere, and if so, the boxes are distributed evenly across the entire sphere (longitude range 0 to 2*PI). If the data spans less than half the sphere in longitude, the actual extent of the data is determined (including correctly handling the cyclic boundary at 0) and the boxes are distributed only within the data extent. This simplifies the actual distance calculations since the distance from the minimum longitude box to the maximum latitude box cannot be shorter going the other way around the sphere. In practice, for a global model the boxes are evenly distributed across the entire surface of the sphere. For local or regional models, the boxes are distributed only across the the extent of the local grid.

This namelist is read in a file called input.nml .

Contents Type Description
horiz_dist_only logical Only compute horizontal distance component? Default: true
vert_normalization_pressure real(r8) How many pascals should be equivalent distance to one radian? Default: 100000.0
vert_normalization_height real(r8) How many geopotential meters should be equal to one radian? Default: 10000.0
vert_normalization_level real(r8) How many model levels should be equal to one radian? Default: 20.0
approximate_distance logical If true, uses a table lookup for fast approximate computation of distances on sphere. Distance computation can be a first order cost for some spherical problems, so this can increase speed significantly at a loss of some precision. Default: .false.
nlon integer Number of boxes (bins) created in the longitude direction. Must be an odd number. (See discussion above for more information about this item.) Default: 71
nlat integer Number of boxes (bins) created in the latitude direction. (See discussion above for more information about this item.) Default: 36
output_box_info logical If true, print more details about the distribution of locations across the array of boxes. Default: .false.


FILES


REFERENCES


ERROR CODES and CONDITIONS

RoutineMessageComment
initialize_module nlon must be odd Tuning parameter for number of longitude boxes must be odd for algorithm to function.
get_dist Dont know how to compute vertical distance for unlike vertical coordinates Need same which_vert for distances.
set_location longitude (#) is not within range [0,360] Is it really a longitude?
set_location latitude (#) is not within range [-90,90] Is it really a latitude?
set_location which_vert (#) must be one of 1, 2, or 3 Vertical coordinate type restricted to:
-1 = surface value
1 = level
2 = pressure
3 = height
read_location Expected location header "loc3d" in input file, got ___ Vertical coordinate confusion involving NetCDF file.
nc_write_location Various NetCDF-f90 interface error messages From one of the NetCDF calls in nc_write_location


KNOWN BUGS

The Hawaii and Workshop versions of this module had an error in the approximate distance computation. The available values in the lookup table for cosine were insufficient for some cases. This manifested itself as potential errors, most commonly for computing distances near the poles. For relatively small horizontal localizations, this problem only occurred for locations very near the pole. None.


FUTURE PLANS

Need to provide more efficient algorithms for getting close observations and document the nlon and nlat choices and their impact on cost.