FMS  2024.03
Flexible Modeling System
gaussian_topog.F90
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 !> @defgroup gaussian_topog_mod gaussian_topog_mod
20 !> @ingroup topography
21 !> @brief Routines for creating Gaussian-shaped land surface topography
22 !! for latitude-longitude grids.
23 !> @author Bruce Wyman
24 !!
25 !! Interfaces generate simple Gaussian-shaped mountains from
26 !! parameters specified by either argument list or namelist input.
27 !! The mountain shapes are controlled by the height, half-width,
28 !! and ridge-width parameters.
29 
30 !> @addtogroup gaussian_topog_mod
31 !> @{
32 module gaussian_topog_mod
33 
34 use fms_mod, only: check_nml_error, &
36  mpp_pe, mpp_root_pe, &
37  error_mesg, fatal
38 
39 use constants_mod, only: pi
40 
41 use mpp_mod, only: input_nml_file
42 use platform_mod, only: r4_kind, r8_kind
43 
44 implicit none
45 private
46 
48 
50  module procedure gaussian_topog_init_r4
51  module procedure gaussian_topog_init_r8
52 end interface gaussian_topog_init
53 
55  module procedure get_gaussian_topog_r4, get_gaussian_topog_r8
56 end interface get_gaussian_topog
57 
58 !! Namelist information for gaussian_topog_nml
59 !!
60 !! The variables in this namelist are only used when routine
61 !! <TT>gaussian_topog_init</TT> is called. The namelist variables
62 !! are dimensioned (by 10), so that multiple mountains can be generated.
63 !!
64 !! Internal parameter mxmtns = 10. By default no mountains are generated.
65 !! </DATA>
66 !!
67 !! NAMELIST FOR GENERATING GAUSSIAN MOUNTAINS
68 !!
69 !! * multiple mountains can be generated
70 !! * the final mountains are the sum of all
71 !!
72 !! height = height in meters
73 !! olon, olat = longitude,latitude origin (degrees)
74 !! rlon, rlat = longitude,latitude half-width of ridge (degrees)
75 !! wlon, wlat = longitude,latitude half-width of tail (degrees)
76 !!
77 !! Note: For the standard gaussian mountain
78 !! set rlon = rlat = 0 .
79 !!
80 !! <PRE>
81 !!
82 !! height --> ___________________________
83 !! / \
84 !! / | \
85 !! gaussian / | \
86 !! sides --> / | \
87 !! / olon \
88 !! _____/ olat \______
89 !!
90 !! | | |
91 !! |<-->|<----------->|
92 !! |wlon| rlon |
93 !! wlat rlat
94 !!
95  integer, parameter :: maxmts = 10
96 
97  real(kind=r8_kind), dimension(maxmts) :: height = 0.0_r8_kind !< height in meters of the gaussian mountiains
98  real(kind=r8_kind), dimension(maxmts) :: olon = 0.0_r8_kind !< longitude of mountain origins (degrees)
99  real(kind=r8_kind), dimension(maxmts) :: olat = 0.0_r8_kind !< Latitude of mountain origins (degrees)
100  real(kind=r8_kind), dimension(maxmts) :: wlon = 0.0_r8_kind !< Longitude of half-width mountain trails (degrees)
101  real(kind=r8_kind), dimension(maxmts) :: wlat = 0.0_r8_kind !< Latitude of half-width mountain trails (degrees)
102  real(kind=r8_kind), dimension(maxmts) :: rlon = 0.0_r8_kind !< Longitude of half-width mountain ridges (degrees)
103  !! for "standard" gaussian mountain set, rlon/rlat = 0
104  real(kind=r8_kind), dimension(maxmts) :: rlat = 0.0_r8_kind !< Latitude of half-width mountain ridges (degrees)
105  !! for "standard" gaussian mountain set, rlon/rlat = 0
106 
107  namelist /gaussian_topog_nml/ height, olon, olat, wlon, wlat, rlon, rlat
108 ! </NAMELIST>
109 
110 !-----------------------------------------------------------------------
111 
112 ! Include variable "version" to be written to log file.
113 #include<file_version.h>
114 
115 logical :: do_nml = .true.
116 logical :: module_is_initialized = .false.
117 
118 !-----------------------------------------------------------------------
119 
120 contains
121 
122 !#######################################################################
123 
124 subroutine read_namelist
125 
126  integer :: iunit, ierr, io
127 
128 !> read namelist
129 
130  read (input_nml_file, gaussian_topog_nml, iostat=io)
131  ierr = check_nml_error(io,'gaussian_topog_nml')
132 
133 !> write version and namelist to log file
134 
135  if (mpp_pe() == mpp_root_pe()) then
136  iunit = stdlog()
137  write (iunit, nml=gaussian_topog_nml)
138  endif
139 
140  do_nml = .false.
141 
142 end subroutine read_namelist
143 
144 !#######################################################################
145 
146 #include "gaussian_topog_r4.fh"
147 #include "gaussian_topog_r8.fh"
148 
149 end module gaussian_topog_mod
150 
151 !> @}
152 ! close documentation grouping
integer function, public check_nml_error(IOSTAT, NML_NAME)
Checks the iostat argument that is returned after reading a namelist and determines if the error code...
Definition: fms.F90:580
subroutine, public write_version_number(version, tag, unit)
Prints to the log file (or a specified unit) the version id string and tag name.
Definition: fms.F90:758
subroutine, public error_mesg(routine, message, level)
Print notes, warnings and error messages; terminates program for warning and error messages....
Definition: fms.F90:498
real(kind=r8_kind), dimension(maxmts) wlon
Longitude of half-width mountain trails (degrees)
real(kind=r8_kind), dimension(maxmts) rlon
Longitude of half-width mountain ridges (degrees) for "standard" gaussian mountain set,...
real(kind=r8_kind), dimension(maxmts) olon
longitude of mountain origins (degrees)
real(kind=r8_kind), dimension(maxmts) olat
Latitude of mountain origins (degrees)
real(kind=r8_kind), dimension(maxmts) height
height in meters of the gaussian mountiains
real(kind=r8_kind), dimension(maxmts) rlat
Latitude of half-width mountain ridges (degrees) for "standard" gaussian mountain set,...
real(kind=r8_kind), dimension(maxmts) wlat
Latitude of half-width mountain trails (degrees)
integer function stdlog()
This function returns the current standard fortran unit numbers for log messages. Log messages,...
Definition: mpp_util.inc:59
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:407