MOM6
MOM_sponge.F90
1 !> Implements sponge regions in isopycnal mode
2 module mom_sponge
3 
4 ! This file is part of MOM6. See LICENSE.md for the license.
5 
6 use mom_coms, only : sum_across_pes
7 use mom_diag_mediator, only : post_data, query_averaging_enabled, register_diag_field
9 use mom_error_handler, only : mom_error, fatal, note, warning, is_root_pe
11 use mom_grid, only : ocean_grid_type
12 use mom_spatial_means, only : global_i_mean
13 use mom_time_manager, only : time_type
15 
16 ! Planned extension: Support for time varying sponge targets.
17 
18 implicit none ; private
19 
20 #include <MOM_memory.h>
21 
22 public set_up_sponge_field, set_up_sponge_ml_density
23 public initialize_sponge, apply_sponge, sponge_end, init_sponge_diags
24 
25 ! A note on unit descriptions in comments: MOM6 uses units that can be rescaled for dimensional
26 ! consistency testing. These are noted in comments with units like Z, H, L, and T, along with
27 ! their mks counterparts with notation like "a velocity [Z T-1 ~> m s-1]". If the units
28 ! vary with the Boussinesq approximation, the Boussinesq variant is given first.
29 
30 !> A structure for creating arrays of pointers to 3D arrays
31 type, public :: p3d
32  real, dimension(:,:,:), pointer :: p => null() !< A pointer to a 3D array
33 end type p3d
34 !> A structure for creating arrays of pointers to 2D arrays
35 type, public :: p2d
36  real, dimension(:,:), pointer :: p => null() !< A pointer to a 2D array
37 end type p2d
38 
39 !> This control structure holds memory and parameters for the MOM_sponge module
40 type, public :: sponge_cs ; private
41  logical :: bulkmixedlayer !< If true, a refined bulk mixed layer is used with
42  !! nkml sublayers and nkbl buffer layer.
43  integer :: nz !< The total number of layers.
44  integer :: isc !< The starting i-index of the computational domain at h.
45  integer :: iec !< The ending i-index of the computational domain at h.
46  integer :: jsc !< The starting j-index of the computational domain at h.
47  integer :: jec !< The ending j-index of the computational domain at h.
48  integer :: isd !< The starting i-index of the data domain at h.
49  integer :: ied !< The ending i-index of the data domain at h.
50  integer :: jsd !< The starting j-index of the data domain at h.
51  integer :: jed !< The ending j-index of the data domain at h.
52  integer :: num_col !< The number of sponge points within the computational domain.
53  integer :: fldno = 0 !< The number of fields which have already been
54  !! registered by calls to set_up_sponge_field
55  integer, pointer :: col_i(:) => null() !< Array of the i-indicies of each of the columns being damped.
56  integer, pointer :: col_j(:) => null() !< Array of the j-indicies of each of the columns being damped.
57  real, pointer :: iresttime_col(:) => null() !< The inverse restoring time of each column.
58  real, pointer :: rcv_ml_ref(:) => null() !< The value toward which the mixed layer
59  !! coordinate-density is being damped [kg m-3].
60  real, pointer :: ref_eta(:,:) => null() !< The value toward which the interface
61  !! heights are being damped [Z ~> m].
62  type(p3d) :: var(max_fields_) !< Pointers to the fields that are being damped.
63  type(p2d) :: ref_val(max_fields_) !< The values to which the fields are damped.
64 
65  logical :: do_i_mean_sponge !< If true, apply sponges to the i-mean fields.
66  real, pointer :: iresttime_im(:) => null() !< The inverse restoring time of
67  !! each row for i-mean sponges.
68  real, pointer :: rcv_ml_ref_im(:) => null() !! The value toward which the i-mean
69  !< mixed layer coordinate-density is being damped [kg m-3].
70  real, pointer :: ref_eta_im(:,:) => null() !< The value toward which the i-mean
71  !! interface heights are being damped [Z ~> m].
72  type(p2d) :: ref_val_im(max_fields_) !< The values toward which the i-means of
73  !! fields are damped.
74 
75  type(diag_ctrl), pointer :: diag => null() !< A structure that is used to
76  !! regulate the timing of diagnostic output.
77  integer :: id_w_sponge = -1 !< A diagnostic ID
78 end type sponge_cs
79 
80 contains
81 
82 !> This subroutine determines the number of points which are within
83 !! sponges in this computational domain. Only points that have
84 !! positive values of Iresttime and which mask2dT indicates are ocean
85 !! points are included in the sponges. It also stores the target interface
86 !! heights.
87 subroutine initialize_sponge(Iresttime, int_height, G, param_file, CS, GV, &
88  Iresttime_i_mean, int_height_i_mean)
89  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
90  real, dimension(SZI_(G),SZJ_(G)), &
91  intent(in) :: iresttime !< The inverse of the restoring time [s-1].
92  real, dimension(SZI_(G),SZJ_(G),SZK_(G)+1), &
93  intent(in) :: int_height !< The interface heights to damp back toward [Z ~> m].
94  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
95  type(sponge_cs), pointer :: cs !< A pointer that is set to point to the control
96  !! structure for this module
97  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
98  real, dimension(SZJ_(G)), &
99  optional, intent(in) :: iresttime_i_mean !< The inverse of the restoring time for
100  !! the zonal mean properties [s-1].
101  real, dimension(SZJ_(G),SZK_(G)+1), &
102  optional, intent(in) :: int_height_i_mean !< The interface heights toward which to
103  !! damp the zonal mean heights [Z ~> m].
104 
105 
106 ! This include declares and sets the variable "version".
107 #include "version_variable.h"
108  character(len=40) :: mdl = "MOM_sponge" ! This module's name.
109  logical :: use_sponge
110  integer :: i, j, k, col, total_sponge_cols
111 
112  if (associated(cs)) then
113  call mom_error(warning, "initialize_sponge called with an associated "// &
114  "control structure.")
115  return
116  endif
117 
118 ! Set default, read and log parameters
119  call log_version(param_file, mdl, version)
120  call get_param(param_file, mdl, "SPONGE", use_sponge, &
121  "If true, sponges may be applied anywhere in the domain. "//&
122  "The exact location and properties of those sponges are "//&
123  "specified from MOM_initialization.F90.", default=.false.)
124 
125  if (.not.use_sponge) return
126  allocate(cs)
127 
128  if (present(iresttime_i_mean) .neqv. present(int_height_i_mean)) &
129  call mom_error(fatal, "initialize_sponge: The optional arguments \n"//&
130  "Iresttime_i_mean and int_height_i_mean must both be present \n"//&
131  "if either one is.")
132 
133  cs%do_i_mean_sponge = present(iresttime_i_mean)
134 
135  cs%nz = g%ke
136 ! CS%isc = G%isc ; CS%iec = G%iec ; CS%jsc = G%jsc ; CS%jec = G%jec
137 ! CS%isd = G%isd ; CS%ied = G%ied ; CS%jsd = G%jsd ; CS%jed = G%jed
138  ! CS%bulkmixedlayer may be set later via a call to set_up_sponge_ML_density.
139  cs%bulkmixedlayer = .false.
140 
141  cs%num_col = 0 ; cs%fldno = 0
142  do j=g%jsc,g%jec ; do i=g%isc,g%iec
143  if ((iresttime(i,j)>0.0) .and. (g%mask2dT(i,j)>0)) &
144  cs%num_col = cs%num_col + 1
145  enddo ; enddo
146 
147  if (cs%num_col > 0) then
148 
149  allocate(cs%Iresttime_col(cs%num_col)) ; cs%Iresttime_col = 0.0
150  allocate(cs%col_i(cs%num_col)) ; cs%col_i = 0
151  allocate(cs%col_j(cs%num_col)) ; cs%col_j = 0
152 
153  col = 1
154  do j=g%jsc,g%jec ; do i=g%isc,g%iec
155  if ((iresttime(i,j)>0.0) .and. (g%mask2dT(i,j)>0)) then
156  cs%col_i(col) = i ; cs%col_j(col) = j
157  cs%Iresttime_col(col) = iresttime(i,j)
158  col = col +1
159  endif
160  enddo ; enddo
161 
162  allocate(cs%Ref_eta(cs%nz+1,cs%num_col))
163  do col=1,cs%num_col ; do k=1,cs%nz+1
164  cs%Ref_eta(k,col) = int_height(cs%col_i(col),cs%col_j(col),k)
165  enddo ; enddo
166 
167  endif
168 
169  if (cs%do_i_mean_sponge) then
170  allocate(cs%Iresttime_im(g%jsd:g%jed)) ; cs%Iresttime_im(:) = 0.0
171  allocate(cs%Ref_eta_im(g%jsd:g%jed,g%ke+1)) ; cs%Ref_eta_im(:,:) = 0.0
172 
173  do j=g%jsc,g%jec
174  cs%Iresttime_im(j) = iresttime_i_mean(j)
175  enddo
176  do k=1,cs%nz+1 ; do j=g%jsc,g%jec
177  cs%Ref_eta_im(j,k) = int_height_i_mean(j,k)
178  enddo ; enddo
179  endif
180 
181  total_sponge_cols = cs%num_col
182  call sum_across_pes(total_sponge_cols)
183 
184  call log_param(param_file, mdl, "!Total sponge columns", total_sponge_cols, &
185  "The total number of columns where sponges are applied.")
186 
187 end subroutine initialize_sponge
188 
189 !> This subroutine sets up diagnostics for the sponges. It is separate
190 !! from initialize_sponge because it requires fields that are not readily
191 !! available where initialize_sponge is called.
192 subroutine init_sponge_diags(Time, G, diag, CS)
193  type(time_type), target, intent(in) :: time !< The current model time
194  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
195  type(diag_ctrl), target, intent(inout) :: diag !< A structure that is used to regulate diagnostic output
196  type(sponge_cs), pointer :: cs !< A pointer to the control structure for this module that
197  !! is set by a previous call to initialize_sponge.
198 
199  if (.not.associated(cs)) return
200 
201  cs%diag => diag
202  cs%id_w_sponge = register_diag_field('ocean_model', 'w_sponge', diag%axesTi, &
203  time, 'The diapycnal motion due to the sponges', 'm s-1')
204 
205 end subroutine init_sponge_diags
206 
207 !> This subroutine stores the reference profile for the variable
208 !! whose address is given by f_ptr. nlay is the number of layers in
209 !! this variable.
210 subroutine set_up_sponge_field(sp_val, f_ptr, G, nlay, CS, sp_val_i_mean)
211  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
212  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
213  intent(in) :: sp_val !< The reference profiles of the quantity being registered.
214  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
215  target, intent(in) :: f_ptr !< a pointer to the field which will be damped
216  integer, intent(in) :: nlay !< the number of layers in this quantity
217  type(sponge_cs), pointer :: cs !< A pointer to the control structure for this module that
218  !! is set by a previous call to initialize_sponge.
219  real, dimension(SZJ_(G),SZK_(G)),&
220  optional, intent(in) :: sp_val_i_mean !< The i-mean reference value for
221  !! this field with i-mean sponges.
222 
223  integer :: j, k, col
224  character(len=256) :: mesg ! String for error messages
225 
226  if (.not.associated(cs)) return
227 
228  cs%fldno = cs%fldno + 1
229 
230  if (cs%fldno > max_fields_) then
231  write(mesg,'("Increase MAX_FIELDS_ to at least ",I3," in MOM_memory.h or decrease &
232  &the number of fields to be damped in the call to &
233  &initialize_sponge." )') cs%fldno
234  call mom_error(fatal,"set_up_sponge_field: "//mesg)
235  endif
236 
237  allocate(cs%Ref_val(cs%fldno)%p(cs%nz,cs%num_col))
238  cs%Ref_val(cs%fldno)%p(:,:) = 0.0
239  do col=1,cs%num_col
240  do k=1,nlay
241  cs%Ref_val(cs%fldno)%p(k,col) = sp_val(cs%col_i(col),cs%col_j(col),k)
242  enddo
243  do k=nlay+1,cs%nz
244  cs%Ref_val(cs%fldno)%p(k,col) = 0.0
245  enddo
246  enddo
247 
248  cs%var(cs%fldno)%p => f_ptr
249 
250  if (nlay/=cs%nz) then
251  write(mesg,'("Danger: Sponge reference fields require nz (",I3,") layers.&
252  & A field with ",I3," layers was passed to set_up_sponge_field.")') &
253  cs%nz, nlay
254  if (is_root_pe()) call mom_error(warning, "set_up_sponge_field: "//mesg)
255  endif
256 
257  if (cs%do_i_mean_sponge) then
258  if (.not.present(sp_val_i_mean)) call mom_error(fatal, &
259  "set_up_sponge_field: sp_val_i_mean must be present with i-mean sponges.")
260 
261  allocate(cs%Ref_val_im(cs%fldno)%p(cs%jsd:cs%jed,cs%nz))
262  cs%Ref_val(cs%fldno)%p(:,:) = 0.0
263  do k=1,cs%nz ; do j=cs%jsc,cs%jec
264  cs%Ref_val_im(cs%fldno)%p(j,k) = sp_val_i_mean(j,k)
265  enddo ; enddo
266  endif
267 
268 end subroutine set_up_sponge_field
269 
270 
271 !> This subroutine stores the reference value for mixed layer density. It is handled differently
272 !! from other values because it is only used in determining which layers can be inflated.
273 subroutine set_up_sponge_ml_density(sp_val, G, CS, sp_val_i_mean)
274  type(ocean_grid_type), intent(in) :: g !< The ocean's grid structure
275  real, dimension(SZI_(G),SZJ_(G)), &
276  intent(in) :: sp_val !< The reference values of the mixed layer density [kg m-3]
277  type(sponge_cs), pointer :: cs !< A pointer to the control structure for this module that is
278  !! set by a previous call to initialize_sponge.
279  real, dimension(SZJ_(G)), &
280  optional, intent(in) :: sp_val_i_mean !< the reference values of the zonal mean mixed
281  !! layer density [kg m-3], for use if Iresttime_i_mean > 0.
282 ! This subroutine stores the reference value for mixed layer density. It is
283 ! handled differently from other values because it is only used in determining
284 ! which layers can be inflated.
285 
286 ! Arguments: sp_val - The reference values of the mixed layer density.
287 ! (in/out) CS - A pointer to the control structure for this module that is
288 ! set by a previous call to initialize_sponge.
289 
290  integer :: j, col
291  character(len=256) :: mesg ! String for error messages
292 
293  if (.not.associated(cs)) return
294 
295  if (associated(cs%Rcv_ml_ref)) then
296  call mom_error(fatal, "set_up_sponge_ML_density appears to have been "//&
297  "called twice.")
298  endif
299 
300  cs%bulkmixedlayer = .true.
301  allocate(cs%Rcv_ml_ref(cs%num_col)) ; cs%Rcv_ml_ref(:) = 0.0
302  do col=1,cs%num_col
303  cs%Rcv_ml_ref(col) = sp_val(cs%col_i(col),cs%col_j(col))
304  enddo
305 
306  if (cs%do_i_mean_sponge) then
307  if (.not.present(sp_val_i_mean)) call mom_error(fatal, &
308  "set_up_sponge_field: sp_val_i_mean must be present with i-mean sponges.")
309 
310  allocate(cs%Rcv_ml_ref_im(cs%jsd:cs%jed)) ; cs%Rcv_ml_ref_im(:) = 0.0
311  do j=cs%jsc,cs%jec
312  cs%Rcv_ml_ref_im(j) = sp_val_i_mean(j)
313  enddo
314  endif
315 
316 end subroutine set_up_sponge_ml_density
317 
318 !> This subroutine applies damping to the layers thicknesses, mixed layer buoyancy, and a variety of
319 !! tracers for every column where there is damping.
320 subroutine apply_sponge(h, dt, G, GV, ea, eb, CS, Rcv_ml)
321  type(ocean_grid_type), intent(inout) :: g !< The ocean's grid structure
322  type(verticalgrid_type), intent(in) :: gv !< The ocean's vertical grid structure
323  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
324  intent(inout) :: h !< Layer thicknesses [H ~> m or kg m-2]
325  real, intent(in) :: dt !< The amount of time covered by this call [s].
326  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
327  intent(inout) :: ea !< An array to which the amount of fluid entrained
328  !! from the layer above during this call will be
329  !! added [H ~> m or kg m-2].
330  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
331  intent(inout) :: eb !< An array to which the amount of fluid entrained
332  !! from the layer below during this call will be
333  !! added [H ~> m or kg m-2].
334  type(sponge_cs), pointer :: cs !< A pointer to the control structure for this module
335  !! that is set by a previous call to initialize_sponge.
336  real, dimension(SZI_(G),SZJ_(G)), &
337  optional, intent(inout) :: rcv_ml !< The coordinate density of the mixed layer [kg m-3].
338 
339 ! This subroutine applies damping to the layers thicknesses, mixed
340 ! layer buoyancy, and a variety of tracers for every column where
341 ! there is damping.
342 
343  ! Local variables
344  real, dimension(SZI_(G), SZJ_(G), SZK_(G)+1) :: &
345  w_int, & ! Water moved upward across an interface within a timestep,
346  ! [H ~> m or kg m-2].
347  e_d ! Interface heights that are dilated to have a value of 0
348  ! at the surface [Z ~> m].
349  real, dimension(SZI_(G), SZJ_(G)) :: &
350  eta_anom, & ! Anomalies in the interface height, relative to the i-mean
351  ! target value [Z ~> m].
352  fld_anom ! Anomalies in a tracer concentration, relative to the
353  ! i-mean target value.
354  real, dimension(SZJ_(G), SZK_(G)+1) :: &
355  eta_mean_anom ! The i-mean interface height anomalies [Z ~> m].
356  real, allocatable, dimension(:,:,:) :: &
357  fld_mean_anom ! THe i-mean tracer concentration anomalies.
358  real, dimension(SZI_(G), SZK_(G)+1) :: &
359  h_above, & ! The total thickness above an interface [H ~> m or kg m-2].
360  h_below ! The total thickness below an interface [H ~> m or kg m-2].
361  real, dimension(SZI_(G)) :: &
362  dilate ! A nondimensional factor by which to dilate layers to
363  ! give 0 at the surface [nondim].
364 
365  real :: e(szk_(g)+1) ! The interface heights [Z ~> m], usually negative.
366  real :: e0 ! The height of the free surface [Z ~> m].
367  real :: e_str ! A nondimensional amount by which the reference
368  ! profile must be stretched for the free surfaces
369  ! heights in the two profiles to agree.
370  real :: w ! The thickness of water moving upward through an
371  ! interface within 1 timestep [H ~> m or kg m-2].
372  real :: wm ! wm is w if w is negative and 0 otherwise [H ~> m or kg m-2].
373  real :: wb ! w at the interface below a layer [H ~> m or kg m-2].
374  real :: wpb ! wpb is wb if wb is positive and 0 otherwise [H ~> m or kg m-2].
375  real :: ea_k, eb_k ! [H ~> m or kg m-2]
376  real :: damp ! The timestep times the local damping coefficient [nondim].
377  real :: i1pdamp ! I1pdamp is 1/(1 + damp). [nondim]
378  real :: damp_1pdamp ! damp_1pdamp is damp/(1 + damp). [nondim]
379  real :: idt ! 1.0/dt [s-1].
380  integer :: c, m, nkmb, i, j, k, is, ie, js, je, nz
381  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
382 
383  if (.not.associated(cs)) return
384  if (cs%bulkmixedlayer) nkmb = gv%nk_rho_varies
385  if (cs%bulkmixedlayer .and. (.not.present(rcv_ml))) &
386  call mom_error(fatal, "Rml must be provided to apply_sponge when using "//&
387  "a bulk mixed layer.")
388 
389  if ((cs%id_w_sponge > 0) .or. cs%do_i_mean_sponge) then
390  do k=1,nz+1 ; do j=js,je ; do i=is,ie
391  w_int(i,j,k) = 0.0
392  enddo ; enddo ; enddo
393  endif
394 
395  if (cs%do_i_mean_sponge) then
396  ! Apply forcing to restore the zonal-mean properties to prescribed values.
397 
398  if (cs%bulkmixedlayer) call mom_error(fatal, "apply_sponge is not yet set up to "//&
399  "work properly with i-mean sponges and a bulk mixed layer.")
400 
401  do j=js,je ; do i=is,ie ; e_d(i,j,nz+1) = -g%bathyT(i,j) ; enddo ; enddo
402  do k=nz,1,-1 ; do j=js,je ; do i=is,ie
403  e_d(i,j,k) = e_d(i,j,k+1) + h(i,j,k)*gv%H_to_Z
404  enddo ; enddo ; enddo
405  do j=js,je
406  do i=is,ie
407  dilate(i) = g%bathyT(i,j) / (e_d(i,j,1) + g%bathyT(i,j))
408  enddo
409  do k=1,nz+1 ; do i=is,ie
410  e_d(i,j,k) = dilate(i) * (e_d(i,j,k) + g%bathyT(i,j)) - g%bathyT(i,j)
411  enddo ; enddo
412  enddo
413 
414  do k=2,nz
415  do j=js,je ; do i=is,ie
416  eta_anom(i,j) = e_d(i,j,k) - cs%Ref_eta_im(j,k)
417  if (cs%Ref_eta_im(j,k) < -g%bathyT(i,j)) eta_anom(i,j) = 0.0
418  enddo ; enddo
419  call global_i_mean(eta_anom(:,:), eta_mean_anom(:,k), g)
420  enddo
421 
422  if (cs%fldno > 0) allocate(fld_mean_anom(g%isd:g%ied,nz,cs%fldno))
423  do m=1,cs%fldno
424  do j=js,je ; do i=is,ie
425  fld_anom(i,j) = cs%var(m)%p(i,j,k) - cs%Ref_val_im(m)%p(j,k)
426  enddo ; enddo
427  call global_i_mean(fld_anom(:,:), fld_mean_anom(:,k,m), g, h(:,:,k))
428  enddo
429 
430  do j=js,je ; if (cs%Iresttime_im(j) > 0.0) then
431  damp = dt*cs%Iresttime_im(j) ; damp_1pdamp = damp / (1.0 + damp)
432 
433  do i=is,ie
434  h_above(i,1) = 0.0 ; h_below(i,nz+1) = 0.0
435  enddo
436  do k=nz,1,-1 ; do i=is,ie
437  h_below(i,k) = h_below(i,k+1) + max(h(i,j,k)-gv%Angstrom_H, 0.0)
438  enddo ; enddo
439  do k=2,nz+1 ; do i=is,ie
440  h_above(i,k) = h_above(i,k-1) + max(h(i,j,k-1)-gv%Angstrom_H, 0.0)
441  enddo ; enddo
442  do k=2,nz
443  ! w is positive for an upward (lightward) flux of mass, resulting
444  ! in the downward movement of an interface.
445  w = damp_1pdamp * eta_mean_anom(j,k) * gv%Z_to_H
446  do i=is,ie
447  if (w > 0.0) then
448  w_int(i,j,k) = min(w, h_below(i,k))
449  eb(i,j,k-1) = eb(i,j,k-1) + w_int(i,j,k)
450  else
451  w_int(i,j,k) = max(w, -h_above(i,k))
452  ea(i,j,k) = ea(i,j,k) - w_int(i,j,k)
453  endif
454  enddo
455  enddo
456  do k=1,nz ; do i=is,ie
457  ea_k = max(0.0, -w_int(i,j,k))
458  eb_k = max(0.0, w_int(i,j,k+1))
459  do m=1,cs%fldno
460  cs%var(m)%p(i,j,k) = (h(i,j,k)*cs%var(m)%p(i,j,k) + &
461  cs%Ref_val_im(m)%p(j,k) * (ea_k + eb_k)) / &
462  (h(i,j,k) + (ea_k + eb_k)) - &
463  damp_1pdamp * fld_mean_anom(j,k,m)
464  enddo
465 
466  h(i,j,k) = max(h(i,j,k) + (w_int(i,j,k+1) - w_int(i,j,k)), &
467  min(h(i,j,k), gv%Angstrom_H))
468  enddo ; enddo
469  endif ; enddo
470 
471  if (cs%fldno > 0) deallocate(fld_mean_anom)
472 
473  endif
474 
475  do c=1,cs%num_col
476 ! c is an index for the next 3 lines but a multiplier for the rest of the loop
477 ! Therefore we use c as per C code and increment the index where necessary.
478  i = cs%col_i(c) ; j = cs%col_j(c)
479  damp = dt*cs%Iresttime_col(c)
480 
481  e(1) = 0.0 ; e0 = 0.0
482  do k=1,nz
483  e(k+1) = e(k) - h(i,j,k)*gv%H_to_Z
484  enddo
485  e_str = e(nz+1) / cs%Ref_eta(nz+1,c)
486 
487  if ( cs%bulkmixedlayer ) then
488  i1pdamp = 1.0 / (1.0 + damp)
489  if (associated(cs%Rcv_ml_ref)) &
490  rcv_ml(i,j) = i1pdamp * (rcv_ml(i,j) + cs%Rcv_ml_ref(c)*damp)
491  do k=1,nkmb
492  do m=1,cs%fldno
493  cs%var(m)%p(i,j,k) = i1pdamp * &
494  (cs%var(m)%p(i,j,k) + cs%Ref_val(m)%p(k,c)*damp)
495  enddo
496  enddo
497 
498  wpb = 0.0; wb = 0.0
499  do k=nz,nkmb+1,-1
500  if (gv%Rlay(k) > rcv_ml(i,j)) then
501  w = min((((e(k)-e0) - e_str*cs%Ref_eta(k,c)) * damp)*gv%Z_to_H, &
502  ((wb + h(i,j,k)) - gv%Angstrom_H))
503  wm = 0.5*(w-abs(w))
504  do m=1,cs%fldno
505  cs%var(m)%p(i,j,k) = (h(i,j,k)*cs%var(m)%p(i,j,k) + &
506  cs%Ref_val(m)%p(k,c)*(damp*h(i,j,k) + (wpb - wm))) / &
507  (h(i,j,k)*(1.0 + damp) + (wpb - wm))
508  enddo
509  else
510  do m=1,cs%fldno
511  cs%var(m)%p(i,j,k) = i1pdamp * &
512  (cs%var(m)%p(i,j,k) + cs%Ref_val(m)%p(k,c)*damp)
513  enddo
514  w = wb + (h(i,j,k) - gv%Angstrom_H)
515  wm = 0.5*(w-abs(w))
516  endif
517  eb(i,j,k) = eb(i,j,k) + wpb
518  ea(i,j,k) = ea(i,j,k) - wm
519  h(i,j,k) = h(i,j,k) + (wb - w)
520  wb = w
521  wpb = w - wm
522  enddo
523 
524  if (wb < 0) then
525  do k=nkmb,1,-1
526  w = min((wb + (h(i,j,k) - gv%Angstrom_H)),0.0)
527  h(i,j,k) = h(i,j,k) + (wb - w)
528  ea(i,j,k) = ea(i,j,k) - w
529  wb = w
530  enddo
531  else
532  w = wb
533  do k=gv%nkml,nkmb
534  eb(i,j,k) = eb(i,j,k) + w
535  enddo
536 
537  k = gv%nkml
538  h(i,j,k) = h(i,j,k) + w
539  do m=1,cs%fldno
540  cs%var(m)%p(i,j,k) = (cs%var(m)%p(i,j,k)*h(i,j,k) + &
541  cs%Ref_val(m)%p(k,c)*w) / (h(i,j,k) + w)
542  enddo
543  endif
544 
545  do k=1,nkmb
546  do m=1,cs%fldno
547  cs%var(m)%p(i,j,k) = i1pdamp * &
548  (cs%var(m)%p(i,j,k) + cs%Ref_val(m)%p(gv%nkml,c)*damp)
549  enddo
550  enddo
551 
552  else ! not BULKMIXEDLAYER
553 
554  wpb = 0.0
555  wb = 0.0
556  do k=nz,1,-1
557  w = min((((e(k)-e0) - e_str*cs%Ref_eta(k,c)) * damp)*gv%Z_to_H, &
558  ((wb + h(i,j,k)) - gv%Angstrom_H))
559  wm = 0.5*(w - abs(w))
560  do m=1,cs%fldno
561  cs%var(m)%p(i,j,k) = (h(i,j,k)*cs%var(m)%p(i,j,k) + &
562  cs%Ref_val(m)%p(k,c) * (damp*h(i,j,k) + (wpb - wm))) / &
563  (h(i,j,k)*(1.0 + damp) + (wpb - wm))
564  enddo
565  eb(i,j,k) = eb(i,j,k) + wpb
566  ea(i,j,k) = ea(i,j,k) - wm
567  h(i,j,k) = h(i,j,k) + (wb - w)
568  wb = w
569  wpb = w - wm
570  enddo
571 
572  endif ! end BULKMIXEDLAYER
573  enddo ! end of c loop
574 
575  if (associated(cs%diag)) then ; if (query_averaging_enabled(cs%diag)) then
576  if (cs%id_w_sponge > 0) then
577  idt = gv%H_to_m / dt
578  do k=1,nz+1 ; do j=js,je ; do i=is,ie
579  w_int(i,j,k) = w_int(i,j,k) * idt ! Scale values by clobbering array since it is local
580  enddo ; enddo ; enddo
581  call post_data(cs%id_w_sponge, w_int(:,:,:), cs%diag)
582  endif
583  endif ; endif
584 
585 end subroutine apply_sponge
586 
587 !> This call deallocates any memory in the sponge control structure.
588 subroutine sponge_end(CS)
589  type(sponge_cs), pointer :: cs !< A pointer to the control structure for this module
590  !! that is set by a previous call to initialize_sponge.
591  integer :: m
592 
593  if (.not.associated(cs)) return
594 
595  if (associated(cs%col_i)) deallocate(cs%col_i)
596  if (associated(cs%col_j)) deallocate(cs%col_j)
597 
598  if (associated(cs%Iresttime_col)) deallocate(cs%Iresttime_col)
599  if (associated(cs%Rcv_ml_ref)) deallocate(cs%Rcv_ml_ref)
600  if (associated(cs%Ref_eta)) deallocate(cs%Ref_eta)
601 
602  if (associated(cs%Iresttime_im)) deallocate(cs%Iresttime_im)
603  if (associated(cs%Rcv_ml_ref_im)) deallocate(cs%Rcv_ml_ref_im)
604  if (associated(cs%Ref_eta_im)) deallocate(cs%Ref_eta_im)
605 
606  do m=1,cs%fldno
607  if (associated(cs%Ref_val(cs%fldno)%p)) deallocate(cs%Ref_val(cs%fldno)%p)
608  if (associated(cs%Ref_val_im(cs%fldno)%p)) &
609  deallocate(cs%Ref_val_im(cs%fldno)%p)
610  enddo
611 
612  deallocate(cs)
613 
614 end subroutine sponge_end
615 
616 !> \namespace mom_sponge
617 !!
618 !! By Robert Hallberg, March 1999-June 2000
619 !!
620 !! This program contains the subroutines that implement sponge
621 !! regions, in which the stratification and water mass properties
622 !! are damped toward some profiles. There are three externally
623 !! callable subroutines in this file.
624 !!
625 !! initialize_sponge determines the mapping from the model
626 !! variables into the arrays of damped columns. This remapping is
627 !! done for efficiency and to conserve memory. Only columns which
628 !! have positive inverse damping times and which are deeper than a
629 !! supplied depth are placed in sponges. The inverse damping
630 !! time is also stored in this subroutine, and memory is allocated
631 !! for all of the reference profiles which will subsequently be
632 !! provided through calls to set_up_sponge_field. The first two
633 !! arguments are a two-dimensional array containing the damping
634 !! rates, and the interface heights to damp towards.
635 !!
636 !! set_up_sponge_field is called to provide a reference profile
637 !! and the location of the field that will be damped back toward
638 !! that reference profile. A third argument, the number of layers
639 !! in the field is also provided, but this should always be nz.
640 !!
641 !! Apply_sponge damps all of the fields that have been registered
642 !! with set_up_sponge_field toward their reference profiles. The
643 !! four arguments are the thickness to be damped, the amount of time
644 !! over which the damping occurs, and arrays to which the movement
645 !! of fluid into a layer from above and below will be added. The
646 !! effect on momentum of the sponge may be accounted for later using
647 !! the movement of water recorded in these later arrays.
648 
649 end module mom_sponge
mom_sponge::p2d
A structure for creating arrays of pointers to 2D arrays.
Definition: MOM_sponge.F90:35
mom_time_manager
Wraps the FMS time manager functions.
Definition: MOM_time_manager.F90:2
mom_spatial_means
Functions and routines to take area, volume, mass-weighted, layerwise, zonal or meridional means.
Definition: MOM_spatial_means.F90:2
mom_verticalgrid
Provides a transparent vertical ocean grid type and supporting routines.
Definition: MOM_verticalGrid.F90:2
mom_file_parser::log_version
An overloaded interface to log version information about modules.
Definition: MOM_file_parser.F90:109
mom_diag_mediator
The subroutines here provide convenient wrappers to the fms diag_manager interfaces with additional d...
Definition: MOM_diag_mediator.F90:3
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_coms
Interfaces to non-domain-oriented communication subroutines, including the MOM6 reproducing sums faci...
Definition: MOM_coms.F90:3
mom_diag_mediator::post_data
Make a diagnostic available for averaging or output.
Definition: MOM_diag_mediator.F90:70
mom_sponge::p3d
A structure for creating arrays of pointers to 3D arrays.
Definition: MOM_sponge.F90:31
mom_verticalgrid::verticalgrid_type
Describes the vertical ocean grid, including unit conversion factors.
Definition: MOM_verticalGrid.F90:24
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:40
mom_file_parser::log_param
An overloaded interface to log the values of various types of parameters.
Definition: MOM_file_parser.F90:96
mom_error_handler
Routines for error handling and I/O management.
Definition: MOM_error_handler.F90:2
mom_grid::ocean_grid_type
Ocean grid type. See mom_grid for details.
Definition: MOM_grid.F90:25
mom_diag_mediator::diag_ctrl
The following data type a list of diagnostic fields an their variants, as well as variables that cont...
Definition: MOM_diag_mediator.F90:239