18 use random_numbers_mod,
only: initializerandomnumberstream, getrandomnumbers, randomnumberstream
20 implicit none ;
private 22 #include <MOM_memory.h> 24 public neverworld_initialize_topography
25 public neverworld_initialize_thickness
35 subroutine neverworld_initialize_topography(D, G, param_file, max_depth)
37 real,
dimension(G%isd:G%ied,G%jsd:G%jed), &
40 real,
intent(in) :: max_depth
48 # include "version_variable.h" 49 character(len=40) :: mdl =
"Neverworld_initialize_topography" 50 integer :: i, j, is, ie, js, je, isd, ied, jsd, jed
51 real :: nl_roughness_amp, nl_top_amp
52 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
53 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
55 call mom_mesg(
" Neverworld_initialization.F90, Neverworld_initialize_topography: setting topography", 5)
58 call get_param(param_file, mdl,
"NL_ROUGHNESS_AMP", nl_roughness_amp, &
59 "Amplitude of wavy signal in bathymetry.", default=0.05)
60 call get_param(param_file, mdl,
"NL_CONTINENT_AMP", nl_top_amp, &
61 "Scale factor for topography - 0.0 for no continents.", default=1.0)
66 do j=js,je ;
do i=is,ie
67 x = (g%geoLonT(i,j)-g%west_lon) / g%len_lon
68 y =( g%geoLatT(i,j)-g%south_lat) / g%len_lat
70 d(i,j) = 1.0 - 1.1 * spike(y-1,0.12) - 1.1 * spike(y,0.12) - &
72 (1.2 * spike(x,0.2) + 1.2 * spike(x-1.0,0.2)) * spike(min(0.0,y-0.3),0.2) &
73 + 1.2 * spike(x-0.5,0.2) * spike(min(0.0,y-0.55),0.2) &
74 + 1.2 * (spike(x,0.12) + spike(x-1,0.12)) * spike(max(0.0,y-0.06),0.12) &
75 + 0.1 * (cosbell(x,0.1) + cosbell(x-1,0.1)) &
76 + 0.5 * cosbell(x-0.16,0.05) * (cosbell(y-0.18,0.13)**0.4) &
77 + 0.4 * (cosbell(x-0.09,0.08)**0.4) * cosbell(y-0.26,0.05) &
78 + 0.4 * (cosbell(x-0.08,0.08)**0.4) * cosbell(y-0.1,0.05)) &
79 - nl_roughness_amp * cos(14*pi*x) * sin(14*pi*y) &
80 - nl_roughness_amp * cos(20*pi*x) * cos(20*pi*y)
81 if (d(i,j) < 0.0) d(i,j) = 0.0
82 d(i,j) = d(i,j) * max_depth
85 end subroutine neverworld_initialize_topography
88 real function cosbell(x, L)
89 real ,
intent(in) :: x
90 real ,
intent(in) :: l
94 cosbell = 0.5 * (1 + cos(pi*min(abs(x/l),1.0)))
98 real function spike(x, L)
100 real ,
intent(in) :: x
101 real ,
intent(in) :: l
105 spike = (1 - sin(pi*min(abs(x/l),0.5)))
112 real function cone(x, x0, L, clip)
113 real,
intent(in) :: x
114 real,
intent(in) :: x0
115 real,
intent(in) :: l
116 real,
optional,
intent(in) :: clip
118 cone = max( 0., 1. - abs(x - x0) / l )
119 if (
present(clip)) cone = min(clip, cone)
123 real function scurve(x, x0, L)
124 real,
intent(in) :: x
125 real,
intent(in) :: x0
126 real,
intent(in) :: l
129 s = max( 0., min( 1.,( x - x0 ) / l ) )
130 scurve = ( 3. - 2.*s ) * ( s * s )
134 real function cstprof(x, x0, L, lf, bf, sf, sh)
135 real,
intent(in) :: x
136 real,
intent(in) :: x0
137 real,
intent(in) :: l
138 real,
intent(in) :: lf
139 real,
intent(in) :: bf
140 real,
intent(in) :: sf
141 real,
intent(in) :: sh
144 s = max( 0., min( 1.,( x - x0 ) / l ) )
145 cstprof = sh * scurve(s-lf,0.,bf) + (1.-sh) * scurve(s - (1.-sf),0.,sf)
149 real function dist_line_fixed_x(x, y, x0, y0, y1)
150 real,
intent(in) :: x
151 real,
intent(in) :: y
152 real,
intent(in) :: x0
153 real,
intent(in) :: y0
154 real,
intent(in) :: y1
158 yr = min( max(y0,y1), max( min(y0,y1), y ) )
160 dist_line_fixed_x = sqrt( dx*dx + dy*dy )
161 end function dist_line_fixed_x
164 real function dist_line_fixed_y(x, y, x0, x1, y0)
165 real,
intent(in) :: x
166 real,
intent(in) :: y
167 real,
intent(in) :: x0
168 real,
intent(in) :: x1
169 real,
intent(in) :: y0
172 dist_line_fixed_y = dist_line_fixed_x(y, x, y0, x0, x1)
173 end function dist_line_fixed_y
176 real function ns_coast(lon, lat, lon0, lat0, lat1, dlon, sh)
177 real,
intent(in) :: lon
178 real,
intent(in) :: lat
179 real,
intent(in) :: lon0
180 real,
intent(in) :: lat0
181 real,
intent(in) :: lat1
182 real,
intent(in) :: dlon
183 real,
intent(in) :: sh
186 r = dist_line_fixed_x( lon, lat, lon0, lat0, lat1 )
187 ns_coast = cstprof(r, 0., dlon, 0.125, 0.125, 0.5, sh)
188 end function ns_coast
191 real function ew_coast(lon, lat, lon0, lon1, lat0, dlat, sh)
192 real,
intent(in) :: lon
193 real,
intent(in) :: lat
194 real,
intent(in) :: lon0
195 real,
intent(in) :: lon1
196 real,
intent(in) :: lat0
197 real,
intent(in) :: dlat
198 real,
intent(in) :: sh
201 r = dist_line_fixed_y( lon, lat, lon0, lon1, lat0 )
202 ew_coast = cstprof(r, 0., dlat, 0.125, 0.125, 0.5, sh)
203 end function ew_coast
206 real function ns_ridge(lon, lat, lon0, lat0, lat1, dlon, rh)
207 real,
intent(in) :: lon
208 real,
intent(in) :: lat
209 real,
intent(in) :: lon0
210 real,
intent(in) :: lat0
211 real,
intent(in) :: lat1
212 real,
intent(in) :: dlon
213 real,
intent(in) :: rh
216 r = dist_line_fixed_x( lon, lat, lon0, lat0, lat1 )
217 ns_ridge = 1. - rh * cone(r, 0., dlon)
218 end function ns_ridge
222 real function circ_ridge(lon, lat, lon0, lat0, ring_radius, ring_thickness, ridge_height)
223 real,
intent(in) :: lon
224 real,
intent(in) :: lat
225 real,
intent(in) :: lon0
226 real,
intent(in) :: lat0
227 real,
intent(in) :: ring_radius
228 real,
intent(in) :: ring_thickness
229 real,
intent(in) :: ridge_height
232 r = sqrt( (lon - lon0)**2 + (lat - lat0)**2 )
233 r = abs( r - ring_radius)
234 r = cone(r, 0., ring_thickness, ridge_height)
236 end function circ_ridge
242 subroutine neverworld_initialize_thickness(h, G, GV, US, param_file, eqn_of_state, P_ref)
246 real,
intent(out),
dimension(SZI_(G),SZJ_(G),SZK_(GV)) :: h
251 type(
eos_type),
pointer :: eqn_of_state
252 real,
intent(in) :: p_ref
255 real :: e0(szk_(g)+1)
257 real,
dimension(SZK_(G)) :: h_profile
263 type(randomnumberstream) :: rns
264 character(len=40) :: mdl =
"Neverworld_initialize_thickness" 265 integer :: i, j, k, k1, is, ie, js, je, nz, itt
267 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
269 call mom_mesg(
" Neverworld_initialization.F90, Neverworld_initialize_thickness: setting thickness", 5)
270 call get_param(param_file, mdl,
"INIT_THICKNESS_PROFILE", h_profile, &
271 "Profile of initial layer thicknesses.", units=
"m", scale=us%m_to_Z, &
272 fail_if_missing=.true.)
273 call get_param(param_file, mdl,
"NL_THICKNESS_PERT_AMP", pert_amp, &
274 "Amplitude of finite scale perturbations as fraction of depth.", &
275 units=
"nondim", default=0.)
276 call get_param(param_file, mdl,
"NL_THICKNESS_NOISE_AMP", h_noise, &
277 "Amplitude of noise to scale layer by.", units=
"nondim", default=0.)
282 e0(k+1) = e0(k) - h_profile(k)
285 do j=js,je ;
do i=is,ie
286 e_interface = -g%bathyT(i,j)
288 h(i,j,k) = gv%Z_to_H * (e0(k) - e_interface)
289 x=(g%geoLonT(i,j)-g%west_lon)/g%len_lon
290 y=(g%geoLatT(i,j)-g%south_lat)/g%len_lat
291 r1=sqrt((x-0.7)**2+(y-0.2)**2)
292 r2=sqrt((x-0.3)**2+(y-0.25)**2)
293 h(i,j,k) = h(i,j,k) + pert_amp * (e0(k) - e0(nz+1)) * gv%Z_to_H * &
294 (spike(r1,0.15)-spike(r2,0.15))
295 if (h_noise /= 0.)
then 296 rns = initializerandomnumberstream( int( 4096*(x + (y+1.)) ) )
297 call getrandomnumbers(rns, noise)
298 noise = h_noise * 2. * ( noise - 0.5 )
299 h(i,j,k) = ( 1. + noise ) * h(i,j,k)
301 h(i,j,k) = max( gv%Angstrom_H, h(i,j,k) )
302 e_interface = e_interface + gv%H_to_Z * h(i,j,k)
304 h(i,j,1) = gv%Z_to_H * (e0(1) - e_interface)
305 h(i,j,1) = max( gv%Angstrom_H, h(i,j,1) )
308 end subroutine neverworld_initialize_thickness
Ocean grid type. See mom_grid for details.
Initialization for the "Neverworld" configuration.
Calculates density of sea water from T, S and P.
A structure that can be parsed to read and document run-time parameters.
Provides the ocean grid type.
Describes the horizontal ocean grid with only dynamic memory arrays.
The MOM6 facility to parse input files for runtime parameters.
Describes various unit conversion factors.
This module contains the tracer_registry_type and the subroutines that handle registration of tracers...
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Type to carry basic tracer information.
Routines for error handling and I/O management.
This control structure holds memory and parameters for the MOM_sponge module.
Implements sponge regions in isopycnal mode.
Calculate the derivatives of density with temperature and salinity from T, S, and P...
Provides subroutines for quantities specific to the equation of state.
An overloaded interface to log version information about modules.
Describes the vertical ocean grid, including unit conversion factors.
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Contains a shareable dynamic type for describing horizontal grids and metric data and utilty routines...
Provides a transparent vertical ocean grid type and supporting routines.
Provides transparent structures with groups of MOM6 variables and supporting routines.
An overloaded interface to read and log the values of various types of parameters.
A control structure for the equation of state.