FV3 Bundle
gaussian_topog.F90
Go to the documentation of this file.
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
10 !*
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14 !* for more details.
15 !*
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 
21 
22 ! <CONTACT EMAIL="Bruce.Wyman@noaa.gov">
23 ! Bruce Wyman
24 ! </CONTACT>
25 
26 ! <HISTORY SRC="http://www.gfdl.noaa.gov/fms-cgi-bin/cvsweb.cgi/FMS/"/>
27 
28 ! <OVERVIEW>
29 ! Routines for creating Gaussian-shaped land surface topography
30 ! for latitude-longitude grids.
31 ! </OVERVIEW>
32 
33 ! <DESCRIPTION>
34 ! Interfaces generate simple Gaussian-shaped mountains from
35 ! parameters specified by either argument list or namelist input.
36 ! The mountain shapes are controlled by the height, half-width,
37 ! and ridge-width parameters.
38 ! </DESCRIPTION>
39 
40 use fms_mod, only: file_exist, open_namelist_file, &
41  check_nml_error, close_file, &
42  stdlog, write_version_number, &
43  mpp_pe, mpp_root_pe, &
44  error_mesg, fatal
45 
46 use constants_mod, only: pi
47 
48 use mpp_mod, only: input_nml_file
49 
50 implicit none
51 private
52 
54 
55 !-----------------------------------------------------------------------
56 ! <NAMELIST NAME="gaussian_topog_nml">
57 ! <DATA NAME="height" UNITS="meter" TYPE="real" DIM="(mxmtns)" DEFAULT="0.">
58 ! Height in meters of the Gaussian mountains.
59 ! </DATA>
60 ! <DATA NAME="olon, olat" UNITS="degree" TYPE="real" DIM="(mxmtns)" DEFAULT="0.">
61 ! The longitude and latitude of mountain origins (in degrees).
62 ! </DATA>
63 ! <DATA NAME="wlon, wlat" UNITS="degree" TYPE="real" DIM="(mxmtns)" DEFAULT="0.">
64 ! The longitude and latitude half-width of mountain tails (in degrees).
65 ! </DATA>
66 ! <DATA NAME="rlon, rlat" UNITS="degree" TYPE="real" DIM="(mxmtns)" DEFAULT="0.">
67 ! The longitude and latitude half-width of mountain ridges (in degrees). For a
68 ! "standard" Gaussian mountain set rlon=rlat=0.
69 ! </DATA>
70 !
71 ! <DATA NAME="NOTE">
72 ! The variables in this namelist are only used when routine
73 ! <TT>gaussian_topog_init</TT> is called. The namelist variables
74 ! are dimensioned (by 10), so that multiple mountains can be generated.
75 !
76 ! Internal parameter mxmtns = 10. By default no mountains are generated.
77 ! </DATA>
78 
79  integer, parameter :: maxmts = 10
80 
81  real, dimension(maxmts) :: height = 0.
82  real, dimension(maxmts) :: olon = 0.
83  real, dimension(maxmts) :: olat = 0.
84  real, dimension(maxmts) :: wlon = 0.
85  real, dimension(maxmts) :: wlat = 0.
86  real, dimension(maxmts) :: rlon = 0.
87  real, dimension(maxmts) :: rlat = 0.
88 
89  namelist /gaussian_topog_nml/ height, olon, olat, wlon, wlat, rlon, rlat
90 ! </NAMELIST>
91 
92 !-----------------------------------------------------------------------
93 
94 ! Include variable "version" to be written to log file.
95 #include<file_version.h>
96 
97 logical :: do_nml = .true.
98 logical :: module_is_initialized = .false.
99 
100 !-----------------------------------------------------------------------
101 
102 contains
103 
104 !#######################################################################
105 
106 ! <SUBROUTINE NAME="gaussian_topog_init">
107 
108 ! <OVERVIEW>
109 ! Returns a surface height field that consists
110 ! of the sum of one or more Gaussian-shaped mountains.
111 ! </OVERVIEW>
112 ! <DESCRIPTION>
113 ! Returns a land surface topography that consists of a "set" of
114 ! simple Gaussian-shaped mountains. The height, position,
115 ! width, and elongation of the mountains can be controlled
116 ! by variables in namelist <LINK SRC="#NAMELIST">&#38;gaussian_topog_nml</LINK>.
117 ! </DESCRIPTION>
118 ! <TEMPLATE>
119 ! <B>call gaussian_topog_init</B> ( lon, lat, zsurf )
120 ! </TEMPLATE>
121 
122 ! <IN NAME="lon" UNITS="radians" TYPE="real" DIM="(:)">
123 ! The mean grid box longitude in radians.
124 ! </IN>
125 ! <IN NAME="lat" UNITS="radians" TYPE="real" DIM="(:)">
126 ! The mean grid box latitude in radians.
127 ! </IN>
128 ! <OUT NAME="zsurf" UNITS="meter" TYPE="real" DIM="(:,:)">
129 ! The surface height (in meters).
130 ! The size of this field must be size(lon) by size(lat).
131 ! </OUT>
132 
133 subroutine gaussian_topog_init ( lon, lat, zsurf )
135 real, intent(in) :: lon(:), lat(:)
136 real, intent(out) :: zsurf(:,:)
137 
138 integer :: n
139 
140  if (.not.module_is_initialized) then
141  call write_version_number("GAUSSIAN_TOPOG_MOD", version)
142  endif
143 
144  if(any(shape(zsurf) /= (/size(lon(:)),size(lat(:))/))) then
145  call error_mesg ('get_gaussian_topog in topography_mod', &
146  'shape(zsurf) is not equal to (/size(lon),size(lat)/)', fatal)
147  endif
148 
149  if (do_nml) call read_namelist
150 
151 ! compute sum of all non-zero mountains
152  zsurf(:,:) = 0.
153  do n = 1, maxmts
154  if ( height(n) == 0. ) cycle
155  zsurf = zsurf + get_gaussian_topog( lon, lat, height(n), &
156  olon(n), olat(n), wlon(n), wlat(n), rlon(n), rlat(n))
157  enddo
158  module_is_initialized = .true.
159 
160 end subroutine gaussian_topog_init
161 ! </SUBROUTINE>
162 
163 !#######################################################################
164 
165 ! <FUNCTION NAME="get_gaussian_topog">
166 
167 ! <OVERVIEW>
168 ! Returns a simple surface height field that consists of a single
169 ! Gaussian-shaped mountain.
170 ! </OVERVIEW>
171 ! <DESCRIPTION>
172 ! Returns a single Gaussian-shaped mountain.
173 ! The height, position, width, and elongation of the mountain
174 ! is controlled by optional arguments.
175 ! </DESCRIPTION>
176 ! <TEMPLATE>
177 ! zsurf = <B>get_gaussian_topog</B> ( lon, lat, height
178 ! [, olond, olatd, wlond, wlatd, rlond, rlatd ] )
179 ! </TEMPLATE>
180 
181 ! <IN NAME="lon" UNITS="radians" TYPE="real" DIM="(:)">
182 ! The mean grid box longitude in radians.
183 ! </IN>
184 ! <IN NAME="lat" UNITS="radians" TYPE="real" DIM="(:)">
185 ! The mean grid box latitude in radians.
186 ! </IN>
187 ! <IN NAME="height" UNITS="meter" TYPE="real" DIM="(scalar)">
188 ! Maximum surface height in meters.
189 ! </IN>
190 ! <IN NAME="olond, olatd" UNITS="degrees" TYPE="real" DIM="(scalar)">
191 ! Position/origin of mountain in degrees longitude and latitude.
192 ! This is the location of the maximum height.
193 ! </IN>
194 ! <IN NAME="wlond, wlatd" UNITS="degrees" TYPE="real" DIM="(scalar)" DEFAULT="15.">
195 ! Gaussian half-width of mountain in degrees longitude and latitude.
196 ! </IN>
197 ! <IN NAME="rlond, rlatd" UNITS="degrees" TYPE="real" DIM="(scalar)" DEFAULT="0.">
198 ! Ridge half-width of mountain in degrees longitude and latitude.
199 ! This is the elongation of the maximum height.
200 ! </IN>
201 ! <OUT NAME="zsurf" UNITS="meter" TYPE="real" DIM="(:,:)">
202 ! The surface height (in meters).
203 ! The size of the returned field is size(lon) by size(lat).
204 ! </OUT>
205 ! <ERROR MSG="shape(zsurf) is not equal to (/size(lon),size(lat)/)" STATUS="FATAL">
206 ! Check the input grid size and output field size.
207 ! The input grid is defined at the midpoint of grid boxes.
208 ! </ERROR>
209 ! <NOTE>
210 ! Mountains do not wrap around the poles.
211 ! </NOTE>
212 
213 function get_gaussian_topog ( lon, lat, height, &
214  olond, olatd, wlond, wlatd, rlond, rlatd ) &
215  result( zsurf )
217 real, intent(in) :: lon(:), lat(:)
218 real, intent(in) :: height
219 real, intent(in), optional :: olond, olatd, wlond, wlatd, rlond, rlatd
220 real :: zsurf(size(lon,1),size(lat,1))
221 
222 integer :: i, j
223 real :: olon, olat, wlon, wlat, rlon, rlat
224 real :: tpi, dtr, dx, dy, xx, yy
225 
226  if (do_nml) call read_namelist
227 
228 ! no need to compute mountain if height=0
229  if ( height == 0. ) then
230  zsurf(:,:) = 0.
231  return
232  endif
233 
234  tpi = 2.0*pi
235  dtr = tpi/360.
236 
237 ! defaults and convert degrees to radians (dtr)
238  olon = 90.*dtr; if (present(olond)) olon=olond*dtr
239  olat = 45.*dtr; if (present(olatd)) olat=olatd*dtr
240  wlon = 15.*dtr; if (present(wlond)) wlon=wlond*dtr
241  wlat = 15.*dtr; if (present(wlatd)) wlat=wlatd*dtr
242  rlon = 0. ; if (present(rlond)) rlon=rlond*dtr
243  rlat = 0. ; if (present(rlatd)) rlat=rlatd*dtr
244 
245 ! compute gaussian-shaped mountain
246  do j=1,size(lat(:))
247  dy = abs(lat(j) - olat) ! dist from y origin
248  yy = max(0., dy-rlat)/wlat
249  do i=1,size(lon(:))
250  dx = abs(lon(i) - olon) ! dist from x origin
251  dx = min(dx, abs(dx-tpi)) ! To ensure that: -pi <= dx <= pi
252  xx = max(0., dx-rlon)/wlon
253  zsurf(i,j) = height*exp(-xx**2 - yy**2)
254  enddo
255  enddo
256 
257 end function get_gaussian_topog
258 ! </FUNCTION>
259 
260 !#######################################################################
261 
262 subroutine read_namelist
264  integer :: unit, ierr, io
265  real :: dtr
266 
267 ! read namelist
268 
269 #ifdef INTERNAL_FILE_NML
270  read (input_nml_file, gaussian_topog_nml, iostat=io)
271  ierr = check_nml_error(io,'gaussian_topog_nml')
272 #else
273  if ( file_exist('input.nml')) then
274  unit = open_namelist_file( )
275  ierr=1; do while (ierr /= 0)
276  read (unit, nml=gaussian_topog_nml, iostat=io, end=10)
277  ierr = check_nml_error(io,'gaussian_topog_nml')
278  enddo
279  10 call close_file (unit)
280  endif
281 #endif
282 
283 ! write version and namelist to log file
284 
285  if (mpp_pe() == mpp_root_pe()) then
286  unit = stdlog()
287  write (unit, nml=gaussian_topog_nml)
288  endif
289 
290  do_nml = .false.
291 
292 end subroutine read_namelist
293 
294 !#######################################################################
295 
296 end module gaussian_topog_mod
297 
298 ! <INFO>
299 ! <NOTE>
300 ! NAMELIST FOR GENERATING GAUSSIAN MOUNTAINS
301 !
302 ! * multiple mountains can be generated
303 ! * the final mountains are the sum of all
304 !
305 ! height = height in meters
306 ! olon, olat = longitude,latitude origin (degrees)
307 ! rlon, rlat = longitude,latitude half-width of ridge (degrees)
308 ! wlon, wlat = longitude,latitude half-width of tail (degrees)
309 !
310 ! Note: For the standard gaussian mountain
311 ! set rlon = rlat = 0 .
312 !
313 ! <PRE>
314 !
315 ! height --> ___________________________
316 ! / \
317 ! / | \
318 ! gaussian / | \
319 ! sides --> / | \
320 ! / olon \
321 ! _____/ olat \______
322 !
323 ! | | |
324 ! |<-->|<----------->|
325 ! |wlon| rlon |
326 ! wlat rlat
327 !
328 ! </PRE>
329 !
330 !See the <LINK SRC="topography.html#TEST PROGRAM">topography </LINK>module documentation for a test program.
331 ! </NOTE>
332 ! </INFO>
Definition: fms.F90:20
real, dimension(maxmts) olat
Definition: mpp.F90:39
integer function, public check_nml_error(IOSTAT, NML_NAME)
Definition: fms.F90:658
character(len=input_str_length), dimension(:), allocatable, target, public input_nml_file
Definition: mpp.F90:1378
subroutine, public gaussian_topog_init(lon, lat, zsurf)
real, dimension(maxmts) height
real, dimension(maxmts) rlon
real, dimension(maxmts) olon
real, dimension(maxmts) wlon
integer, parameter maxmts
real function, dimension(size(lon, 1), size(lat, 1)), public get_gaussian_topog(lon, lat, height, olond, olatd, wlond, wlatd, rlond, rlatd)
************************************************************************GNU Lesser General Public License **This file is part of the GFDL Flexible Modeling System(FMS). ! *! *FMS is free software without even the implied warranty of MERCHANTABILITY or *FITNESS FOR A PARTICULAR PURPOSE See the GNU General Public License *for more details **You should have received a copy of the GNU Lesser General Public *License along with FMS If see< http:! ***********************************************************************subroutine READ_RECORD_CORE_(unit, field, nwords, data, start, axsiz) integer, intent(in) ::unit type(fieldtype), intent(in) ::field integer, intent(in) ::nwords MPP_TYPE_, intent(inout) ::data(nwords) integer, intent(in) ::start(:), axsiz(:) integer(SHORT_KIND) ::i2vals(nwords)!rab used in conjunction with transfer intrinsic to determine size of a variable integer(KIND=1) ::one_byte(8) integer ::word_sz!#ifdef __sgi integer(INT_KIND) ::ivals(nwords) real(FLOAT_KIND) ::rvals(nwords)!#else! integer ::ivals(nwords)! real ::rvals(nwords)!#endif real(DOUBLE_KIND) ::r8vals(nwords) pointer(ptr1, i2vals) pointer(ptr2, ivals) pointer(ptr3, rvals) pointer(ptr4, r8vals) if(mpp_io_stack_size< nwords) call mpp_io_set_stack_size(nwords) call mpp_error(FATAL, 'MPP_READ currently requires use_netCDF option') end subroutine READ_RECORD_CORE_ subroutine READ_RECORD_(unit, field, nwords, data, time_level, domain, position, tile_count, start_in, axsiz_in)!routine that is finally called by all mpp_read routines to perform the read!a non-netCDF record contains:! field ID! a set of 4 coordinates(is:ie, js:je) giving the data subdomain! a timelevel and a timestamp(=NULLTIME if field is static)! 3D real data(stored as 1D)!if you are using direct access I/O, the RECL argument to OPEN must be large enough for the above!in a global direct access file, record position on PE is given by %record.!Treatment of timestamp:! We assume that static fields have been passed without a timestamp.! Here that is converted into a timestamp of NULLTIME.! For non-netCDF fields, field is treated no differently, but is written! with a timestamp of NULLTIME. There is no check in the code to prevent! the user from repeatedly writing a static field. integer, intent(in) ::unit, nwords type(fieldtype), intent(in) ::field MPP_TYPE_, intent(inout) ::data(nwords) integer, intent(in), optional ::time_level type(domain2D), intent(in), optional ::domain integer, intent(in), optional ::position, tile_count integer, intent(in), optional ::start_in(:), axsiz_in(:) integer, dimension(size(field%axes(:))) ::start, axsiz integer ::tlevel !, subdomain(4) integer ::i, error, is, ie, js, je, isg, ieg, jsg, jeg type(domain2d), pointer ::io_domain=> tlevel if(PRESENT(start_in) .AND. PRESENT(axsiz_in)) then if(size(start(! the data domain and compute domain must refer to the subdomain being passed ! In this ! since that attempts to gather all data on PE size(field%axes(:)) axsiz(i)
#define max(a, b)
Definition: mosaic_util.h:33
real, dimension(maxmts) rlat
#define min(a, b)
Definition: mosaic_util.h:32
subroutine, public error_mesg(routine, message, level)
Definition: fms.F90:529
real(fp), parameter, public pi
real, dimension(maxmts) wlat