Routines for creating land surface topography fields and land-water masks for latitude-longitude grids.
idealized_topog_mod
horiz_interp_mod
fms_mod
fms_io_mod
get_topog_mean
| blon | The longitude (in radians) at grid box boundaries. [real, dimension(:)] |
| blat | The latitude (in radians) at grid box boundaries. [real, dimension(:)] |
| zmean | The mean surface height (meters).
The size of this field must be size(blon)-1 by size(blat)-1. [real, dimension(:,:)] |
| get_topog_mean | A logical value of TRUE is returned if the surface height field
was successfully created. A value of FALSE may be returned if the
input topography data set was not readable. [logical] |
get_topog_stdev
| blon | The longitude (in radians) at grid box boundaries. [real, dimension(:)] |
| blat | The latitude (in radians) at grid box boundaries. [real, dimension(:)] |
| stdev | The standard deviation of surface height (in meters) within
given input model grid boxes.
The size of this field must be size(blon)-1 by size(blat)-1. [real, dimension(:,:)] |
| get_topog_stdev | A logical value of TRUE is returned if the output field was
successfully created. A value of FALSE may be returned if the
input topography data set was not readable. [logical] |
get_ocean_frac
| blon | The longitude (in radians) at grid box boundaries. [real, dimension(:)] |
| blat | The latitude (in radians) at grid box boundaries. [real, dimension(:)] |
| ocean_frac | The fractional amount (0 to 1) of ocean in a grid box.
The size of this field must be size(blon)-1 by size(blat)-1. [real, dimension(:,:)] |
| get_ocean_frac | A logical value of TRUE is returned if the output field
was successfully created. A value of FALSE may be returned
if the Navy 1/6 degree percent water data set was not readable. [logical] |
get_ocean_mask
| blon | The longitude (in radians) at grid box boundaries. [real, dimension(:)] |
| blat | The latitude (in radians) at grid box boundaries. [real, dimension(:)] |
| ocean_frac | The fractional amount (0 to 1) of ocean in a grid box.
The size of this field must be size(blon)-1 by size(blat)-1. [real, dimension(:,:)] |
| get_ocean_mask | A logical value of TRUE is returned if the output field
was successfully created. A value of FALSE may be returned
if the Navy 1/6 degree percent water data set was not readable. [logical] |
get_water_frac
| blon | The longitude (in radians) at grid box boundaries. [real, dimension(:)] |
| blat | The latitude (in radians) at grid box boundaries. [real, dimension(:)] |
| water_frac | The fractional amount (0 to 1) of water in a grid box.
The size of this field must be size(blon)-1 by size(blat)-1. [real, dimension(:,:)] |
| get_water_frac | A logical value of TRUE is returned if the output field
was successfully created. A value of FALSE may be returned
if the Navy 1/6 degree percent water data set was not readable. [logical] |
get_water_mask
| blon | The longitude (in radians) at grid box boundaries. [real, dimension(:)] |
| blat | The latitude (in radians) at grid box boundaries. [real, dimension(:)] |
| water_mask | A binary mask for water (true) or land (false).
The size of this field must be size(blon)-1 by size(blat)-1. [real, dimension(:,:)] |
| get_water_mask | A logical value of TRUE is returned if the output field
was successfully created. A value of FALSE may be returned
if the Navy 1/6 degree percent water data set was not readable. [logical] |
record = 1 nlon, nlat
record = 2 blon, blat
record = 3 data where: nlon, nlat = The number of longitude and latitude points
in the horizontal grid. For the 1/6 degree
data sets this is 2160 x 1080. [integer]
blon, blat = The longitude and latitude grid box boundaries in degrees.
[real :: blon(nlon+1), blat(nlat+1)]
data = The topography or percent water data.
[real :: data(nlon,nlat)]
use topography_mod
implicit none
integer, parameter :: nlon=24, nlat=18
real :: x(nlon), y(nlat), xb(nlon+1), yb(nlat+1), z(nlon,nlat)
real :: hpi, rtd
integer :: i,j
logical :: a
gaussian mountain parameters
real, parameter :: ht=4000.
real, parameter :: x0=90., y0=45. ! origin in degrees
real, parameter :: xw=15., yw=15. ! half-width in degees
real, parameter :: xr=30., yr= 0. ! ridge-width in degrees
create lat/lon grid in radians
hpi = acos(0.0)
rtd = 90./hpi ! rad to deg
do i=1,nlon
xb(i) = 4.*hpi*real(i-1)/real(nlon)
enddo
xb(nlon+1) = xb(1)+4.*hpi
yb(1) = -hpi
do j=2,nlat
yb(j) = yb(j-1) + 2.*hpi/real(nlat)
enddo
yb(nlat+1) = hpi
mid-point of grid boxes
x(1:nlon) = 0.5*(xb(1:nlon)+xb(2:nlon+1))
y(1:nlat) = 0.5*(yb(1:nlat)+yb(2:nlat+1))
test topography_mod routines
a = get_topog_mean(xb,yb,z)
call printz ('get_topog_mean')
a = get_water_frac(xb,yb,z)
z = z*100. ! in percent
call printz ('get_water_frac')
a = get_ocean_frac(xb,yb,z)
z = z*100. ! in percent
call printz ('get_ocean_frac')
test idealized_topog_mod routines
a = .true.
z = get_gaussian_topog(x,y,ht,x0,y0,xw,yw,xr,yr)
call printz ('get_gaussian_topog')
call gaussian_topog_init (x,y,z)
call printz ('gaussian_topog_init')
contains
simple printout of topog/water array
subroutine printz (lab)
character(len=*), intent(in) :: lab
if (a) then
print '(/a)', trim(lab)
else
print '(/a)', 'no data available: '//trim(lab)
return
endif
print full grid
print '(3x,25i5)', (nint(x(i)*rtd),i=1,nlon)
do j=nlat,1,-1
print '(i3,25i5)', nint(y(j)*rtd), (nint(z(i,j)),i=1,nlon)
enddo
end subroutine printz
end program test
Water mask produces some possible erroneous water points along the coast of Antarctic (at about 90W).
Use of netcdf data sets.
Incorporate other topography and ocean data sets.