Routines for creating Gaussian-shaped land surface topography for latitude-longitude grids.
More...
|
subroutine | gaussian_topog_init_ (lon, lat, zsurf) |
| Returns a simple surface height field that consists of a single Gaussian-shaped mountain.
|
|
| gaussian_topog_init_r4 |
|
| gaussian_topog_init_r8 |
|
real(kind=fms_top_kind_) function, dimension(size(lon, 1), size(lat, 1)) | get_gaussian_topog_ (lon, lat, height, olond, olatd, wlond, wlatd, rlond, rlatd) |
| The height, position, width, and elongation of the mountain is controlled by optional arguments.
|
|
| get_gaussian_topog_r4 |
|
| get_gaussian_topog_r8 |
|
subroutine | read_namelist |
|
|
logical | do_nml = .true. |
|
real(kind=r8_kind), dimension(maxmts) | height = 0.0_r8_kind |
| height in meters of the gaussian mountiains
|
|
integer, parameter | maxmts = 10 |
|
logical | module_is_initialized = .FALSE. |
|
real(kind=r8_kind), dimension(maxmts) | olat = 0.0_r8_kind |
| Latitude of mountain origins (degrees)
|
|
real(kind=r8_kind), dimension(maxmts) | olon = 0.0_r8_kind |
| longitude of mountain origins (degrees)
|
|
real(kind=r8_kind), dimension(maxmts) | rlat = 0.0_r8_kind |
| Latitude of half-width mountain ridges (degrees) for "standard" gaussian mountain set, rlon/rlat = 0.
|
|
real(kind=r8_kind), dimension(maxmts) | rlon = 0.0_r8_kind |
| Longitude of half-width mountain ridges (degrees) for "standard" gaussian mountain set, rlon/rlat = 0.
|
|
real(kind=r8_kind), dimension(maxmts) | wlat = 0.0_r8_kind |
| Latitude of half-width mountain trails (degrees)
|
|
real(kind=r8_kind), dimension(maxmts) | wlon = 0.0_r8_kind |
| Longitude of half-width mountain trails (degrees)
|
|
Routines for creating Gaussian-shaped land surface topography for latitude-longitude grids.
- Author
- Bruce Wyman
Interfaces generate simple Gaussian-shaped mountains from parameters specified by either argument list or namelist input. The mountain shapes are controlled by the height, half-width, and ridge-width parameters.
◆ gaussian_topog_mod::get_gaussian_topog
interface gaussian_topog_mod::get_gaussian_topog |
Definition at line 54 of file gaussian_topog.F90.
Public Member Functions |
| get_gaussian_topog_r4 |
|
| get_gaussian_topog_r8 |
|
◆ gaussian_topog_init_()
subroutine gaussian_topog_init_ |
( |
real(kind=fms_top_kind_), dimension(:), intent(in) |
lon, |
|
|
real(kind=fms_top_kind_), dimension(:), intent(in) |
lat, |
|
|
real(kind=fms_top_kind_), dimension(:,:), intent(out) |
zsurf |
|
) |
| |
Returns a simple surface height field that consists of a single Gaussian-shaped mountain.
Returns a surface height field that consists of the sum of one or more Gaussian-shaped mountains.
- Parameters
-
[in] | lon | The mean grid box longitude in radians |
[in] | lat | The mean grid box latitude in radians |
[out] | zsurf | The surface height (meters). Size must be size(lon) by size(lat) |
Definition at line 40 of file gaussian_topog.inc.
◆ get_gaussian_topog_()
real(kind=fms_top_kind_) function, dimension(size(lon,1),size(lat,1)) get_gaussian_topog_ |
( |
real(kind=fms_top_kind_), dimension(:), intent(in) |
lon, |
|
|
real(kind=fms_top_kind_), dimension(:), intent(in) |
lat, |
|
|
real(kind=fms_top_kind_), intent(in) |
height, |
|
|
real(kind=fms_top_kind_), intent(in), optional |
olond, |
|
|
real(kind=fms_top_kind_), intent(in), optional |
olatd, |
|
|
real(kind=fms_top_kind_), intent(in), optional |
wlond, |
|
|
real(kind=fms_top_kind_), intent(in), optional |
wlatd, |
|
|
real(kind=fms_top_kind_), intent(in), optional |
rlond, |
|
|
real(kind=fms_top_kind_), intent(in), optional |
rlatd |
|
) |
| |
The height, position, width, and elongation of the mountain is controlled by optional arguments.
- Parameters
-
real | lon The mean grid box longitude in radians. |
real | lat The mean grid box latitude in radians. |
real | height Maximum surface height in meters. |
real | olond, olatd Position/origin of mountain in degrees longitude and latitude. This is the location of the maximum height. |
real | wlond, wlatd Gaussian half-width of mountain in degrees longitude and latitude. |
real | rlond, rlatd Ridge half-width of mountain in degrees longitude and latitude. This is the elongation of the maximum height. |
real | zsurf The surface height (in meters). The size of the returned field is size(lon) by size(lat). </OUT> |
- Exceptions
-
FATAL | shape(zsurf) is not equal to (/size(lon),size(lat)/) Check the input grid size and output field size. The input grid is defined at the midpoint of grid boxes. |
- Note
- Mountains do not wrap around the poles. Returns a land surface topography that consists of a "set" of simple Gaussian-shaped mountains. The height, position, width, and elongation of the mountains can be controlled by variables in the namelist.
Definition at line 100 of file gaussian_topog.inc.
◆ read_namelist()
◆ do_nml
◆ height
real(kind=r8_kind), dimension(maxmts) height = 0.0_r8_kind |
|
private |
◆ maxmts
integer, parameter maxmts = 10 |
|
private |
◆ module_is_initialized
logical module_is_initialized = .FALSE. |
|
private |
◆ olat
real(kind=r8_kind), dimension(maxmts) olat = 0.0_r8_kind |
|
private |
◆ olon
real(kind=r8_kind), dimension(maxmts) olon = 0.0_r8_kind |
|
private |
◆ rlat
real(kind=r8_kind), dimension(maxmts) rlat = 0.0_r8_kind |
|
private |
Latitude of half-width mountain ridges (degrees) for "standard" gaussian mountain set, rlon/rlat = 0.
Definition at line 104 of file gaussian_topog.F90.
◆ rlon
real(kind=r8_kind), dimension(maxmts) rlon = 0.0_r8_kind |
|
private |
Longitude of half-width mountain ridges (degrees) for "standard" gaussian mountain set, rlon/rlat = 0.
Definition at line 102 of file gaussian_topog.F90.
◆ wlat
real(kind=r8_kind), dimension(maxmts) wlat = 0.0_r8_kind |
|
private |
◆ wlon
real(kind=r8_kind), dimension(maxmts) wlon = 0.0_r8_kind |
|
private |