MOM6
mom_sponge Module Reference

Detailed Description

Implements sponge regions in isopycnal mode.

By Robert Hallberg, March 1999-June 2000

This program contains the subroutines that implement sponge regions, in which the stratification and water mass properties are damped toward some profiles. There are three externally callable subroutines in this file.

initialize_sponge determines the mapping from the model variables into the arrays of damped columns. This remapping is done for efficiency and to conserve memory. Only columns which have positive inverse damping times and which are deeper than a supplied depth are placed in sponges. The inverse damping time is also stored in this subroutine, and memory is allocated for all of the reference profiles which will subsequently be provided through calls to set_up_sponge_field. The first two arguments are a two-dimensional array containing the damping rates, and the interface heights to damp towards.

set_up_sponge_field is called to provide a reference profile and the location of the field that will be damped back toward that reference profile. A third argument, the number of layers in the field is also provided, but this should always be nz.

Apply_sponge damps all of the fields that have been registered with set_up_sponge_field toward their reference profiles. The four arguments are the thickness to be damped, the amount of time over which the damping occurs, and arrays to which the movement of fluid into a layer from above and below will be added. The effect on momentum of the sponge may be accounted for later using the movement of water recorded in these later arrays.

Data Types

type  p2d
 A structure for creating arrays of pointers to 2D arrays. More...
 
type  p3d
 A structure for creating arrays of pointers to 3D arrays. More...
 
type  sponge_cs
 This control structure holds memory and parameters for the MOM_sponge module. More...
 

Functions/Subroutines

subroutine, public initialize_sponge (Iresttime, int_height, G, param_file, CS, GV, Iresttime_i_mean, int_height_i_mean)
 This subroutine determines the number of points which are within sponges in this computational domain. Only points that have positive values of Iresttime and which mask2dT indicates are ocean points are included in the sponges. It also stores the target interface heights. More...
 
subroutine, public init_sponge_diags (Time, G, GV, US, diag, CS)
 This subroutine sets up diagnostics for the sponges. It is separate from initialize_sponge because it requires fields that are not readily available where initialize_sponge is called. More...
 
subroutine, public set_up_sponge_field (sp_val, f_ptr, G, nlay, CS, sp_val_i_mean)
 This subroutine stores the reference profile for the variable whose address is given by f_ptr. nlay is the number of layers in this variable. More...
 
subroutine, public set_up_sponge_ml_density (sp_val, G, CS, sp_val_i_mean)
 This subroutine stores the reference value for mixed layer density. It is handled differently from other values because it is only used in determining which layers can be inflated. More...
 
subroutine, public apply_sponge (h, dt, G, GV, US, ea, eb, CS, Rcv_ml)
 This subroutine applies damping to the layers thicknesses, mixed layer buoyancy, and a variety of tracers for every column where there is damping. More...
 
subroutine, public sponge_end (CS)
 This call deallocates any memory in the sponge control structure. More...
 

Function/Subroutine Documentation

◆ apply_sponge()

subroutine, public mom_sponge::apply_sponge ( real, dimension(szi_(g),szj_(g),szk_(g)), intent(inout)  h,
real, intent(in)  dt,
type(ocean_grid_type), intent(inout)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(inout)  ea,
real, dimension(szi_(g),szj_(g),szk_(g)), intent(inout)  eb,
type(sponge_cs), pointer  CS,
real, dimension(szi_(g),szj_(g)), intent(inout), optional  Rcv_ml 
)

This subroutine applies damping to the layers thicknesses, mixed layer buoyancy, and a variety of tracers for every column where there is damping.

Parameters
[in,out]gThe ocean's grid structure
[in]gvThe ocean's vertical grid structure
[in]usA dimensional unit scaling type
[in,out]hLayer thicknesses [H ~> m or kg m-2]
[in]dtThe amount of time covered by this call [T ~> s].
[in,out]eaAn array to which the amount of fluid entrained
[in,out]ebAn array to which the amount of fluid entrained
csA pointer to the control structure for this module that is set by a previous call to initialize_sponge.
[in,out]rcv_mlThe coordinate density of the mixed layer [R ~> kg m-3].

Definition at line 324 of file MOM_sponge.F90.

324  type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
325  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
326  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
327  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
328  intent(inout) :: h !< Layer thicknesses [H ~> m or kg m-2]
329  real, intent(in) :: dt !< The amount of time covered by this call [T ~> s].
330  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
331  intent(inout) :: ea !< An array to which the amount of fluid entrained
332  !! from the layer above during this call will be
333  !! added [H ~> m or kg m-2].
334  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
335  intent(inout) :: eb !< An array to which the amount of fluid entrained
336  !! from the layer below during this call will be
337  !! added [H ~> m or kg m-2].
338  type(sponge_CS), pointer :: CS !< A pointer to the control structure for this module
339  !! that is set by a previous call to initialize_sponge.
340  real, dimension(SZI_(G),SZJ_(G)), &
341  optional, intent(inout) :: Rcv_ml !< The coordinate density of the mixed layer [R ~> kg m-3].
342 
343 ! This subroutine applies damping to the layers thicknesses, mixed
344 ! layer buoyancy, and a variety of tracers for every column where
345 ! there is damping.
346 
347  ! Local variables
348  real, dimension(SZI_(G), SZJ_(G), SZK_(G)+1) :: &
349  w_int, & ! Water moved upward across an interface within a timestep,
350  ! [H ~> m or kg m-2].
351  e_d ! Interface heights that are dilated to have a value of 0
352  ! at the surface [Z ~> m].
353  real, dimension(SZI_(G), SZJ_(G)) :: &
354  eta_anom, & ! Anomalies in the interface height, relative to the i-mean
355  ! target value [Z ~> m].
356  fld_anom ! Anomalies in a tracer concentration, relative to the
357  ! i-mean target value.
358  real, dimension(SZJ_(G), SZK_(G)+1) :: &
359  eta_mean_anom ! The i-mean interface height anomalies [Z ~> m].
360  real, allocatable, dimension(:,:,:) :: &
361  fld_mean_anom ! THe i-mean tracer concentration anomalies.
362  real, dimension(SZI_(G), SZK_(G)+1) :: &
363  h_above, & ! The total thickness above an interface [H ~> m or kg m-2].
364  h_below ! The total thickness below an interface [H ~> m or kg m-2].
365  real, dimension(SZI_(G)) :: &
366  dilate ! A nondimensional factor by which to dilate layers to
367  ! give 0 at the surface [nondim].
368 
369  real :: e(SZK_(G)+1) ! The interface heights [Z ~> m], usually negative.
370  real :: e0 ! The height of the free surface [Z ~> m].
371  real :: e_str ! A nondimensional amount by which the reference
372  ! profile must be stretched for the free surfaces
373  ! heights in the two profiles to agree.
374  real :: w ! The thickness of water moving upward through an
375  ! interface within 1 timestep [H ~> m or kg m-2].
376  real :: wm ! wm is w if w is negative and 0 otherwise [H ~> m or kg m-2].
377  real :: wb ! w at the interface below a layer [H ~> m or kg m-2].
378  real :: wpb ! wpb is wb if wb is positive and 0 otherwise [H ~> m or kg m-2].
379  real :: ea_k, eb_k ! [H ~> m or kg m-2]
380  real :: damp ! The timestep times the local damping coefficient [nondim].
381  real :: I1pdamp ! I1pdamp is 1/(1 + damp). [nondim]
382  real :: damp_1pdamp ! damp_1pdamp is damp/(1 + damp). [nondim]
383  real :: Idt ! 1.0/dt times a height unit conversion factor [m H-1 T-1 ~> s-1 or m3 kg-1 s-1].
384  integer :: c, m, nkmb, i, j, k, is, ie, js, je, nz
385  is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
386 
387  if (.not.associated(cs)) return
388  if (cs%bulkmixedlayer) nkmb = gv%nk_rho_varies
389  if (cs%bulkmixedlayer .and. (.not.present(rcv_ml))) &
390  call mom_error(fatal, "Rml must be provided to apply_sponge when using "//&
391  "a bulk mixed layer.")
392 
393  if ((cs%id_w_sponge > 0) .or. cs%do_i_mean_sponge) then
394  do k=1,nz+1 ; do j=js,je ; do i=is,ie
395  w_int(i,j,k) = 0.0
396  enddo ; enddo ; enddo
397  endif
398 
399  if (cs%do_i_mean_sponge) then
400  ! Apply forcing to restore the zonal-mean properties to prescribed values.
401 
402  if (cs%bulkmixedlayer) call mom_error(fatal, "apply_sponge is not yet set up to "//&
403  "work properly with i-mean sponges and a bulk mixed layer.")
404 
405  do j=js,je ; do i=is,ie ; e_d(i,j,nz+1) = -g%bathyT(i,j) ; enddo ; enddo
406  do k=nz,1,-1 ; do j=js,je ; do i=is,ie
407  e_d(i,j,k) = e_d(i,j,k+1) + h(i,j,k)*gv%H_to_Z
408  enddo ; enddo ; enddo
409  do j=js,je
410  do i=is,ie
411  dilate(i) = g%bathyT(i,j) / (e_d(i,j,1) + g%bathyT(i,j))
412  enddo
413  do k=1,nz+1 ; do i=is,ie
414  e_d(i,j,k) = dilate(i) * (e_d(i,j,k) + g%bathyT(i,j)) - g%bathyT(i,j)
415  enddo ; enddo
416  enddo
417 
418  do k=2,nz
419  do j=js,je ; do i=is,ie
420  eta_anom(i,j) = e_d(i,j,k) - cs%Ref_eta_im(j,k)
421  if (cs%Ref_eta_im(j,k) < -g%bathyT(i,j)) eta_anom(i,j) = 0.0
422  enddo ; enddo
423  call global_i_mean(eta_anom(:,:), eta_mean_anom(:,k), g, tmp_scale=us%Z_to_m)
424  enddo
425 
426  if (cs%fldno > 0) allocate(fld_mean_anom(g%isd:g%ied,nz,cs%fldno))
427  do m=1,cs%fldno
428  do j=js,je ; do i=is,ie
429  fld_anom(i,j) = cs%var(m)%p(i,j,k) - cs%Ref_val_im(m)%p(j,k)
430  enddo ; enddo
431  call global_i_mean(fld_anom(:,:), fld_mean_anom(:,k,m), g, h(:,:,k))
432  enddo
433 
434  do j=js,je ; if (cs%Iresttime_im(j) > 0.0) then
435  damp = dt * cs%Iresttime_im(j) ; damp_1pdamp = damp / (1.0 + damp)
436 
437  do i=is,ie
438  h_above(i,1) = 0.0 ; h_below(i,nz+1) = 0.0
439  enddo
440  do k=nz,1,-1 ; do i=is,ie
441  h_below(i,k) = h_below(i,k+1) + max(h(i,j,k)-gv%Angstrom_H, 0.0)
442  enddo ; enddo
443  do k=2,nz+1 ; do i=is,ie
444  h_above(i,k) = h_above(i,k-1) + max(h(i,j,k-1)-gv%Angstrom_H, 0.0)
445  enddo ; enddo
446  do k=2,nz
447  ! w is positive for an upward (lightward) flux of mass, resulting
448  ! in the downward movement of an interface.
449  w = damp_1pdamp * eta_mean_anom(j,k) * gv%Z_to_H
450  do i=is,ie
451  if (w > 0.0) then
452  w_int(i,j,k) = min(w, h_below(i,k))
453  eb(i,j,k-1) = eb(i,j,k-1) + w_int(i,j,k)
454  else
455  w_int(i,j,k) = max(w, -h_above(i,k))
456  ea(i,j,k) = ea(i,j,k) - w_int(i,j,k)
457  endif
458  enddo
459  enddo
460  do k=1,nz ; do i=is,ie
461  ea_k = max(0.0, -w_int(i,j,k))
462  eb_k = max(0.0, w_int(i,j,k+1))
463  do m=1,cs%fldno
464  cs%var(m)%p(i,j,k) = (h(i,j,k)*cs%var(m)%p(i,j,k) + &
465  cs%Ref_val_im(m)%p(j,k) * (ea_k + eb_k)) / &
466  (h(i,j,k) + (ea_k + eb_k)) - &
467  damp_1pdamp * fld_mean_anom(j,k,m)
468  enddo
469 
470  h(i,j,k) = max(h(i,j,k) + (w_int(i,j,k+1) - w_int(i,j,k)), &
471  min(h(i,j,k), gv%Angstrom_H))
472  enddo ; enddo
473  endif ; enddo
474 
475  if (cs%fldno > 0) deallocate(fld_mean_anom)
476 
477  endif
478 
479  do c=1,cs%num_col
480 ! c is an index for the next 3 lines but a multiplier for the rest of the loop
481 ! Therefore we use c as per C code and increment the index where necessary.
482  i = cs%col_i(c) ; j = cs%col_j(c)
483  damp = dt * cs%Iresttime_col(c)
484 
485  e(1) = 0.0 ; e0 = 0.0
486  do k=1,nz
487  e(k+1) = e(k) - h(i,j,k)*gv%H_to_Z
488  enddo
489  e_str = e(nz+1) / cs%Ref_eta(nz+1,c)
490 
491  if ( cs%bulkmixedlayer ) then
492  i1pdamp = 1.0 / (1.0 + damp)
493  if (associated(cs%Rcv_ml_ref)) &
494  rcv_ml(i,j) = i1pdamp * (rcv_ml(i,j) + cs%Rcv_ml_ref(c)*damp)
495  do k=1,nkmb
496  do m=1,cs%fldno
497  cs%var(m)%p(i,j,k) = i1pdamp * &
498  (cs%var(m)%p(i,j,k) + cs%Ref_val(m)%p(k,c)*damp)
499  enddo
500  enddo
501 
502  wpb = 0.0; wb = 0.0
503  do k=nz,nkmb+1,-1
504  if (gv%Rlay(k) > rcv_ml(i,j)) then
505  w = min((((e(k)-e0) - e_str*cs%Ref_eta(k,c)) * damp)*gv%Z_to_H, &
506  ((wb + h(i,j,k)) - gv%Angstrom_H))
507  wm = 0.5*(w-abs(w))
508  do m=1,cs%fldno
509  cs%var(m)%p(i,j,k) = (h(i,j,k)*cs%var(m)%p(i,j,k) + &
510  cs%Ref_val(m)%p(k,c)*(damp*h(i,j,k) + (wpb - wm))) / &
511  (h(i,j,k)*(1.0 + damp) + (wpb - wm))
512  enddo
513  else
514  do m=1,cs%fldno
515  cs%var(m)%p(i,j,k) = i1pdamp * &
516  (cs%var(m)%p(i,j,k) + cs%Ref_val(m)%p(k,c)*damp)
517  enddo
518  w = wb + (h(i,j,k) - gv%Angstrom_H)
519  wm = 0.5*(w-abs(w))
520  endif
521  eb(i,j,k) = eb(i,j,k) + wpb
522  ea(i,j,k) = ea(i,j,k) - wm
523  h(i,j,k) = h(i,j,k) + (wb - w)
524  wb = w
525  wpb = w - wm
526  enddo
527 
528  if (wb < 0) then
529  do k=nkmb,1,-1
530  w = min((wb + (h(i,j,k) - gv%Angstrom_H)),0.0)
531  h(i,j,k) = h(i,j,k) + (wb - w)
532  ea(i,j,k) = ea(i,j,k) - w
533  wb = w
534  enddo
535  else
536  w = wb
537  do k=gv%nkml,nkmb
538  eb(i,j,k) = eb(i,j,k) + w
539  enddo
540 
541  k = gv%nkml
542  h(i,j,k) = h(i,j,k) + w
543  do m=1,cs%fldno
544  cs%var(m)%p(i,j,k) = (cs%var(m)%p(i,j,k)*h(i,j,k) + &
545  cs%Ref_val(m)%p(k,c)*w) / (h(i,j,k) + w)
546  enddo
547  endif
548 
549  do k=1,nkmb
550  do m=1,cs%fldno
551  cs%var(m)%p(i,j,k) = i1pdamp * &
552  (cs%var(m)%p(i,j,k) + cs%Ref_val(m)%p(gv%nkml,c)*damp)
553  enddo
554  enddo
555 
556  else ! not BULKMIXEDLAYER
557 
558  wpb = 0.0
559  wb = 0.0
560  do k=nz,1,-1
561  w = min((((e(k)-e0) - e_str*cs%Ref_eta(k,c)) * damp)*gv%Z_to_H, &
562  ((wb + h(i,j,k)) - gv%Angstrom_H))
563  wm = 0.5*(w - abs(w))
564  do m=1,cs%fldno
565  cs%var(m)%p(i,j,k) = (h(i,j,k)*cs%var(m)%p(i,j,k) + &
566  cs%Ref_val(m)%p(k,c) * (damp*h(i,j,k) + (wpb - wm))) / &
567  (h(i,j,k)*(1.0 + damp) + (wpb - wm))
568  enddo
569  eb(i,j,k) = eb(i,j,k) + wpb
570  ea(i,j,k) = ea(i,j,k) - wm
571  h(i,j,k) = h(i,j,k) + (wb - w)
572  wb = w
573  wpb = w - wm
574  enddo
575 
576  endif ! end BULKMIXEDLAYER
577  enddo ! end of c loop
578 
579  if (associated(cs%diag)) then ; if (query_averaging_enabled(cs%diag)) then
580  if (cs%id_w_sponge > 0) then
581  idt = gv%H_to_m / dt ! Do any height unit conversion here for efficiency.
582  do k=1,nz+1 ; do j=js,je ; do i=is,ie
583  w_int(i,j,k) = w_int(i,j,k) * idt ! Scale values by clobbering array since it is local
584  enddo ; enddo ; enddo
585  call post_data(cs%id_w_sponge, w_int(:,:,:), cs%diag)
586  endif
587  endif ; endif
588 

References mom_spatial_means::global_i_mean(), and mom_error_handler::mom_error().

Referenced by mom_diabatic_driver::layered_diabatic().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ init_sponge_diags()

subroutine, public mom_sponge::init_sponge_diags ( type(time_type), intent(in), target  Time,
type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(unit_scale_type), intent(in)  US,
type(diag_ctrl), intent(inout), target  diag,
type(sponge_cs), pointer  CS 
)

This subroutine sets up diagnostics for the sponges. It is separate from initialize_sponge because it requires fields that are not readily available where initialize_sponge is called.

Parameters
[in]timeThe current model time
[in]gThe ocean's grid structure
[in]gvThe ocean's vertical grid structure
[in]usA dimensional unit scaling type
[in,out]diagA structure that is used to regulate diagnostic output
csA pointer to the control structure for this module that is set by a previous call to initialize_sponge.

Definition at line 194 of file MOM_sponge.F90.

194  type(time_type), target, intent(in) :: Time !< The current model time
195  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
196  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
197  type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
198  type(diag_ctrl), target, intent(inout) :: diag !< A structure that is used to regulate diagnostic output
199  type(sponge_CS), pointer :: CS !< A pointer to the control structure for this module that
200  !! is set by a previous call to initialize_sponge.
201 
202  if (.not.associated(cs)) return
203 
204  cs%diag => diag
205  cs%id_w_sponge = register_diag_field('ocean_model', 'w_sponge', diag%axesTi, &
206  time, 'The diapycnal motion due to the sponges', 'm s-1', conversion=us%s_to_T)
207 

Referenced by mom::initialize_mom().

Here is the caller graph for this function:

◆ initialize_sponge()

subroutine, public mom_sponge::initialize_sponge ( real, dimension(szi_(g),szj_(g)), intent(in)  Iresttime,
real, dimension(szi_(g),szj_(g),szk_(g)+1), intent(in)  int_height,
type(ocean_grid_type), intent(in)  G,
type(param_file_type), intent(in)  param_file,
type(sponge_cs), pointer  CS,
type(verticalgrid_type), intent(in)  GV,
real, dimension( g %jsd: g %jed), intent(in), optional  Iresttime_i_mean,
real, dimension( g %jsd: g %jed, g %ke+1), intent(in), optional  int_height_i_mean 
)

This subroutine determines the number of points which are within sponges in this computational domain. Only points that have positive values of Iresttime and which mask2dT indicates are ocean points are included in the sponges. It also stores the target interface heights.

Parameters
[in]gThe ocean's grid structure
[in]iresttimeThe inverse of the restoring time [s-1].
[in]int_heightThe interface heights to damp back toward [Z ~> m].
[in]param_fileA structure to parse for run-time parameters
csA pointer that is set to point to the control structure for this module
[in]gvThe ocean's vertical grid structure
[in]iresttime_i_meanThe inverse of the restoring time for
[in]int_height_i_meanThe interface heights toward which to

Definition at line 90 of file MOM_sponge.F90.

90  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
91  real, dimension(SZI_(G),SZJ_(G)), &
92  intent(in) :: Iresttime !< The inverse of the restoring time [s-1].
93  real, dimension(SZI_(G),SZJ_(G),SZK_(G)+1), &
94  intent(in) :: int_height !< The interface heights to damp back toward [Z ~> m].
95  type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
96  type(sponge_CS), pointer :: CS !< A pointer that is set to point to the control
97  !! structure for this module
98  type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure
99  real, dimension(SZJ_(G)), &
100  optional, intent(in) :: Iresttime_i_mean !< The inverse of the restoring time for
101  !! the zonal mean properties [s-1].
102  real, dimension(SZJ_(G),SZK_(G)+1), &
103  optional, intent(in) :: int_height_i_mean !< The interface heights toward which to
104  !! damp the zonal mean heights [Z ~> m].
105 
106 
107 ! This include declares and sets the variable "version".
108 #include "version_variable.h"
109  character(len=40) :: mdl = "MOM_sponge" ! This module's name.
110  logical :: use_sponge
111  integer :: i, j, k, col, total_sponge_cols
112 
113  if (associated(cs)) then
114  call mom_error(warning, "initialize_sponge called with an associated "// &
115  "control structure.")
116  return
117  endif
118 
119 ! Set default, read and log parameters
120  call log_version(param_file, mdl, version)
121  call get_param(param_file, mdl, "SPONGE", use_sponge, &
122  "If true, sponges may be applied anywhere in the domain. "//&
123  "The exact location and properties of those sponges are "//&
124  "specified from MOM_initialization.F90.", default=.false.)
125 
126  if (.not.use_sponge) return
127  allocate(cs)
128 
129  if (present(iresttime_i_mean) .neqv. present(int_height_i_mean)) &
130  call mom_error(fatal, "initialize_sponge: The optional arguments \n"//&
131  "Iresttime_i_mean and int_height_i_mean must both be present \n"//&
132  "if either one is.")
133 
134  cs%do_i_mean_sponge = present(iresttime_i_mean)
135 
136  cs%nz = g%ke
137 ! CS%isc = G%isc ; CS%iec = G%iec ; CS%jsc = G%jsc ; CS%jec = G%jec
138 ! CS%isd = G%isd ; CS%ied = G%ied ; CS%jsd = G%jsd ; CS%jed = G%jed
139  ! CS%bulkmixedlayer may be set later via a call to set_up_sponge_ML_density.
140  cs%bulkmixedlayer = .false.
141 
142  cs%num_col = 0 ; cs%fldno = 0
143  do j=g%jsc,g%jec ; do i=g%isc,g%iec
144  if ((iresttime(i,j)>0.0) .and. (g%mask2dT(i,j)>0)) &
145  cs%num_col = cs%num_col + 1
146  enddo ; enddo
147 
148  if (cs%num_col > 0) then
149 
150  allocate(cs%Iresttime_col(cs%num_col)) ; cs%Iresttime_col = 0.0
151  allocate(cs%col_i(cs%num_col)) ; cs%col_i = 0
152  allocate(cs%col_j(cs%num_col)) ; cs%col_j = 0
153 
154  col = 1
155  do j=g%jsc,g%jec ; do i=g%isc,g%iec
156  if ((iresttime(i,j)>0.0) .and. (g%mask2dT(i,j)>0)) then
157  cs%col_i(col) = i ; cs%col_j(col) = j
158  cs%Iresttime_col(col) = g%US%T_to_s*iresttime(i,j)
159  col = col +1
160  endif
161  enddo ; enddo
162 
163  allocate(cs%Ref_eta(cs%nz+1,cs%num_col))
164  do col=1,cs%num_col ; do k=1,cs%nz+1
165  cs%Ref_eta(k,col) = int_height(cs%col_i(col),cs%col_j(col),k)
166  enddo ; enddo
167 
168  endif
169 
170  if (cs%do_i_mean_sponge) then
171  allocate(cs%Iresttime_im(g%jsd:g%jed)) ; cs%Iresttime_im(:) = 0.0
172  allocate(cs%Ref_eta_im(g%jsd:g%jed,g%ke+1)) ; cs%Ref_eta_im(:,:) = 0.0
173 
174  do j=g%jsc,g%jec
175  cs%Iresttime_im(j) = g%US%T_to_s*iresttime_i_mean(j)
176  enddo
177  do k=1,cs%nz+1 ; do j=g%jsc,g%jec
178  cs%Ref_eta_im(j,k) = int_height_i_mean(j,k)
179  enddo ; enddo
180  endif
181 
182  total_sponge_cols = cs%num_col
183  call sum_across_pes(total_sponge_cols)
184 
185  call log_param(param_file, mdl, "!Total sponge columns", total_sponge_cols, &
186  "The total number of columns where sponges are applied.")
187 

References mom_error_handler::mom_error().

Referenced by bfb_initialization::bfb_initialize_sponges_southonly(), dome2d_initialization::dome2d_initialize_sponges(), dome_initialization::dome_initialize_sponges(), isomip_initialization::isomip_initialize_sponges(), and phillips_initialization::phillips_initialize_sponges().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_up_sponge_field()

subroutine, public mom_sponge::set_up_sponge_field ( real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in)  sp_val,
real, dimension( g %isd: g %ied, g %jsd: g %jed, g %ke), intent(in), target  f_ptr,
type(ocean_grid_type), intent(in)  G,
integer, intent(in)  nlay,
type(sponge_cs), pointer  CS,
real, dimension( g %jsd: g %jed, g %ke), intent(in), optional  sp_val_i_mean 
)

This subroutine stores the reference profile for the variable whose address is given by f_ptr. nlay is the number of layers in this variable.

Parameters
[in]gThe ocean's grid structure
[in]sp_valThe reference profiles of the quantity being registered.
[in]f_ptra pointer to the field which will be damped
[in]nlaythe number of layers in this quantity
csA pointer to the control structure for this module that is set by a previous call to initialize_sponge.
[in]sp_val_i_meanThe i-mean reference value for

Definition at line 214 of file MOM_sponge.F90.

214  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
215  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
216  intent(in) :: sp_val !< The reference profiles of the quantity being registered.
217  real, dimension(SZI_(G),SZJ_(G),SZK_(G)), &
218  target, intent(in) :: f_ptr !< a pointer to the field which will be damped
219  integer, intent(in) :: nlay !< the number of layers in this quantity
220  type(sponge_CS), pointer :: CS !< A pointer to the control structure for this module that
221  !! is set by a previous call to initialize_sponge.
222  real, dimension(SZJ_(G),SZK_(G)),&
223  optional, intent(in) :: sp_val_i_mean !< The i-mean reference value for
224  !! this field with i-mean sponges.
225 
226  integer :: j, k, col
227  character(len=256) :: mesg ! String for error messages
228 
229  if (.not.associated(cs)) return
230 
231  cs%fldno = cs%fldno + 1
232 
233  if (cs%fldno > max_fields_) then
234  write(mesg,'("Increase MAX_FIELDS_ to at least ",I3," in MOM_memory.h or decrease &
235  &the number of fields to be damped in the call to &
236  &initialize_sponge." )') cs%fldno
237  call mom_error(fatal,"set_up_sponge_field: "//mesg)
238  endif
239 
240  allocate(cs%Ref_val(cs%fldno)%p(cs%nz,cs%num_col))
241  cs%Ref_val(cs%fldno)%p(:,:) = 0.0
242  do col=1,cs%num_col
243  do k=1,nlay
244  cs%Ref_val(cs%fldno)%p(k,col) = sp_val(cs%col_i(col),cs%col_j(col),k)
245  enddo
246  do k=nlay+1,cs%nz
247  cs%Ref_val(cs%fldno)%p(k,col) = 0.0
248  enddo
249  enddo
250 
251  cs%var(cs%fldno)%p => f_ptr
252 
253  if (nlay/=cs%nz) then
254  write(mesg,'("Danger: Sponge reference fields require nz (",I3,") layers.&
255  & A field with ",I3," layers was passed to set_up_sponge_field.")') &
256  cs%nz, nlay
257  if (is_root_pe()) call mom_error(warning, "set_up_sponge_field: "//mesg)
258  endif
259 
260  if (cs%do_i_mean_sponge) then
261  if (.not.present(sp_val_i_mean)) call mom_error(fatal, &
262  "set_up_sponge_field: sp_val_i_mean must be present with i-mean sponges.")
263 
264  allocate(cs%Ref_val_im(cs%fldno)%p(cs%jsd:cs%jed,cs%nz))
265  cs%Ref_val(cs%fldno)%p(:,:) = 0.0
266  do k=1,cs%nz ; do j=cs%jsc,cs%jec
267  cs%Ref_val_im(cs%fldno)%p(j,k) = sp_val_i_mean(j,k)
268  enddo ; enddo
269  endif
270 

References mom_error_handler::is_root_pe(), and mom_error_handler::mom_error().

Referenced by dome_initialization::dome_initialize_sponges(), dome_tracer::initialize_dome_tracer(), rgc_tracer::initialize_rgc_tracer(), isomip_initialization::isomip_initialize_sponges(), and user_tracer_example::user_initialize_tracer().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ set_up_sponge_ml_density()

subroutine, public mom_sponge::set_up_sponge_ml_density ( real, dimension(szi_(g),szj_(g)), intent(in)  sp_val,
type(ocean_grid_type), intent(in)  G,
type(sponge_cs), pointer  CS,
real, dimension(szj_(g)), intent(in), optional  sp_val_i_mean 
)

This subroutine stores the reference value for mixed layer density. It is handled differently from other values because it is only used in determining which layers can be inflated.

Parameters
[in]gThe ocean's grid structure
[in]sp_valThe reference values of the mixed layer density [R ~> kg m-3]
csA pointer to the control structure for this module that is set by a previous call to initialize_sponge.
[in]sp_val_i_meanthe reference values of the zonal mean mixed

Definition at line 277 of file MOM_sponge.F90.

277  type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
278  real, dimension(SZI_(G),SZJ_(G)), &
279  intent(in) :: sp_val !< The reference values of the mixed layer density [R ~> kg m-3]
280  type(sponge_CS), pointer :: CS !< A pointer to the control structure for this module that is
281  !! set by a previous call to initialize_sponge.
282  real, dimension(SZJ_(G)), &
283  optional, intent(in) :: sp_val_i_mean !< the reference values of the zonal mean mixed
284  !! layer density [R ~> kg m-3], for use if Iresttime_i_mean > 0.
285 ! This subroutine stores the reference value for mixed layer density. It is
286 ! handled differently from other values because it is only used in determining
287 ! which layers can be inflated.
288 
289 ! Arguments: sp_val - The reference values of the mixed layer density.
290 ! (in/out) CS - A pointer to the control structure for this module that is
291 ! set by a previous call to initialize_sponge.
292 
293  integer :: j, col
294  character(len=256) :: mesg ! String for error messages
295 
296  if (.not.associated(cs)) return
297 
298  if (associated(cs%Rcv_ml_ref)) then
299  call mom_error(fatal, "set_up_sponge_ML_density appears to have been "//&
300  "called twice.")
301  endif
302 
303  cs%bulkmixedlayer = .true.
304  allocate(cs%Rcv_ml_ref(cs%num_col)) ; cs%Rcv_ml_ref(:) = 0.0
305  do col=1,cs%num_col
306  cs%Rcv_ml_ref(col) = sp_val(cs%col_i(col),cs%col_j(col))
307  enddo
308 
309  if (cs%do_i_mean_sponge) then
310  if (.not.present(sp_val_i_mean)) call mom_error(fatal, &
311  "set_up_sponge_field: sp_val_i_mean must be present with i-mean sponges.")
312 
313  allocate(cs%Rcv_ml_ref_im(cs%jsd:cs%jed)) ; cs%Rcv_ml_ref_im(:) = 0.0
314  do j=cs%jsc,cs%jec
315  cs%Rcv_ml_ref_im(j) = sp_val_i_mean(j)
316  enddo
317  endif
318 

References mom_error_handler::mom_error().

Referenced by rgc_initialization::rgc_initialize_sponges().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ sponge_end()

subroutine, public mom_sponge::sponge_end ( type(sponge_cs), pointer  CS)

This call deallocates any memory in the sponge control structure.

Parameters
csA pointer to the control structure for this module that is set by a previous call to initialize_sponge.

Definition at line 593 of file MOM_sponge.F90.

593  type(sponge_CS), pointer :: CS !< A pointer to the control structure for this module
594  !! that is set by a previous call to initialize_sponge.
595  integer :: m
596 
597  if (.not.associated(cs)) return
598 
599  if (associated(cs%col_i)) deallocate(cs%col_i)
600  if (associated(cs%col_j)) deallocate(cs%col_j)
601 
602  if (associated(cs%Iresttime_col)) deallocate(cs%Iresttime_col)
603  if (associated(cs%Rcv_ml_ref)) deallocate(cs%Rcv_ml_ref)
604  if (associated(cs%Ref_eta)) deallocate(cs%Ref_eta)
605 
606  if (associated(cs%Iresttime_im)) deallocate(cs%Iresttime_im)
607  if (associated(cs%Rcv_ml_ref_im)) deallocate(cs%Rcv_ml_ref_im)
608  if (associated(cs%Ref_eta_im)) deallocate(cs%Ref_eta_im)
609 
610  do m=1,cs%fldno
611  if (associated(cs%Ref_val(cs%fldno)%p)) deallocate(cs%Ref_val(cs%fldno)%p)
612  if (associated(cs%Ref_val_im(cs%fldno)%p)) &
613  deallocate(cs%Ref_val_im(cs%fldno)%p)
614  enddo
615 
616  deallocate(cs)
617