15 implicit none ;
private 17 #include <MOM_memory.h> 23 character(len=40) ::
mdl =
"dense_water_initialization" 34 real,
dimension(SZI_(G),SZJ_(G)),
intent(out) :: D
36 real,
intent(in) :: max_depth
38 real,
dimension(5) :: domain_params
39 real :: sill_frac, shelf_frac
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.", &
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.", &
61 domain_params(i) = domain_params(i-1) + domain_params(i)
67 x = (g%geoLonT(i,j) - g%west_lon) / g%len_lon
69 if (x <= domain_params(1))
then 72 else if (x <= domain_params(2))
then 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 78 d(i,j) = sill_frac * max_depth
79 else if (x <= domain_params(4))
then 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))
85 d(i,j) = shelf_frac * max_depth
96 type(
eos_type),
pointer :: eqn_of_state
97 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)),
intent(out) :: T, S
98 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)),
intent(in) :: h
99 logical,
optional,
intent(in) :: just_read_params
102 real :: mld, S_ref, S_range, T_ref
105 integer :: i, j, k, nz
109 just_read = .false. ;
if (
present(just_read_params)) just_read = just_read_params
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)
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.)
119 if (just_read)
return 129 zmid = zi + 0.5 * h(i,j,k) / g%max_depth
136 s(i,j,k) = s_ref + s_range * (zmid - mld) / (1.0 - mld)
139 zi = zi + h(i,j,k) / g%max_depth
151 logical,
intent(in) :: use_ALE
156 real :: west_sponge_time_scale, west_sponge_width
157 real :: east_sponge_time_scale, east_sponge_width
159 real,
dimension(SZI_(G),SZJ_(G)) :: Idamp
160 real,
dimension(SZI_(G),SZJ_(G),SZK_(GV)) :: h, T, S
161 real,
dimension(SZK_(GV)+1) :: e0, eta1D
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
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)
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)
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.)
194 if (west_sponge_time_scale <= 0. .and. east_sponge_time_scale <= 0.)
return 201 if (g%mask2dT(i,j) > 0.)
then 203 x = (g%geoLonT(i,j) - g%west_lon) / g%len_lon
205 if (west_sponge_time_scale > 0. .and. x < west_sponge_width)
then 206 dist = 1. - x / west_sponge_width
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))
220 e0(k) = -gv%max_depth * (
real(k - 1) /
real(nz))
222 e0(nz+1) = -gv%max_depth
226 eta1d(nz+1) = -g%bathyT(i,j)
230 if (eta1d(k) < (eta1d(k+1) + gv%Angstrom_z))
then 232 eta1d(k) = eta1d(k+1) + gv%Angstrom_z
233 h(i,j,k) = gv%Angstrom_z
235 h(i,j,k) = eta1d(k) - eta1d(k+1)
251 x = (g%geoLonT(i,j) - g%west_lon) / g%len_lon
254 zmid = zi + 0.5 * h(i,j,k) / gv%max_depth
256 if (x > (1. - east_sponge_width))
then 258 s(i,j,k) = s_ref + s_dense
262 s(i,j,k) = s_ref + s_range * (zmid - mld) / (1.0 - mld)
265 zi = zi + h(i,j,k) / gv%max_depth
273 call mom_error(fatal,
"dense_water_initialize_sponges: trying to use non ALE sponge")
real, parameter default_mld
Default depth of the mixed layer [nondim].
subroutine, public dense_water_initialize_topography(D, G, param_file, max_depth)
Initialize the topography field for the dense water experiment.
Ocean grid type. See mom_grid for details.
Provides the ocean grid type.
subroutine, public set_up_ale_sponge_field(sp_val, G, f_ptr, CS)
This subroutine stores the reference profile at h points for the variable.
This module contains the routines used to apply sponge layers when using the ALE mode. Applying sponges requires the following: (1) initialize_ALE_sponge (2) set_up_ALE_sponge_field (tracers) and set_up_ALE_sponge_vel_field (vel) (3) apply_ALE_sponge (4) init_ALE_sponge_diags (not being used for now) (5) ALE_sponge_end (not being used for now)
SPONGE control structure.
real, parameter default_sill
Default depth of the sill [nondim].
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.
subroutine, public initialize_ale_sponge(Iresttime, data_h, nz_data, G, param_file, CS)
This subroutine determines the number of points which are within.
Provides subroutines for quantities specific to the equation of state.
The thermo_var_ptrs structure contains pointers to an assortment of thermodynamic fields that may be ...
Initialization routines for the dense water formation and overflow experiment.
subroutine, public dense_water_initialize_sponges(G, GV, tv, param_file, use_ALE, CSp, ACSp)
Initialize the restoring sponges for the dense water experiment.
subroutine, public mom_error(level, message, all_print)
A control structure for the equation of state.
real, parameter default_shelf
Default depth of the shelf [nondim].