MOM6
dense_water_initialization Module Reference

Detailed Description

Initialization routines for the dense water formation and overflow experiment.

Functions/Subroutines

subroutine, public dense_water_initialize_topography (D, G, param_file, max_depth)
 Initialize the topography field for the dense water experiment. More...
 
subroutine, public dense_water_initialize_ts (G, GV, param_file, eqn_of_state, T, S, h, just_read_params)
 Initialize the temperature and salinity for the dense water experiment. More...
 
subroutine, public dense_water_initialize_sponges (G, GV, tv, param_file, use_ALE, CSp, ACSp)
 Initialize the restoring sponges for the dense water experiment. More...
 

Variables

character(len=40) mdl = "dense_water_initialization"
 
real, parameter default_sill = 0.2
 Default depth of the sill [nondim]. More...
 
real, parameter default_shelf = 0.4
 Default depth of the shelf [nondim]. More...
 
real, parameter default_mld = 0.25
 Default depth of the mixed layer [nondim]. More...
 

Function/Subroutine Documentation

◆ dense_water_initialize_sponges()

subroutine, public dense_water_initialization::dense_water_initialize_sponges ( type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(thermo_var_ptrs), intent(in)  tv,
type(param_file_type), intent(in)  param_file,
logical, intent(in)  use_ALE,
type(sponge_cs), pointer  CSp,
type(ale_sponge_cs), pointer  ACSp 
)

Initialize the restoring sponges for the dense water experiment.

Parameters
[in]gHorizontal grid control structure
[in]gvVertical grid control structure
[in]tvThermodynamic variables
[in]param_fileParameter file structure
[in]use_aleALE flag
cspLayered sponge control structure pointer
acspALE sponge control structure pointer

Definition at line 147 of file dense_water_initialization.F90.

References default_mld, default_sill, mom_ale_sponge::initialize_ale_sponge(), mdl, mom_error_handler::mom_error(), and mom_ale_sponge::set_up_ale_sponge_field().

Referenced by mom_state_initialization::mom_initialize_state().

147  type(ocean_grid_type), intent(in) :: g !< Horizontal grid control structure
148  type(verticalgrid_type), intent(in) :: gv !< Vertical grid control structure
149  type(thermo_var_ptrs), intent(in) :: tv !< Thermodynamic variables
150  type(param_file_type), intent(in) :: param_file !< Parameter file structure
151  logical, intent(in) :: use_ale !< ALE flag
152  type(sponge_cs), pointer :: csp !< Layered sponge control structure pointer
153  type(ale_sponge_cs), pointer :: acsp !< ALE sponge control structure pointer
154 
155  ! Local variables
156  real :: west_sponge_time_scale, west_sponge_width
157  real :: east_sponge_time_scale, east_sponge_width
158 
159  real, dimension(SZI_(G),SZJ_(G)) :: idamp ! inverse damping timescale
160  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)) :: h, t, s ! sponge thicknesses, temp and salt
161  real, dimension(SZK_(GV)+1) :: e0, eta1d ! interface positions for ALE sponge
162 
163  integer :: i, j, k, nz
164  real :: x, zi, zmid, dist
165  real :: mld, s_ref, s_range, s_dense, t_ref, sill_height
166 
167  nz = gv%ke
168 
169  call get_param(param_file, mdl, "DENSE_WATER_WEST_SPONGE_TIME_SCALE", west_sponge_time_scale, &
170  "The time scale on the west (outflow) of the domain for restoring. If zero, the sponge is disabled.", &
171  units="s", default=0.)
172  call get_param(param_file, mdl, "DENSE_WATER_WEST_SPONGE_WIDTH", west_sponge_width, &
173  "The fraction of the domain in which the western (outflow) sponge is active.", &
174  units="nondim", default=0.1)
175  call get_param(param_file, mdl, "DENSE_WATER_EAST_SPONGE_TIME_SCALE", east_sponge_time_scale, &
176  "The time scale on the east (outflow) of the domain for restoring. If zero, the sponge is disabled.", &
177  units="s", default=0.)
178  call get_param(param_file, mdl, "DENSE_WATER_EAST_SPONGE_WIDTH", east_sponge_width, &
179  "The fraction of the domain in which the eastern (outflow) sponge is active.", &
180  units="nondim", default=0.1)
181 
182  call get_param(param_file, mdl, "DENSE_WATER_EAST_SPONGE_SALT", s_dense, &
183  "Salt anomaly of the dense water being formed in the overflow region.", &
184  units="1e-3", default=4.0)
185 
186  call get_param(param_file, mdl, "DENSE_WATER_MLD", mld, default=default_mld, do_not_log=.true.)
187  call get_param(param_file, mdl, "DENSE_WATER_SILL_HEIGHT", sill_height, default=default_sill, do_not_log=.true.)
188 
189  call get_param(param_file, mdl, "S_REF", s_ref, do_not_log=.true.)
190  call get_param(param_file, mdl, "S_RANGE", s_range, do_not_log=.true.)
191  call get_param(param_file, mdl, "T_REF", t_ref, do_not_log=.true.)
192 
193  ! no active sponges
194  if (west_sponge_time_scale <= 0. .and. east_sponge_time_scale <= 0.) return
195 
196  ! everywhere is initially unsponged
197  idamp(:,:) = 0.0
198 
199  do j = g%jsc, g%jec
200  do i = g%isc,g%iec
201  if (g%mask2dT(i,j) > 0.) then
202  ! nondimensional x position
203  x = (g%geoLonT(i,j) - g%west_lon) / g%len_lon
204 
205  if (west_sponge_time_scale > 0. .and. x < west_sponge_width) then
206  dist = 1. - x / west_sponge_width
207  ! scale restoring by depth into sponge
208  idamp(i,j) = 1. / west_sponge_time_scale * max(0., min(1., dist))
209  else if (east_sponge_time_scale > 0. .and. x > (1. - east_sponge_width)) then
210  dist = 1. - (1. - x) / east_sponge_width
211  idamp(i,j) = 1. / east_sponge_time_scale * max(0., min(1., dist))
212  endif
213  endif
214  enddo
215  enddo
216 
217  if (use_ale) then
218  ! construct a uniform grid for the sponge
219  do k = 1,nz
220  e0(k) = -gv%max_depth * (real(k - 1) / real(nz))
221  enddo
222  e0(nz+1) = -gv%max_depth
223 
224  do j = g%jsc,g%jec
225  do i = g%isc,g%iec
226  eta1d(nz+1) = -g%bathyT(i,j)
227  do k = nz,1,-1
228  eta1d(k) = e0(k)
229 
230  if (eta1d(k) < (eta1d(k+1) + gv%Angstrom_z)) then
231  ! is this layer vanished?
232  eta1d(k) = eta1d(k+1) + gv%Angstrom_z
233  h(i,j,k) = gv%Angstrom_z
234  else
235  h(i,j,k) = eta1d(k) - eta1d(k+1)
236  endif
237  enddo
238  enddo
239  enddo
240 
241  call initialize_ale_sponge(idamp, h, nz, g, param_file, acsp)
242 
243  ! construct temperature and salinity for the sponge
244  ! start with initial condition
245  t(:,:,:) = t_ref
246  s(:,:,:) = s_ref
247 
248  do j = g%jsc,g%jec
249  do i = g%isc,g%iec
250  zi = 0.
251  x = (g%geoLonT(i,j) - g%west_lon) / g%len_lon
252  do k = 1,nz
253  ! nondimensional middle of layer
254  zmid = zi + 0.5 * h(i,j,k) / gv%max_depth
255 
256  if (x > (1. - east_sponge_width)) then
257  !if (zmid >= 0.9 * sill_height) &
258  s(i,j,k) = s_ref + s_dense
259  else
260  ! linear between bottom of mixed layer and bottom
261  if (zmid >= mld) &
262  s(i,j,k) = s_ref + s_range * (zmid - mld) / (1.0 - mld)
263  endif
264 
265  zi = zi + h(i,j,k) / gv%max_depth
266  enddo
267  enddo
268  enddo
269 
270  if (associated(tv%T)) call set_up_ale_sponge_field(t, g, tv%T, acsp)
271  if (associated(tv%S)) call set_up_ale_sponge_field(s, g, tv%S, acsp)
272  else
273  call mom_error(fatal, "dense_water_initialize_sponges: trying to use non ALE sponge")
274  endif
Here is the call graph for this function:
Here is the caller graph for this function:

◆ dense_water_initialize_topography()

subroutine, public dense_water_initialization::dense_water_initialize_topography ( real, dimension( g %isd: g %ied, g %jsd: g %jed), intent(out)  D,
type(dyn_horgrid_type), intent(in)  G,
type(param_file_type), intent(in)  param_file,
real, intent(in)  max_depth 
)

Initialize the topography field for the dense water experiment.

Parameters
[in]gGrid control structure
[out]dOutput topography field
[in]param_fileParameter file structure
[in]max_depthMaximum depth of the model

Definition at line 33 of file dense_water_initialization.F90.

References default_shelf, default_sill, and mdl.

Referenced by mom_fixed_initialization::mom_initialize_topography().

33  type(dyn_horgrid_type), intent(in) :: g !< Grid control structure
34  real, dimension(SZI_(G),SZJ_(G)), intent(out) :: d !< Output topography field
35  type(param_file_type), intent(in) :: param_file !< Parameter file structure
36  real, intent(in) :: max_depth !< Maximum depth of the model
37 
38  real, dimension(5) :: domain_params ! nondimensional widths of all domain sections
39  real :: sill_frac, shelf_frac
40  integer :: i, j
41  real :: x
42 
43  call get_param(param_file, mdl, "DENSE_WATER_DOMAIN_PARAMS", domain_params, &
44  "Fractional widths of all the domain sections for the dense water experiment.\n"//&
45  "As a 5-element vector:\n"//&
46  " - open ocean, the section at maximum depth\n"//&
47  " - downslope, the downward overflow slope\n"//&
48  " - sill separating downslope from upslope\n"//&
49  " - upslope, the upward slope accumulating dense water\n"//&
50  " - the shelf in the dense formation region.", &
51  units="nondim", fail_if_missing=.true.)
52  call get_param(param_file, mdl, "DENSE_WATER_SILL_DEPTH", sill_frac, &
53  "Depth of the sill separating downslope from upslope, as fraction of basin depth.", &
54  units="nondim", default=default_sill)
55  call get_param(param_file, mdl, "DENSE_WATER_SHELF_DEPTH", shelf_frac, &
56  "Depth of the shelf region accumulating dense water for overflow, as fraction of basin depth.", &
57  units="nondim", default=default_shelf)
58 
59  do i = 2, 5
60  ! turn widths into positions
61  domain_params(i) = domain_params(i-1) + domain_params(i)
62  enddo
63 
64  do i = g%isc,g%iec
65  do j = g%jsc,g%jec
66  ! compute normalised zonal coordinate
67  x = (g%geoLonT(i,j) - g%west_lon) / g%len_lon
68 
69  if (x <= domain_params(1)) then
70  ! open ocean region
71  d(i,j) = max_depth
72  else if (x <= domain_params(2)) then
73  ! downslope region, linear
74  d(i,j) = max_depth - (1.0 - sill_frac) * max_depth * &
75  (x - domain_params(1)) / (domain_params(2) - domain_params(1))
76  else if (x <= domain_params(3)) then
77  ! sill region
78  d(i,j) = sill_frac * max_depth
79  else if (x <= domain_params(4)) then
80  ! upslope region
81  d(i,j) = sill_frac * max_depth + (shelf_frac - sill_frac) * max_depth * &
82  (x - domain_params(3)) / (domain_params(4) - domain_params(3))
83  else
84  ! shelf region
85  d(i,j) = shelf_frac * max_depth
86  endif
87  enddo
88  enddo
Here is the caller graph for this function:

◆ dense_water_initialize_ts()

subroutine, public dense_water_initialization::dense_water_initialize_ts ( type(ocean_grid_type), intent(in)  G,
type(verticalgrid_type), intent(in)  GV,
type(param_file_type), intent(in)  param_file,
type(eos_type), pointer  eqn_of_state,
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(out)  T,
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(out)  S,
real, dimension( g %isd: g %ied, g %jsd: g %jed, gv %ke), intent(in)  h,
logical, intent(in), optional  just_read_params 
)

Initialize the temperature and salinity for the dense water experiment.

Parameters
[in]gHorizontal grid control structure
[in]gvVertical grid control structure
[in]param_fileParameter file structure
eqn_of_stateEOS structure
[out]sOutput state
[in]hLayer thicknesses
[in]just_read_paramsIf present and true, this call will only read parameters without changing h.

Definition at line 93 of file dense_water_initialization.F90.

References default_mld, and mdl.

93  type(ocean_grid_type), intent(in) :: g !< Horizontal grid control structure
94  type(verticalgrid_type), intent(in) :: gv !< Vertical grid control structure
95  type(param_file_type), intent(in) :: param_file !< Parameter file structure
96  type(eos_type), pointer :: eqn_of_state !< EOS structure
97  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(out) :: t, s !< Output state
98  real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: h !< Layer thicknesses
99  logical, optional, intent(in) :: just_read_params !< If present and true, this call will
100  !! only read parameters without changing h.
101 
102  real :: mld, s_ref, s_range, t_ref
103  real :: zi, zmid
104  logical :: just_read ! If true, just read parameters but set nothing.
105  integer :: i, j, k, nz
106 
107  nz = gv%ke
108 
109  just_read = .false. ; if (present(just_read_params)) just_read = just_read_params
110 
111  call get_param(param_file, mdl, "DENSE_WATER_MLD", mld, &
112  "Depth of unstratified mixed layer as a fraction of the water column.", &
113  units="nondim", default=default_mld, do_not_log=just_read)
114 
115  call get_param(param_file, mdl, "S_REF", s_ref, do_not_log=.true.)
116  call get_param(param_file, mdl, "S_RANGE", s_range, do_not_log=.true.)
117  call get_param(param_file, mdl, "T_REF", t_ref, do_not_log=.true.)
118 
119  if (just_read) return ! All run-time parameters have been read, so return.
120 
121  ! uniform temperature everywhere
122  t(:,:,:) = t_ref
123 
124  do j = g%jsc,g%jec
125  do i = g%isc,g%iec
126  zi = 0.
127  do k = 1,nz
128  ! nondimensional middle of layer
129  zmid = zi + 0.5 * h(i,j,k) / g%max_depth
130 
131  if (zmid < mld) then
132  ! use reference salinity in the mixed layer
133  s(i,j,k) = s_ref
134  else
135  ! linear between bottom of mixed layer and bottom
136  s(i,j,k) = s_ref + s_range * (zmid - mld) / (1.0 - mld)
137  endif
138 
139  zi = zi + h(i,j,k) / g%max_depth
140  enddo
141  enddo
142  enddo

Variable Documentation

◆ default_mld

real, parameter dense_water_initialization::default_mld = 0.25
private

Default depth of the mixed layer [nondim].

Definition at line 27 of file dense_water_initialization.F90.

Referenced by dense_water_initialize_sponges(), and dense_water_initialize_ts().

27 real, parameter :: default_mld = 0.25 !< Default depth of the mixed layer [nondim]

◆ default_shelf

real, parameter dense_water_initialization::default_shelf = 0.4
private

Default depth of the shelf [nondim].

Definition at line 26 of file dense_water_initialization.F90.

Referenced by dense_water_initialize_topography().

26 real, parameter :: default_shelf = 0.4 !< Default depth of the shelf [nondim]

◆ default_sill

real, parameter dense_water_initialization::default_sill = 0.2
private

Default depth of the sill [nondim].

Definition at line 25 of file dense_water_initialization.F90.

Referenced by dense_water_initialize_sponges(), and dense_water_initialize_topography().

25 real, parameter :: default_sill = 0.2 !< Default depth of the sill [nondim]

◆ mdl

character(len=40) dense_water_initialization::mdl = "dense_water_initialization"
private

Definition at line 23 of file dense_water_initialization.F90.

Referenced by dense_water_initialize_sponges(), dense_water_initialize_topography(), and dense_water_initialize_ts().

23 character(len=40) :: mdl = "dense_water_initialization"