FMS  2025.04
Flexible Modeling System
gaussian_topog.F90
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 module gaussian_topog_mod
32 
33 use fms_mod, only: check_nml_error, &
35  mpp_pe, mpp_root_pe, &
36  error_mesg, fatal
37 
38 use constants_mod, only: pi
39 
40 use mpp_mod, only: input_nml_file
41 use platform_mod, only: r4_kind, r8_kind
42 
43 implicit none
44 private
45 
47 
49  module procedure gaussian_topog_init_r4
50  module procedure gaussian_topog_init_r8
51 end interface gaussian_topog_init
52 
54  module procedure get_gaussian_topog_r4, get_gaussian_topog_r8
55 end interface get_gaussian_topog
56 
57 !! Namelist information for gaussian_topog_nml
58 !!
59 !! The variables in this namelist are only used when routine
60 !! <TT>gaussian_topog_init</TT> is called. The namelist variables
61 !! are dimensioned (by 10), so that multiple mountains can be generated.
62 !!
63 !! Internal parameter mxmtns = 10. By default no mountains are generated.
64 !! </DATA>
65 !!
66 !! NAMELIST FOR GENERATING GAUSSIAN MOUNTAINS
67 !!
68 !! * multiple mountains can be generated
69 !! * the final mountains are the sum of all
70 !!
71 !! height = height in meters
72 !! olon, olat = longitude,latitude origin (degrees)
73 !! rlon, rlat = longitude,latitude half-width of ridge (degrees)
74 !! wlon, wlat = longitude,latitude half-width of tail (degrees)
75 !!
76 !! Note: For the standard gaussian mountain
77 !! set rlon = rlat = 0 .
78 !!
79 !! <PRE>
80 !!
81 !! height --> ___________________________
82 !! / \
83 !! / | \
84 !! gaussian / | \
85 !! sides --> / | \
86 !! / olon \
87 !! _____/ olat \______
88 !!
89 !! | | |
90 !! |<-->|<----------->|
91 !! |wlon| rlon |
92 !! wlat rlat
93 !!
94  integer, parameter :: maxmts = 10
95 
96  real(kind=r8_kind), dimension(maxmts) :: height = 0.0_r8_kind !< height in meters of the gaussian mountiains
97  real(kind=r8_kind), dimension(maxmts) :: olon = 0.0_r8_kind !< longitude of mountain origins (degrees)
98  real(kind=r8_kind), dimension(maxmts) :: olat = 0.0_r8_kind !< Latitude of mountain origins (degrees)
99  real(kind=r8_kind), dimension(maxmts) :: wlon = 0.0_r8_kind !< Longitude of half-width mountain trails (degrees)
100  real(kind=r8_kind), dimension(maxmts) :: wlat = 0.0_r8_kind !< Latitude of half-width mountain trails (degrees)
101  real(kind=r8_kind), dimension(maxmts) :: rlon = 0.0_r8_kind !< Longitude of half-width mountain ridges (degrees)
102  !! for "standard" gaussian mountain set, rlon/rlat = 0
103  real(kind=r8_kind), dimension(maxmts) :: rlat = 0.0_r8_kind !< Latitude of half-width mountain ridges (degrees)
104  !! for "standard" gaussian mountain set, rlon/rlat = 0
105 
106  namelist /gaussian_topog_nml/ height, olon, olat, wlon, wlat, rlon, rlat
107 ! </NAMELIST>
108 
109 !-----------------------------------------------------------------------
110 
111 ! Include variable "version" to be written to log file.
112 #include<file_version.h>
113 
114 logical :: do_nml = .true.
115 logical :: module_is_initialized = .false.
116 
117 !-----------------------------------------------------------------------
118 
119 contains
120 
121 !#######################################################################
122 
123 subroutine read_namelist
124 
125  integer :: iunit, ierr, io
126 
127 !> read namelist
128 
129  read (input_nml_file, gaussian_topog_nml, iostat=io)
130  ierr = check_nml_error(io,'gaussian_topog_nml')
131 
132 !> write version and namelist to log file
133 
134  if (mpp_pe() == mpp_root_pe()) then
135  iunit = stdlog()
136  write (iunit, nml=gaussian_topog_nml)
137  endif
138 
139  do_nml = .false.
140 
141 end subroutine read_namelist
142 
143 !#######################################################################
144 
145 #include "gaussian_topog_r4.fh"
146 #include "gaussian_topog_r8.fh"
147 
148 end module gaussian_topog_mod
149 
150 !> @}
151 ! 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:523
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:701
subroutine, public error_mesg(routine, message, level)
Print notes, warnings and error messages; terminates program for warning and error messages....
Definition: fms.F90:441
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:58
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:406