MOM6
dense_water_initialization.F90
1 !> Initialization routines for the dense water formation
2 !! and overflow experiment.
4 
5 ! This file is part of MOM6. See LICENSE.md for the license.
6 
9 use mom_eos, only : eos_type
10 use mom_error_handler, only : mom_error, fatal
12 use mom_grid, only : ocean_grid_type
13 use mom_sponge, only : sponge_cs
17 
18 implicit none ; private
19 
20 #include <MOM_memory.h>
21 
22 public dense_water_initialize_topography
23 public dense_water_initialize_ts
24 public dense_water_initialize_sponges
25 
26 character(len=40) :: mdl = "dense_water_initialization" !< Module name
27 
28 real, parameter :: default_sill = 0.2 !< Default depth of the sill [nondim]
29 real, parameter :: default_shelf = 0.4 !< Default depth of the shelf [nondim]
30 real, parameter :: default_mld = 0.25 !< Default depth of the mixed layer [nondim]
31 
32 contains
33 
34 !> Initialize the topography field for the dense water experiment
35 subroutine dense_water_initialize_topography(D, G, param_file, max_depth)
36  type(dyn_horgrid_type), intent(in) :: g !< The dynamic horizontal grid type
37  real, dimension(G%isd:G%ied,G%jsd:G%jed), &
38  intent(out) :: d !< Ocean bottom depth in the units of depth_max
39  type(param_file_type), intent(in) :: param_file !< Parameter file structure
40  real, intent(in) :: max_depth !< Maximum ocean depth in arbitrary units
41 
42  ! Local variables
43  real, dimension(5) :: domain_params ! nondimensional widths of all domain sections
44  real :: sill_frac, shelf_frac
45  integer :: i, j
46  real :: x
47 
48  call get_param(param_file, mdl, "DENSE_WATER_DOMAIN_PARAMS", domain_params, &
49  "Fractional widths of all the domain sections for the dense water experiment.\n"//&
50  "As a 5-element vector:\n"//&
51  " - open ocean, the section at maximum depth\n"//&
52  " - downslope, the downward overflow slope\n"//&
53  " - sill separating downslope from upslope\n"//&
54  " - upslope, the upward slope accumulating dense water\n"//&
55  " - the shelf in the dense formation region.", &
56  units="nondim", fail_if_missing=.true.)
57  call get_param(param_file, mdl, "DENSE_WATER_SILL_DEPTH", sill_frac, &
58  "Depth of the sill separating downslope from upslope, as fraction of basin depth.", &
59  units="nondim", default=default_sill)
60  call get_param(param_file, mdl, "DENSE_WATER_SHELF_DEPTH", shelf_frac, &
61  "Depth of the shelf region accumulating dense water for overflow, as fraction of basin depth.", &
62  units="nondim", default=default_shelf)
63 
64  do i = 2, 5
65  ! turn widths into positions
66  domain_params(i) = domain_params(i-1) + domain_params(i)
67  enddo
68 
69  do j = g%jsc,g%jec
70  do i = g%isc,g%iec
71  ! compute normalised zonal coordinate
72  x = (g%geoLonT(i,j) - g%west_lon) / g%len_lon
73 
74  if (x <= domain_params(1)) then
75  ! open ocean region
76  d(i,j) = max_depth
77  elseif (x <= domain_params(2)) then
78  ! downslope region, linear
79  d(i,j) = max_depth - (1.0 - sill_frac) * max_depth * &
80  (x - domain_params(1)) / (domain_params(2) - domain_params(1))
81  elseif (x <= domain_params(3)) then
82  ! sill region
83  d(i,j) = sill_frac * max_depth
84  elseif (x <= domain_params(4)) then
85  ! upslope region
86  d(i,j) = sill_frac * max_depth + (shelf_frac - sill_frac) * max_depth * &
87  (x - domain_params(3)) / (domain_params(4) - domain_params(3))
88  else
89  ! shelf region
90  d(i,j) = shelf_frac * max_depth
91  endif
92  enddo
93  enddo
94 
95 end subroutine dense_water_initialize_topography
96 
97 !> Initialize the temperature and salinity for the dense water experiment
98 subroutine dense_water_initialize_ts(G, GV, param_file, eqn_of_state, T, S, h, just_read_params)
99  type(ocean_grid_type), intent(in) :: g !< Horizontal grid control structure
100  type(verticalgrid_type), intent(in) :: gv !< Vertical grid control structure
101  type(param_file_type), intent(in) :: param_file !< Parameter file structure
102  type(eos_type), pointer :: eqn_of_state !< EOS structure
103  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(out) :: t !< Output temperature [degC]
104  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(out) :: s !< Output salinity [ppt]
105  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2]
106  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
107  !! only read parameters without changing h.
108  ! Local variables
109  real :: mld, s_ref, s_range, t_ref
110  real :: zi, zmid
111  logical :: just_read ! If true, just read parameters but set nothing.
112  integer :: i, j, k, nz
113 
114  nz = gv%ke
115 
116  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
117 
118  call get_param(param_file, mdl, "DENSE_WATER_MLD", mld, &
119  "Depth of unstratified mixed layer as a fraction of the water column.", &
120  units="nondim", default=default_mld, do_not_log=just_read)
121  call get_param(param_file, mdl, "S_REF", s_ref, 'Reference salinity', &
122  default=35.0, units='1e-3', do_not_log=just_read)
123  call get_param(param_file, mdl,"T_REF", t_ref, 'Reference temperature', units='degC', &
124  fail_if_missing=.not.just_read, do_not_log=just_read)
125  call get_param(param_file, mdl,"S_RANGE", s_range, 'Initial salinity range', &
126  units='1e-3', default=2.0, do_not_log=just_read)
127 
128  if (just_read) return ! All run-time parameters have been read, so return.
129 
130  ! uniform temperature everywhere
131  t(:,:,:) = t_ref
132 
133  do j = g%jsc,g%jec
134  do i = g%isc,g%iec
135  zi = 0.
136  do k = 1,nz
137  ! nondimensional middle of layer
138  zmid = zi + 0.5 * h(i,j,k) / (gv%Z_to_H * g%max_depth)
139 
140  if (zmid < mld) then
141  ! use reference salinity in the mixed layer
142  s(i,j,k) = s_ref
143  else
144  ! linear between bottom of mixed layer and bottom
145  s(i,j,k) = s_ref + s_range * (zmid - mld) / (1.0 - mld)
146  endif
147 
148  zi = zi + h(i,j,k) / (gv%Z_to_H * g%max_depth)
149  enddo
150  enddo
151  enddo
152 end subroutine dense_water_initialize_ts
153 
154 !> Initialize the restoring sponges for the dense water experiment
155 subroutine dense_water_initialize_sponges(G, GV, US, tv, param_file, use_ALE, CSp, ACSp)
156  type(ocean_grid_type), intent(in) :: g !< Horizontal grid control structure
157  type(verticalgrid_type), intent(in) :: gv !< Vertical grid control structure
158  type(unit_scale_type), intent(in) :: us !< A dimensional unit scaling type
159  type(thermo_var_ptrs), intent(in) :: tv !< Thermodynamic variables
160  type(param_file_type), intent(in) :: param_file !< Parameter file structure
161  logical, intent(in) :: use_ale !< ALE flag
162  type(sponge_cs), pointer :: csp !< Layered sponge control structure pointer
163  type(ale_sponge_cs), pointer :: acsp !< ALE sponge control structure pointer
164  ! Local variables
165  real :: west_sponge_time_scale, east_sponge_time_scale ! Sponge timescales [T ~> s]
166  real :: west_sponge_width, east_sponge_width
167 
168  real, dimension(SZI_(G),SZJ_(G)) :: idamp ! inverse damping timescale [T-1 ~> s-1]
169  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)) :: h ! sponge thicknesses [H ~> m or kg m-2]
170  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)) :: t ! sponge temperature [degC]
171  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)) :: s ! sponge salinity [ppt]
172  real, dimension(SZK_(GV)+1) :: e0, eta1d ! interface positions for ALE sponge [Z ~> m]
173 
174  integer :: i, j, k, nz
175  real :: x, zi, zmid, dist
176  real :: mld, s_ref, s_range, s_dense, t_ref, sill_height
177 
178  nz = gv%ke
179 
180  call get_param(param_file, mdl, "DENSE_WATER_WEST_SPONGE_TIME_SCALE", west_sponge_time_scale, &
181  "The time scale on the west (outflow) of the domain for restoring. If zero, the sponge is disabled.", &
182  units="s", default=0., scale=us%s_to_T)
183  call get_param(param_file, mdl, "DENSE_WATER_WEST_SPONGE_WIDTH", west_sponge_width, &
184  "The fraction of the domain in which the western (outflow) sponge is active.", &
185  units="nondim", default=0.1)
186  call get_param(param_file, mdl, "DENSE_WATER_EAST_SPONGE_TIME_SCALE", east_sponge_time_scale, &
187  "The time scale on the east (outflow) of the domain for restoring. If zero, the sponge is disabled.", &
188  units="s", default=0., scale=us%s_to_T)
189  call get_param(param_file, mdl, "DENSE_WATER_EAST_SPONGE_WIDTH", east_sponge_width, &
190  "The fraction of the domain in which the eastern (outflow) sponge is active.", &
191  units="nondim", default=0.1)
192 
193  call get_param(param_file, mdl, "DENSE_WATER_EAST_SPONGE_SALT", s_dense, &
194  "Salt anomaly of the dense water being formed in the overflow region.", &
195  units="1e-3", default=4.0)
196 
197  call get_param(param_file, mdl, "DENSE_WATER_MLD", mld, default=default_mld, do_not_log=.true.)
198  call get_param(param_file, mdl, "DENSE_WATER_SILL_HEIGHT", sill_height, default=default_sill, do_not_log=.true.)
199 
200  call get_param(param_file, mdl, "S_REF", s_ref, default=35.0, do_not_log=.true.)
201  call get_param(param_file, mdl, "S_RANGE", s_range, do_not_log=.true.)
202  call get_param(param_file, mdl, "T_REF", t_ref, do_not_log=.true.)
203 
204  ! no active sponges
205  if (west_sponge_time_scale <= 0. .and. east_sponge_time_scale <= 0.) return
206 
207  ! everywhere is initially unsponged
208  idamp(:,:) = 0.0
209 
210  do j = g%jsc, g%jec
211  do i = g%isc,g%iec
212  if (g%mask2dT(i,j) > 0.) then
213  ! nondimensional x position
214  x = (g%geoLonT(i,j) - g%west_lon) / g%len_lon
215 
216  if (west_sponge_time_scale > 0. .and. x < west_sponge_width) then
217  dist = 1. - x / west_sponge_width
218  ! scale restoring by depth into sponge
219  idamp(i,j) = 1. / west_sponge_time_scale * max(0., min(1., dist))
220  elseif (east_sponge_time_scale > 0. .and. x > (1. - east_sponge_width)) then
221  dist = 1. - (1. - x) / east_sponge_width
222  idamp(i,j) = 1. / east_sponge_time_scale * max(0., min(1., dist))
223  endif
224  endif
225  enddo
226  enddo
227 
228  if (use_ale) then
229  ! construct a uniform grid for the sponge
230  do k = 1,nz
231  e0(k) = -g%max_depth * (real(k - 1) / real(nz))
232  enddo
233  e0(nz+1) = -g%max_depth
234 
235  do j = g%jsc,g%jec
236  do i = g%isc,g%iec
237  eta1d(nz+1) = -g%bathyT(i,j)
238  do k = nz,1,-1
239  eta1d(k) = e0(k)
240 
241  if (eta1d(k) < (eta1d(k+1) + gv%Angstrom_Z)) then
242  ! is this layer vanished?
243  eta1d(k) = eta1d(k+1) + gv%Angstrom_Z
244  h(i,j,k) = gv%Angstrom_H
245  else
246  h(i,j,k) = gv%Z_to_H * (eta1d(k) - eta1d(k+1))
247  endif
248  enddo
249  enddo
250  enddo
251 
252  call initialize_ale_sponge(idamp, g, param_file, acsp, h, nz)
253 
254  ! construct temperature and salinity for the sponge
255  ! start with initial condition
256  t(:,:,:) = t_ref
257  s(:,:,:) = s_ref
258 
259  do j = g%jsc,g%jec
260  do i = g%isc,g%iec
261  zi = 0.
262  x = (g%geoLonT(i,j) - g%west_lon) / g%len_lon
263  do k = 1,nz
264  ! nondimensional middle of layer
265  zmid = zi + 0.5 * h(i,j,k) / (gv%Z_to_H * g%max_depth)
266 
267  if (x > (1. - east_sponge_width)) then
268  !if (zmid >= 0.9 * sill_height) &
269  s(i,j,k) = s_ref + s_dense
270  else
271  ! linear between bottom of mixed layer and bottom
272  if (zmid >= mld) &
273  s(i,j,k) = s_ref + s_range * (zmid - mld) / (1.0 - mld)
274  endif
275 
276  zi = zi + h(i,j,k) / (gv%Z_to_H * g%max_depth)
277  enddo
278  enddo
279  enddo
280 
281  if (associated(tv%T)) call set_up_ale_sponge_field(t, g, tv%T, acsp)
282  if (associated(tv%S)) call set_up_ale_sponge_field(s, g, tv%S, acsp)
283  else
284  call mom_error(fatal, "dense_water_initialize_sponges: trying to use non ALE sponge")
285  endif
286 end subroutine dense_water_initialize_sponges
287 
289 
290 !> \namespace dense_water_initialization
291 !!
292 !! This experiment consists of a shelf accumulating dense water, which spills
293 !! over an upward slope and a sill, before flowing down a slope into an open
294 !! ocean region. It's intended as a test of one of the motivating situations
295 !! for the adaptive coordinate.
296 !!
297 !! The nondimensional widths of the 5 regions are controlled by the
298 !! <code>DENSE_WATER_DOMAIN_PARAMS</code>, and the heights of the sill and shelf
299 !! as a fraction of the total domain depth are controlled by
300 !! <code>DENSE_WATER_SILL_HEIGHT</code> and <code>DENSE_WATER_SHELF_HEIGHT</code>.
301 !!
302 !! The density in the domain is governed by a linear equation of state, and
303 !! is set up with a mixed layer of non-dimensional depth <code>DENSE_WATER_MLD</code>
304 !! below which there is a linear stratification from <code>S_REF</code>, increasing
305 !! by <code>S_RANGE</code> to the bottom.
306 !!
307 !! To force the experiment, there are sponges on the inflow and outflow of the
308 !! domain. The inflow sponge has a salinity anomaly of
309 !! <code>DENSE_WATER_EAST_SPONGE_SALT</code> through the entire depth. The outflow
310 !! sponge simply restores to the initial condition. Both sponges have controllable
311 !! widths and restoring timescales.
mom_verticalgrid
Provides a transparent vertical ocean grid type and supporting routines.
Definition: MOM_verticalGrid.F90:2
mom_variables::thermo_var_ptrs
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Definition: MOM_variables.F90:80
mom_ale_sponge::initialize_ale_sponge
Ddetermine the number of points which are within sponges in this computational domain.
Definition: MOM_ALE_sponge.F90:49
mom_dyn_horgrid
Contains a shareable dynamic type for describing horizontal grids and metric data and utilty routines...
Definition: MOM_dyn_horgrid.F90:3
mom_ale_sponge::ale_sponge_cs
ALE sponge control structure.
Definition: MOM_ALE_sponge.F90:86
mom_eos
Provides subroutines for quantities specific to the equation of state.
Definition: MOM_EOS.F90:2
mom_file_parser::param_file_type
A structure that can be parsed to read and document run-time parameters.
Definition: MOM_file_parser.F90:54
mom_file_parser::get_param
An overloaded interface to read and log the values of various types of parameters.
Definition: MOM_file_parser.F90:102
mom_ale_sponge
This module contains the routines used to apply sponge layers when using the ALE mode.
Definition: MOM_ALE_sponge.F90:11
mom_unit_scaling::unit_scale_type
Describes various unit conversion factors.
Definition: MOM_unit_scaling.F90:14
mom_eos::eos_type
A control structure for the equation of state.
Definition: MOM_EOS.F90:108
dense_water_initialization
Initialization routines for the dense water formation and overflow experiment.
Definition: dense_water_initialization.F90:3
mom_ale_sponge::set_up_ale_sponge_field
Store the reference profile at h points for a variable.
Definition: MOM_ALE_sponge.F90:34
mom_verticalgrid::verticalgrid_type
Describes the vertical ocean grid, including unit conversion factors.
Definition: MOM_verticalGrid.F90:24
mom_variables
Provides transparent structures with groups of MOM6 variables and supporting routines.
Definition: MOM_variables.F90:2
mom_file_parser
The MOM6 facility to parse input files for runtime parameters.
Definition: MOM_file_parser.F90:2
mom_sponge
Implements sponge regions in isopycnal mode.
Definition: MOM_sponge.F90:2
mom_grid
Provides the ocean grid type.
Definition: MOM_grid.F90:2
mom_sponge::sponge_cs
This control structure holds memory and parameters for the MOM_sponge module.
Definition: MOM_sponge.F90:41
mom_unit_scaling
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Definition: MOM_unit_scaling.F90:2
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
mom_dyn_horgrid::dyn_horgrid_type
Describes the horizontal ocean grid with only dynamic memory arrays.
Definition: MOM_dyn_horgrid.F90:23
mom_grid::ocean_grid_type
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:26