FMS  2025.04
Flexible Modeling System
gaussian_topog.inc
1 !***********************************************************************
2 !* Apache License 2.0
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* Licensed under the Apache License, Version 2.0 (the "License");
7 !* you may not use this file except in compliance with the License.
8 !* You may obtain a copy of the License at
9 !*
10 !* http://www.apache.org/licenses/LICENSE-2.0
11 !*
12 !* FMS is distributed in the hope that it will be useful, but WITHOUT
13 !* WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied;
14 !* without even the implied warranty of MERCHANTABILITY or FITNESS FOR A
15 !* PARTICULAR PURPOSE. See the License for the specific language
16 !* governing permissions and limitations under the License.
17 !***********************************************************************
18 !> @defgroup gaussian_topog_mod gaussian_topog_mod
19 !> @ingroup topography
20 !> @brief Routines for creating Gaussian-shaped land surface topography
21 !! for latitude-longitude grids.
22 !> @author Bruce Wyman
23 !!
24 !! Interfaces generate simple Gaussian-shaped mountains from
25 !! parameters specified by either argument list or namelist input.
26 !! The mountain shapes are controlled by the height, half-width,
27 !! and ridge-width parameters.
28 
29 !> @addtogroup gaussian_topog_mod
30 !> @{
31 
32 !#######################################################################
33 !> @brief Returns a simple surface height field that consists of a single
34 !! Gaussian-shaped mountain.
35 !!
36 !> Returns a surface height field that consists
37 !! of the sum of one or more Gaussian-shaped mountains.
38 !!
39 subroutine gaussian_topog_init_ ( lon, lat, zsurf )
40 
41 real(kind=fms_top_kind_), intent(in) :: lon(:) !< The mean grid box longitude in radians
42 real(kind=fms_top_kind_), intent(in) :: lat(:) !< The mean grid box latitude in radians
43 real(kind=fms_top_kind_), intent(out) :: zsurf(:,:) !< The surface height (meters). Size must be size(lon) by size(lat)
44 
45 integer :: n
46 integer, parameter :: lkind=fms_top_kind_ !local FMS_TOP_KIND_ kind
47 
48  if (.not.module_is_initialized) then
49  call write_version_number("GAUSSIAN_TOPOG_MOD", version)
50  endif
51 
52  if(any(shape(zsurf) /= (/size(lon(:)),size(lat(:))/))) then
53  call error_mesg ('get_gaussian_topog in topography_mod', &
54  'shape(zsurf) is not equal to (/size(lon),size(lat)/)', fatal)
55  endif
56 
57  if (do_nml) call read_namelist
58 
59 ! compute sum of all non-zero mountains
60  zsurf(:,:) = 0.0_lkind
61  do n = 1, maxmts
62  if ( height(n) == 0.0_r8_kind ) cycle
63  zsurf = zsurf + get_gaussian_topog( lon, lat, real(height(n),lkind), &
64  real(olon(n),lkind), real(olat(n),lkind), real(wlon(n),lkind), &
65  real(wlat(n),lkind), real(rlon(n),lkind), real(rlat(n),lkind))
66  enddo
67  module_is_initialized = .true.
68 
69 end subroutine gaussian_topog_init_
70 !#######################################################################
71 !> The height, position, width, and elongation of the mountain
72 !! is controlled by optional arguments.
73 !! @param real lon The mean grid box longitude in radians.
74 !! @param real lat The mean grid box latitude in radians.
75 !! @param real height Maximum surface height in meters.
76 !! @param real olond, olatd Position/origin of mountain in degrees longitude and latitude.
77 !! This is the location of the maximum height.
78 !! @param real wlond, wlatd Gaussian half-width of mountain in degrees longitude and latitude.
79 !! @param real rlond, rlatd Ridge half-width of mountain in degrees longitude and latitude.
80 !! This is the elongation of the maximum height.
81 !! @param real zsurf The surface height (in meters).
82 !! The size of the returned field is size(lon) by size(lat).
83 !! </OUT>
84 !!
85 !! @throws FATAL shape(zsurf) is not equal to (/size(lon),size(lat)/)
86 !! Check the input grid size and output field size.
87 !! The input grid is defined at the midpoint of grid boxes.
88 !!
89 !! @note
90 !! Mountains do not wrap around the poles.
91 !
92 !! <br>Example usage:
93 !! @code{.F90} zsurf = <B>get_gaussian_topog</B> ( lon, lat, height
94 !! [, olond, olatd, wlond, wlatd, rlond, rlatd ] )@endcode
95 !> Returns a land surface topography that consists of a "set" of
96 !! simple Gaussian-shaped mountains. The height, position,
97 !! width, and elongation of the mountains can be controlled
98 !! by variables in the namelist.
99 function get_gaussian_topog_(lon, lat, height, &
100  olond, olatd, wlond, wlatd, rlond, rlatd ) &
101  result(zsurf )
102 
103 real(kind=fms_top_kind_), intent(in) :: lon(:), lat(:)
104 real(kind=fms_top_kind_), intent(in) :: height
105 real(kind=fms_top_kind_), intent(in), optional :: olond, olatd, wlond, wlatd, rlond, rlatd
106 real(kind=fms_top_kind_) :: zsurf(size(lon,1),size(lat,1))
107 
108 integer :: i, j
109 real(kind=fms_top_kind_) :: olon, olat, wlon, wlat, rlon, rlat
110 real(kind=fms_top_kind_) :: tpi, dtr, dx, dy, xx, yy
111 integer, parameter :: lkind = fms_top_kind_
112 
113  if (do_nml) call read_namelist
114 
115 ! no need to compute mountain if height=0
116  if ( height == 0.0_lkind) then
117  zsurf(:,:) = 0.0_lkind
118  return
119  endif
120 
121  tpi = 2.0_lkind*real(pi, fms_top_kind_)
122  dtr = tpi/360.0_lkind
123 
124 ! defaults and convert degrees to radians (dtr)
125  olon = 90.0_r8_kind*real(dtr, r8_kind); if (present(olond)) olon=real(olond*dtr,r8_kind)
126  olat = 45.0_r8_kind*real(dtr, r8_kind); if (present(olatd)) olat=real(olatd*dtr,r8_kind)
127  wlon = 15.0_r8_kind*real(dtr, r8_kind); if (present(wlond)) wlon=real(wlond*dtr,r8_kind)
128  wlat = 15.0_r8_kind*real(dtr, r8_kind); if (present(wlatd)) wlat=real(wlatd*dtr,r8_kind)
129  rlon = 0.0_r8_kind ; if (present(rlond)) rlon=real(rlond*dtr,r8_kind)
130  rlat = 0.0_r8_kind ; if (present(rlatd)) rlat=real(rlatd*dtr,r8_kind)
131 
132 ! compute gaussian-shaped mountain
133  do j=1,size(lat(:))
134  dy = abs(lat(j) - real(olat,lkind)) ! dist from y origin
135  yy = max(0.0_lkind, dy-real(rlat,lkind))/real(wlat,lkind)
136  do i=1,size(lon(:))
137  dx = abs(lon(i) - real(olon,lkind)) ! dist from x origin
138  dx = min(dx, abs(dx-tpi)) ! To ensure that: -pi <= dx <= pi
139  xx = max(0.0_lkind, dx-real(rlon,lkind))/real(wlon,lkind)
140  zsurf(i,j) = real(height,lkind)*exp(-xx**2 - yy**2)
141  enddo
142  enddo
143 
144 end function get_gaussian_topog_
145 
146 !> @}
147 ! close documentation grouping
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.
subroutine gaussian_topog_init_(lon, lat, zsurf)
Returns a simple surface height field that consists of a single Gaussian-shaped mountain.