17 use cvmix_shear,
only : cvmix_init_shear, cvmix_coeffs_shear
19 implicit none ;
private
21 #include <MOM_memory.h>
23 public calculate_cvmix_shear, cvmix_shear_init, cvmix_shear_is_used, cvmix_shear_end
39 real,
allocatable,
dimension(:,:,:) :: n2
40 real,
allocatable,
dimension(:,:,:) :: s2
41 real,
allocatable,
dimension(:,:,:) :: ri_grad
42 real,
allocatable,
dimension(:,:,:) :: ri_grad_smooth
44 character(10) :: mix_scheme
48 integer :: id_n2 = -1, id_s2 = -1, id_ri_grad = -1, id_kv = -1, id_kd = -1
49 integer :: id_ri_grad_smooth = -1
54 character(len=40) :: mdl =
"MOM_CVMix_shear"
59 subroutine calculate_cvmix_shear(u_H, v_H, h, tv, kd, kv, G, GV, US, CS )
63 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: u_h
64 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: v_h
66 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: h
68 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1),
intent(out) :: kd
70 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1),
intent(out) :: kv
75 integer :: i, j, k, kk, km1
84 real,
dimension(2*(G%ke)) :: pres_1d
85 real,
dimension(2*(G%ke)) :: temp_1d
86 real,
dimension(2*(G%ke)) :: salt_1d
87 real,
dimension(2*(G%ke)) :: rho_1d
88 real,
dimension(G%ke+1) :: ri_grad
89 real,
dimension(G%ke+1) :: kvisc
90 real,
dimension(G%ke+1) :: kdiff
94 gorho = us%L_to_Z**2 * gv%g_Earth / gv%Rho0
95 epsln = 1.e-10 * gv%m_to_H
101 if (g%mask2dT(i,j)==0.) cycle
104 pref = 0. ;
if (
associated(tv%p_surf)) pref = tv%p_surf(i,j)
114 temp_1d(kk+1) = tv%T(i,j,k)
115 temp_1d(kk+2) = tv%T(i,j,km1)
116 salt_1d(kk+1) = tv%S(i,j,k)
117 salt_1d(kk+2) = tv%S(i,j,km1)
121 pref = pref + (gv%g_Earth * gv%H_to_RZ) * h(i,j,k)
132 du = u_h(i,j,k) - u_h(i,j,km1)
133 dv = v_h(i,j,k) - v_h(i,j,km1)
134 drho = gorho * (rho_1d(kk+1) - rho_1d(kk+2))
135 dz = (0.5*(h(i,j,km1) + h(i,j,k))+gv%H_subroundoff)*gv%H_to_Z
137 s2 = us%L_to_Z**2*(du*du+dv*dv)/(dz*dz)
138 ri_grad(k) = max(0., n2) / max(s2, 1.e-10*us%T_to_s**2)
141 if (cs%id_N2 > 0) cs%N2(i,j,k) = n2
142 if (cs%id_S2 > 0) cs%S2(i,j,k) = s2
146 ri_grad(g%ke+1) = ri_grad(g%ke)
148 if (cs%id_ri_grad > 0) cs%ri_grad(i,j,:) = ri_grad(:)
150 if (cs%smooth_ri)
then
153 if (h(i,j,k) <= epsln) ri_grad(k) = ri_grad(k-1)
156 ri_grad(g%ke+1) = ri_grad(g%ke)
159 dummy = 0.25 * ri_grad(2)
160 ri_grad(g%ke+1) = ri_grad(g%ke)
162 ri_grad(k) = dummy + 0.5 * ri_grad(k) + 0.25 * ri_grad(k+1)
163 dummy = 0.25 * ri_grad(k)
166 if (cs%id_ri_grad_smooth > 0) cs%ri_grad_smooth(i,j,:) = ri_grad(:)
170 kvisc(k) = us%Z2_T_to_m2_s * kv(i,j,k)
171 kdiff(k) = us%Z2_T_to_m2_s * kd(i,j,k)
175 call cvmix_coeffs_shear(mdiff_out=kvisc(:), &
176 tdiff_out=kdiff(:), &
181 kv(i,j,k) = us%m2_s_to_Z2_T * kvisc(k)
182 kd(i,j,k) = us%m2_s_to_Z2_T * kdiff(k)
188 if (cs%id_kd > 0)
call post_data(cs%id_kd, kd, cs%diag)
189 if (cs%id_kv > 0)
call post_data(cs%id_kv, kv, cs%diag)
190 if (cs%id_N2 > 0)
call post_data(cs%id_N2, cs%N2, cs%diag)
191 if (cs%id_S2 > 0)
call post_data(cs%id_S2, cs%S2, cs%diag)
192 if (cs%id_ri_grad > 0)
call post_data(cs%id_ri_grad, cs%ri_grad, cs%diag)
193 if (cs%id_ri_grad_smooth > 0)
call post_data(cs%id_ri_grad_smooth ,cs%ri_grad_smooth, cs%diag)
195 end subroutine calculate_cvmix_shear
203 logical function cvmix_shear_init(Time, G, GV, US, param_file, diag, CS)
204 type(time_type),
intent(in) :: time
209 type(
diag_ctrl),
target,
intent(inout) :: diag
212 integer :: numbertrue=0
215 #include "version_variable.h"
217 if (
associated(cs))
then
218 call mom_error(warning,
"CVMix_shear_init called with an associated "// &
219 "control structure.")
225 call get_param(param_file, mdl,
"USE_LMD94", cs%use_LMD94, default=.false., do_not_log=.true.)
226 call get_param(param_file, mdl,
"USE_PP81", cs%use_PP81, default=.false., do_not_log=.true.)
228 "Parameterization of shear-driven turbulence via CVMix (various options)", &
229 all_default=.not.(cs%use_PP81.or.cs%use_LMD94))
230 call get_param(param_file, mdl,
"USE_LMD94", cs%use_LMD94, &
231 "If true, use the Large-McWilliams-Doney (JGR 1994) "//&
232 "shear mixing parameterization.", default=.false.)
233 if (cs%use_LMD94)
then
234 numbertrue=numbertrue + 1
237 call get_param(param_file, mdl,
"USE_PP81", cs%use_PP81, &
238 "If true, use the Pacanowski and Philander (JPO 1981) "//&
239 "shear mixing parameterization.", default=.false.)
240 if (cs%use_PP81)
then
241 numbertrue = numbertrue + 1
244 use_jhl=kappa_shear_is_used(param_file)
245 if (use_jhl) numbertrue = numbertrue + 1
248 if ((numbertrue) > 1)
then
249 call mom_error(fatal,
'MOM_CVMix_shear_init: '// &
250 'Multiple shear driven internal mixing schemes selected,'//&
251 ' please disable all but one scheme to proceed.')
253 cvmix_shear_init=(cs%use_PP81.or.cs%use_LMD94)
256 if (.not. cvmix_shear_init)
return
257 call get_param(param_file, mdl,
"NU_ZERO", cs%Nu_Zero, &
258 "Leading coefficient in KPP shear mixing.", &
259 units=
"nondim", default=5.e-3)
260 call get_param(param_file, mdl,
"RI_ZERO", cs%Ri_Zero, &
261 "Critical Richardson for KPP shear mixing, "// &
262 "NOTE this the internal mixing and this is "// &
263 "not for setting the boundary layer depth." &
264 ,units=
"nondim", default=0.8)
265 call get_param(param_file, mdl,
"KPP_EXP", cs%KPP_exp, &
266 "Exponent of unitless factor of diffusivities, "// &
267 "for KPP internal shear mixing scheme." &
268 ,units=
"nondim", default=3.0)
269 call get_param(param_file, mdl,
"SMOOTH_RI", cs%smooth_ri, &
270 "If true, vertically smooth the Richardson "// &
271 "number by applying a 1-2-1 filter once.", &
273 call cvmix_init_shear(mix_scheme=cs%Mix_Scheme, &
274 kpp_nu_zero=cs%Nu_Zero, &
275 kpp_ri_zero=cs%Ri_zero, &
281 cs%id_N2 = register_diag_field(
'ocean_model',
'N2_shear', diag%axesTi, time, &
282 'Square of Brunt-Vaisala frequency used by MOM_CVMix_shear module',
'1/s2', conversion=us%s_to_T**2)
283 if (cs%id_N2 > 0)
then
284 allocate( cs%N2( szi_(g), szj_(g), szk_(g)+1 ) ) ; cs%N2(:,:,:) = 0.
287 cs%id_S2 = register_diag_field(
'ocean_model',
'S2_shear', diag%axesTi, time, &
288 'Square of vertical shear used by MOM_CVMix_shear module',
'1/s2', conversion=us%s_to_T**2)
289 if (cs%id_S2 > 0)
then
290 allocate( cs%S2( szi_(g), szj_(g), szk_(g)+1 ) ) ; cs%S2(:,:,:) = 0.
293 cs%id_ri_grad = register_diag_field(
'ocean_model',
'ri_grad_shear', diag%axesTi, time, &
294 'Gradient Richarson number used by MOM_CVMix_shear module',
'nondim')
295 if (cs%id_ri_grad > 0)
then
296 allocate( cs%ri_grad( szi_(g), szj_(g), szk_(g)+1 )) ; cs%ri_grad(:,:,:) = 1.e8
299 cs%id_ri_grad_smooth = register_diag_field(
'ocean_model',
'ri_grad_shear_smooth', &
301 'Smoothed gradient Richarson number used by MOM_CVMix_shear module',
'nondim')
302 if (cs%id_ri_grad_smooth > 0)
then
303 allocate( cs%ri_grad_smooth( szi_(g), szj_(g), szk_(g)+1 )) ; cs%ri_grad_smooth(:,:,:) = 1.e8
306 cs%id_kd = register_diag_field(
'ocean_model',
'kd_shear_CVMix', diag%axesTi, time, &
307 'Vertical diffusivity added by MOM_CVMix_shear module',
'm2/s', conversion=us%Z2_T_to_m2_s)
308 cs%id_kv = register_diag_field(
'ocean_model',
'kv_shear_CVMix', diag%axesTi, time, &
309 'Vertical viscosity added by MOM_CVMix_shear module',
'm2/s', conversion=us%Z2_T_to_m2_s)
311 end function cvmix_shear_init
316 logical function cvmix_shear_is_used(param_file)
319 logical :: lmd94, pp81
320 call get_param(param_file, mdl,
"USE_LMD94", lmd94, &
321 default=.false., do_not_log = .true.)
322 call get_param(param_file, mdl,
"Use_PP81", pp81, &
323 default=.false., do_not_log = .true.)
324 cvmix_shear_is_used = (lmd94 .or. pp81)
325 end function cvmix_shear_is_used
328 subroutine cvmix_shear_end(CS)
332 if (.not.
associated(cs))
return
334 if (cs%id_N2 > 0)
deallocate(cs%N2)
335 if (cs%id_S2 > 0)
deallocate(cs%S2)
336 if (cs%id_ri_grad > 0)
deallocate(cs%ri_grad)
339 end subroutine cvmix_shear_end