20 use mom_cpu_clock,
only : cpu_clock_id, cpu_clock_begin, cpu_clock_end
23 use cvmix_kpp,
only : cvmix_init_kpp, cvmix_put_kpp, cvmix_get_kpp_real
24 use cvmix_kpp,
only : cvmix_coeffs_kpp
25 use cvmix_kpp,
only : cvmix_kpp_compute_obl_depth
26 use cvmix_kpp,
only : cvmix_kpp_compute_turbulent_scales
27 use cvmix_kpp,
only : cvmix_kpp_compute_bulk_richardson
28 use cvmix_kpp,
only : cvmix_kpp_compute_unresolved_shear
29 use cvmix_kpp,
only : cvmix_kpp_params_type
30 use cvmix_kpp,
only : cvmix_kpp_compute_kobl_depth
32 implicit none ;
private 34 #include "MOM_memory.h" 37 public :: kpp_compute_bld
38 public :: kpp_calculate
40 public :: kpp_nonlocaltransport_temp
41 public :: kpp_nonlocaltransport_saln
45 integer,
private,
parameter :: nlt_shape_cvmix = 0
46 integer,
private,
parameter :: nlt_shape_linear = 1
47 integer,
private,
parameter :: nlt_shape_parabolic = 2
48 integer,
private,
parameter :: nlt_shape_cubic = 3
49 integer,
private,
parameter :: nlt_shape_cubic_lmd = 4
52 integer,
private,
parameter :: sw_method_all_sw = 0
53 integer,
private,
parameter :: sw_method_mxl_sw = 1
54 integer,
private,
parameter :: sw_method_lv1_sw = 2
55 integer,
private,
parameter :: lt_k_constant = 1, & !< Constant enhance K through column
57 lt_k_mode_constant = 1, &
62 lt_vt2_mode_constant = 1, &
63 lt_vt2_mode_vr12 = 2, &
65 lt_vt2_mode_rw16 = 3, &
79 logical :: enhance_diffusion
80 character(len=10) :: interptype
81 character(len=10) :: interptype2
82 logical :: computeekman
83 logical :: computemoninobukhov
84 logical :: passivemode
88 real :: surf_layer_ext
90 logical :: fixedobldepth
91 real :: fixedobldepth_value
93 character(len=30) :: matchtechnique
95 logical :: applynonlocaltrans
97 logical :: deepen_only
98 logical :: kppzerodiffusivity
100 logical :: kppisadditive
103 real :: min_thickness
106 logical :: correctsurflayeravg
107 real :: surflayerdepth
110 logical :: lt_k_enhancement
111 integer :: lt_k_shape
112 integer :: lt_k_method
113 real :: kpp_k_enh_fac
114 logical :: lt_vt2_enhancement
115 integer :: lt_vt2_method
116 real :: kpp_vt2_enh_fac
117 logical :: stokes_mixing
121 type(cvmix_kpp_params_type),
pointer :: kpp_params => null()
125 integer :: id_obldepth = -1, id_bulkri = -1
126 integer :: id_n = -1, id_n2 = -1
127 integer :: id_ws = -1, id_vt2 = -1
128 integer :: id_bulkuz2 = -1, id_bulkdrho = -1
129 integer :: id_ustar = -1, id_buoyflux = -1
130 integer :: id_qminussw = -1, id_nets = -1
131 integer :: id_sigma = -1, id_kv_kpp = -1
132 integer :: id_kt_kpp = -1, id_ks_kpp = -1
133 integer :: id_tsurf = -1, id_ssurf = -1
134 integer :: id_usurf = -1, id_vsurf = -1
135 integer :: id_kd_in = -1
136 integer :: id_nltt = -1
137 integer :: id_nlts = -1
138 integer :: id_nlt_dsdt = -1
139 integer :: id_nlt_dtdt = -1
140 integer :: id_nlt_temp_budget = -1
141 integer :: id_nlt_saln_budget = -1
142 integer :: id_enhk = -1, id_enhvt2 = -1
143 integer :: id_enhw = -1
144 integer :: id_la_sl = -1
145 integer :: id_obldepth_original = -1
149 real,
allocatable,
dimension(:,:) :: obldepth
150 real,
allocatable,
dimension(:,:) :: obldepth_original
151 real,
allocatable,
dimension(:,:) :: kobl
152 real,
allocatable,
dimension(:,:) :: obldepthprev
153 real,
allocatable,
dimension(:,:) :: la_sl
154 real,
allocatable,
dimension(:,:,:) :: drho
155 real,
allocatable,
dimension(:,:,:) :: uz2
156 real,
allocatable,
dimension(:,:,:) :: bulkri
157 real,
allocatable,
dimension(:,:,:) :: sigma
158 real,
allocatable,
dimension(:,:,:) :: ws
159 real,
allocatable,
dimension(:,:,:) :: n
160 real,
allocatable,
dimension(:,:,:) :: n2
161 real,
allocatable,
dimension(:,:,:) :: vt2
162 real,
allocatable,
dimension(:,:,:) :: kt_kpp
163 real,
allocatable,
dimension(:,:,:) :: ks_kpp
164 real,
allocatable,
dimension(:,:,:) :: kv_kpp
165 real,
allocatable,
dimension(:,:) :: tsurf
166 real,
allocatable,
dimension(:,:) :: ssurf
167 real,
allocatable,
dimension(:,:) :: usurf
168 real,
allocatable,
dimension(:,:) :: vsurf
169 real,
allocatable,
dimension(:,:,:) :: enhk
170 real,
allocatable,
dimension(:,:,:) :: enhvt2
175 integer :: id_clock_kpp_calc, id_clock_kpp_compute_bld, id_clock_kpp_smoothing
178 #define __DO_SAFETY_CHECKS__ 184 logical function kpp_init(paramFile, G, GV, US, diag, Time, CS, passive, Waves)
191 type(
diag_ctrl),
target,
intent(in) :: diag
192 type(time_type),
intent(in) :: time
193 type(
kpp_cs),
pointer :: cs
194 logical,
optional,
intent(out) :: passive
198 # include "version_variable.h" 199 character(len=40) :: mdl =
'MOM_CVMix_KPP' 200 character(len=20) :: string
201 logical :: cs_is_one=.false.
202 logical :: lnodgat1=.false.
204 if (
associated(cs))
call mom_error(fatal,
'MOM_CVMix_KPP, KPP_init: '// &
205 'Control structure has already been initialized')
208 call get_param(paramfile, mdl,
"USE_KPP", kpp_init, default=.false., do_not_log=.true.)
209 call log_version(paramfile, mdl, version,
'This is the MOM wrapper to CVMix:KPP\n' // &
210 'See http://cvmix.github.io/', all_default=.not.kpp_init)
211 call get_param(paramfile, mdl,
"USE_KPP", kpp_init, &
212 "If true, turns on the [CVMix] KPP scheme of Large et al., 1994, "// &
213 "to calculate diffusivities and non-local transport in the OBL.", &
216 if (.not. kpp_init)
return 219 call openparameterblock(paramfile,
'KPP')
220 call get_param(paramfile, mdl,
'PASSIVE', cs%passiveMode, &
221 'If True, puts KPP into a passive-diagnostic mode.', &
225 if (
present(passive)) passive=cs%passiveMode
227 call get_param(paramfile, mdl,
'APPLY_NONLOCAL_TRANSPORT', cs%applyNonLocalTrans, &
228 'If True, applies the non-local transport to heat and scalars. '// &
229 'If False, calculates the non-local transport and tendencies but '//&
230 'purely for diagnostic purposes.', &
231 default=.not. cs%passiveMode)
232 call get_param(paramfile, mdl,
'N_SMOOTH', cs%n_smooth, &
233 'The number of times the 1-1-4-1-1 Laplacian filter is applied on '// &
236 if (cs%n_smooth > g%domain%nihalo)
then 237 call mom_error(fatal,
'KPP smoothing number (N_SMOOTH) cannot be greater than NIHALO.')
238 elseif (cs%n_smooth > g%domain%njhalo)
then 239 call mom_error(fatal,
'KPP smoothing number (N_SMOOTH) cannot be greater than NJHALO.')
241 if (cs%n_smooth > 0)
then 242 call get_param(paramfile, mdl,
'DEEPEN_ONLY_VIA_SMOOTHING', cs%deepen_only, &
243 'If true, apply OBLdepth smoothing at a cell only if the OBLdepth '// &
244 'gets deeper via smoothing.', &
246 id_clock_kpp_smoothing = cpu_clock_id(
'(Ocean KPP BLD smoothing)', grain=clock_routine)
248 call get_param(paramfile, mdl,
'RI_CRIT', cs%Ri_crit, &
249 'Critical bulk Richardson number used to define depth of the '// &
250 'surface Ocean Boundary Layer (OBL).', &
251 units=
'nondim', default=0.3)
252 call get_param(paramfile, mdl,
'VON_KARMAN', cs%vonKarman, &
253 'von Karman constant.', &
254 units=
'nondim', default=0.40)
255 call get_param(paramfile, mdl,
'ENHANCE_DIFFUSION', cs%enhance_diffusion, &
256 'If True, adds enhanced diffusion at the based of the boundary layer.', &
258 call get_param(paramfile, mdl,
'INTERP_TYPE', cs%interpType, &
259 'Type of interpolation to determine the OBL depth.\n'// &
260 'Allowed types are: linear, quadratic, cubic.', &
262 call get_param(paramfile, mdl,
'INTERP_TYPE2', cs%interpType2, &
263 'Type of interpolation to compute diff and visc at OBL_depth.\n'// &
264 'Allowed types are: linear, quadratic, cubic or LMD94.', &
266 call get_param(paramfile, mdl,
'COMPUTE_EKMAN', cs%computeEkman, &
267 'If True, limit OBL depth to be no deeper than Ekman depth.', &
269 call get_param(paramfile, mdl,
'COMPUTE_MONIN_OBUKHOV', cs%computeMoninObukhov, &
270 'If True, limit the OBL depth to be no deeper than '// &
271 'Monin-Obukhov depth.', &
273 call get_param(paramfile, mdl,
'CS', cs%cs, &
274 'Parameter for computing velocity scale function.', &
275 units=
'nondim', default=98.96)
276 call get_param(paramfile, mdl,
'CS2', cs%cs2, &
277 'Parameter for computing non-local term.', &
278 units=
'nondim', default=6.32739901508)
279 call get_param(paramfile, mdl,
'DEEP_OBL_OFFSET', cs%deepOBLoffset, &
280 'If non-zero, the distance above the bottom to which the OBL is clipped '// &
281 'if it would otherwise reach the bottom. The smaller of this and 0.1D is used.', &
282 units=
'm',default=0.)
283 call get_param(paramfile, mdl,
'FIXED_OBLDEPTH', cs%fixedOBLdepth, &
284 'If True, fix the OBL depth to FIXED_OBLDEPTH_VALUE '// &
285 'rather than using the OBL depth from CVMix. '// &
286 'This option is just for testing purposes.', &
288 call get_param(paramfile, mdl,
'FIXED_OBLDEPTH_VALUE', cs%fixedOBLdepth_value, &
289 'Value for the fixed OBL depth when fixedOBLdepth==True. '// &
290 'This parameter is for just for testing purposes. '// &
291 'It will over-ride the OBLdepth computed from CVMix.', &
292 units=
'm',default=30.0)
293 call get_param(paramfile, mdl,
'SURF_LAYER_EXTENT', cs%surf_layer_ext, &
294 'Fraction of OBL depth considered in the surface layer.', &
295 units=
'nondim',default=0.10)
296 call get_param(paramfile, mdl,
'MINIMUM_OBL_DEPTH', cs%minOBLdepth, &
297 'If non-zero, a minimum depth to use for KPP OBL depth. Independent of '// &
298 'this parameter, the OBL depth is always at least as deep as the first layer.', &
299 units=
'm',default=0.)
300 call get_param(paramfile, mdl,
'MINIMUM_VT2', cs%minVtsqr, &
301 'Min of the unresolved velocity Vt2 used in Rib CVMix calculation.\n'// &
302 'Scaling: MINIMUM_VT2 = const1*d*N*ws, with d=1m, N=1e-5/s, ws=1e-6 m/s.', &
303 units=
'm2/s2',default=1e-10)
306 call get_param(paramfile, mdl,
'CORRECT_SURFACE_LAYER_AVERAGE', cs%correctSurfLayerAvg, &
307 'If true, applies a correction step to the averaging of surface layer '// &
308 'properties. This option is obsolete.', default=.false.)
309 if (cs%correctSurfLayerAvg) &
310 call mom_error(fatal,
'Correct surface layer average disabled in code. To recover \n'// &
311 ' feature will require code intervention.')
312 call get_param(paramfile, mdl,
'FIRST_GUESS_SURFACE_LAYER_DEPTH', cs%surfLayerDepth, &
313 'The first guess at the depth of the surface layer used for averaging '// &
314 'the surface layer properties. If =0, the top model level properties '// &
315 'will be used for the surface layer. If CORRECT_SURFACE_LAYER_AVERAGE=True, a '// &
316 'subsequent correction is applied. This parameter is obsolete', units=
'm', default=0.)
319 call get_param(paramfile, mdl,
'NLT_SHAPE', string, &
320 'MOM6 method to set nonlocal transport profile. '// &
321 'Over-rides the result from CVMix. Allowed values are: \n'// &
322 '\t CVMix - Uses the profiles from CVMix specified by MATCH_TECHNIQUE\n'//&
323 '\t LINEAR - A linear profile, 1-sigma\n'// &
324 '\t PARABOLIC - A parablic profile, (1-sigma)^2\n'// &
325 '\t CUBIC - A cubic profile, (1-sigma)^2(1+2*sigma)\n'// &
326 '\t CUBIC_LMD - The original KPP profile', &
328 select case ( trim(string) )
329 case (
"CVMix") ; cs%NLT_shape = nlt_shape_cvmix
330 case (
"LINEAR") ; cs%NLT_shape = nlt_shape_linear
331 case (
"PARABOLIC") ; cs%NLT_shape = nlt_shape_parabolic
332 case (
"CUBIC") ; cs%NLT_shape = nlt_shape_cubic
333 case (
"CUBIC_LMD") ; cs%NLT_shape = nlt_shape_cubic_lmd
334 case default ;
call mom_error(fatal,
"KPP_init: "// &
335 "Unrecognized NLT_SHAPE option"//trim(string))
337 call get_param(paramfile, mdl,
'MATCH_TECHNIQUE', cs%MatchTechnique, &
338 'CVMix method to set profile function for diffusivity and NLT, '// &
339 'as well as matching across OBL base. Allowed values are: \n'// &
340 '\t SimpleShapes = sigma*(1-sigma)^2 for both diffusivity and NLT\n'// &
341 '\t MatchGradient = sigma*(1-sigma)^2 for NLT; diffusivity profile from matching\n'//&
342 '\t MatchBoth = match gradient for both diffusivity and NLT\n'// &
343 '\t ParabolicNonLocal = sigma*(1-sigma)^2 for diffusivity; (1-sigma)^2 for NLT', &
344 default=
'SimpleShapes')
345 if (cs%MatchTechnique ==
'ParabolicNonLocal')
then 350 if (cs%MatchTechnique ==
'ParabolicNonLocal' .or. cs%MatchTechnique ==
'SimpleShapes')
then 356 if (cs%MatchTechnique ==
'MatchBoth' .and. (cs%interpType2 ==
'cubic' .or. &
357 cs%interpType2 ==
'quadratic'))
then 358 call mom_error(fatal,
"If MATCH_TECHNIQUE=MatchBoth, INTERP_TYPE2 must be set to \n"//&
359 "linear or LMD94 (recommended) to avoid negative viscosity and diffusivity.\n"//&
360 "Please select one of these valid options." )
363 call get_param(paramfile, mdl,
'KPP_ZERO_DIFFUSIVITY', cs%KPPzeroDiffusivity, &
364 'If True, zeroes the KPP diffusivity and viscosity; for testing purpose.',&
366 call get_param(paramfile, mdl,
'KPP_IS_ADDITIVE', cs%KPPisAdditive, &
367 'If true, adds KPP diffusivity to diffusivity from other schemes.\n'//&
368 'If false, KPP is the only diffusivity wherever KPP is non-zero.', &
370 call get_param(paramfile, mdl,
'KPP_SHORTWAVE_METHOD',string, &
371 'Determines contribution of shortwave radiation to KPP surface '// &
372 'buoyancy flux. Options include:\n'// &
373 ' ALL_SW: use total shortwave radiation\n'// &
374 ' MXL_SW: use shortwave radiation absorbed by mixing layer\n'// &
375 ' LV1_SW: use shortwave radiation absorbed by top model layer', &
377 select case ( trim(string) )
378 case (
"ALL_SW") ; cs%SW_METHOD = sw_method_all_sw
379 case (
"MXL_SW") ; cs%SW_METHOD = sw_method_mxl_sw
380 case (
"LV1_SW") ; cs%SW_METHOD = sw_method_lv1_sw
381 case default ;
call mom_error(fatal,
"KPP_init: "// &
382 "Unrecognized KPP_SHORTWAVE_METHOD option"//trim(string))
384 call get_param(paramfile, mdl,
'CVMix_ZERO_H_WORK_AROUND', cs%min_thickness, &
385 'A minimum thickness used to avoid division by small numbers in the vicinity '// &
386 'of vanished layers. This is independent of MIN_THICKNESS used in other parts of MOM.', &
387 units=
'm', default=0.)
391 call get_param(paramfile, mdl,
"USE_KPP_LT_K", cs%LT_K_Enhancement, &
392 'Flag for Langmuir turbulence enhancement of turbulent'//&
393 'mixing coefficient.', units=
"", default=.false.)
394 call get_param(paramfile, mdl,
"STOKES_MIXING", cs%STOKES_MIXING, &
395 'Flag for Langmuir turbulence enhancement of turbulent'//&
396 'mixing coefficient.', units=
"", default=.false.)
397 if (cs%LT_K_Enhancement)
then 398 call get_param(paramfile, mdl,
'KPP_LT_K_SHAPE', string, &
399 'Vertical dependence of LT enhancement of mixing. '// &
400 'Valid options are: \n'// &
401 '\t CONSTANT = Constant value for full OBL\n'// &
402 '\t SCALED = Varies based on normalized shape function.', &
404 select case ( trim(string))
405 case (
"CONSTANT") ; cs%LT_K_SHAPE = lt_k_constant
406 case (
"SCALED") ; cs%LT_K_SHAPE = lt_k_scaled
407 case default ;
call mom_error(fatal,
"KPP_init: "//&
408 "Unrecognized KPP_LT_K_SHAPE option: "//trim(string))
410 call get_param(paramfile, mdl,
"KPP_LT_K_METHOD", string , &
411 'Method to enhance mixing coefficient in KPP. '// &
412 'Valid options are: \n'// &
413 '\t CONSTANT = Constant value (KPP_K_ENH_FAC) \n'// &
414 '\t VR12 = Function of Langmuir number based on VR12\n'// &
415 '\t RW16 = Function of Langmuir number based on RW16', &
417 select case ( trim(string))
418 case (
"CONSTANT") ; cs%LT_K_METHOD = lt_k_mode_constant
419 case (
"VR12") ; cs%LT_K_METHOD = lt_k_mode_vr12
420 case (
"RW16") ; cs%LT_K_METHOD = lt_k_mode_rw16
421 case default ;
call mom_error(fatal,
"KPP_init: "//&
422 "Unrecognized KPP_LT_K_METHOD option: "//trim(string))
424 if (cs%LT_K_METHOD==lt_k_mode_constant)
then 425 call get_param(paramfile, mdl,
"KPP_K_ENH_FAC",cs%KPP_K_ENH_FAC , &
426 'Constant value to enhance mixing coefficient in KPP.', &
431 call get_param(paramfile, mdl,
"USE_KPP_LT_VT2", cs%LT_Vt2_Enhancement, &
432 'Flag for Langmuir turbulence enhancement of Vt2'//&
433 'in Bulk Richardson Number.', units=
"", default=.false.)
434 if (cs%LT_Vt2_Enhancement)
then 435 call get_param(paramfile, mdl,
"KPP_LT_VT2_METHOD",string , &
436 'Method to enhance Vt2 in KPP. '// &
437 'Valid options are: \n'// &
438 '\t CONSTANT = Constant value (KPP_VT2_ENH_FAC) \n'// &
439 '\t VR12 = Function of Langmuir number based on VR12\n'// &
440 '\t RW16 = Function of Langmuir number based on RW16\n'// &
441 '\t LF17 = Function of Langmuir number based on LF17', &
443 select case ( trim(string))
444 case (
"CONSTANT") ; cs%LT_VT2_METHOD = lt_vt2_mode_constant
445 case (
"VR12") ; cs%LT_VT2_METHOD = lt_vt2_mode_vr12
446 case (
"RW16") ; cs%LT_VT2_METHOD = lt_vt2_mode_rw16
447 case (
"LF17") ; cs%LT_VT2_METHOD = lt_vt2_mode_lf17
448 case default ;
call mom_error(fatal,
"KPP_init: "//&
449 "Unrecognized KPP_LT_VT2_METHOD option: "//trim(string))
451 if (cs%LT_VT2_METHOD==lt_vt2_mode_constant)
then 452 call get_param(paramfile, mdl,
"KPP_VT2_ENH_FAC",cs%KPP_VT2_ENH_FAC , &
453 'Constant value to enhance VT2 in KPP.', &
458 call closeparameterblock(paramfile)
459 call get_param(paramfile, mdl,
'DEBUG', cs%debug, default=.false., do_not_log=.true.)
461 call cvmix_init_kpp( ri_crit=cs%Ri_crit, &
462 minobldepth=cs%minOBLdepth, &
463 minvtsqr=cs%minVtsqr, &
464 vonkarman=cs%vonKarman, &
465 surf_layer_ext=cs%surf_layer_ext, &
466 interp_type=cs%interpType, &
467 interp_type2=cs%interpType2, &
468 lekman=cs%computeEkman, &
469 lmonob=cs%computeMoninObukhov, &
470 matchtechnique=cs%MatchTechnique, &
471 lenhanced_diff=cs%enhance_diffusion,&
472 lnonzero_surf_nonlocal=cs_is_one ,&
474 cvmix_kpp_params_user=cs%KPP_params )
478 cs%id_OBLdepth = register_diag_field(
'ocean_model',
'KPP_OBLdepth', diag%axesT1, time, &
479 'Thickness of the surface Ocean Boundary Layer calculated by [CVMix] KPP',
'meter', &
480 cmor_field_name=
'oml', cmor_long_name=
'ocean_mixed_layer_thickness_defined_by_mixing_scheme', &
481 cmor_units=
'm', cmor_standard_name=
'Ocean Mixed Layer Thickness Defined by Mixing Scheme')
485 if (cs%n_smooth > 0)
then 486 cs%id_OBLdepth_original = register_diag_field(
'ocean_model',
'KPP_OBLdepth_original', diag%axesT1, time, &
487 'Thickness of the surface Ocean Boundary Layer without smoothing calculated by [CVMix] KPP',
'meter', &
488 cmor_field_name=
'oml', cmor_long_name=
'ocean_mixed_layer_thickness_defined_by_mixing_scheme', &
489 cmor_units=
'm', cmor_standard_name=
'Ocean Mixed Layer Thickness Defined by Mixing Scheme')
491 cs%id_BulkDrho = register_diag_field(
'ocean_model',
'KPP_BulkDrho', diag%axesTL, time, &
492 'Bulk difference in density used in Bulk Richardson number, as used by [CVMix] KPP', &
493 'kg/m3', conversion=us%R_to_kg_m3)
494 cs%id_BulkUz2 = register_diag_field(
'ocean_model',
'KPP_BulkUz2', diag%axesTL, time, &
495 'Square of bulk difference in resolved velocity used in Bulk Richardson number via [CVMix] KPP',
'm2/s2')
496 cs%id_BulkRi = register_diag_field(
'ocean_model',
'KPP_BulkRi', diag%axesTL, time, &
497 'Bulk Richardson number used to find the OBL depth used by [CVMix] KPP',
'nondim')
498 cs%id_Sigma = register_diag_field(
'ocean_model',
'KPP_sigma', diag%axesTi, time, &
499 'Sigma coordinate used by [CVMix] KPP',
'nondim')
500 cs%id_Ws = register_diag_field(
'ocean_model',
'KPP_Ws', diag%axesTL, time, &
501 'Turbulent vertical velocity scale for scalars used by [CVMix] KPP',
'm/s')
502 cs%id_N = register_diag_field(
'ocean_model',
'KPP_N', diag%axesTi, time, &
503 '(Adjusted) Brunt-Vaisala frequency used by [CVMix] KPP',
'1/s')
504 cs%id_N2 = register_diag_field(
'ocean_model',
'KPP_N2', diag%axesTi, time, &
505 'Square of Brunt-Vaisala frequency used by [CVMix] KPP',
'1/s2')
506 cs%id_Vt2 = register_diag_field(
'ocean_model',
'KPP_Vt2', diag%axesTL, time, &
507 'Unresolved shear turbulence used by [CVMix] KPP',
'm2/s2')
508 cs%id_uStar = register_diag_field(
'ocean_model',
'KPP_uStar', diag%axesT1, time, &
509 'Friction velocity, u*, as used by [CVMix] KPP',
'm/s', conversion=us%Z_to_m*us%s_to_T)
510 cs%id_buoyFlux = register_diag_field(
'ocean_model',
'KPP_buoyFlux', diag%axesTi, time, &
511 'Surface (and penetrating) buoyancy flux, as used by [CVMix] KPP',
'm2/s3', conversion=us%L_to_m**2*us%s_to_T**3)
512 cs%id_QminusSW = register_diag_field(
'ocean_model',
'KPP_QminusSW', diag%axesT1, time, &
513 'Net temperature flux ignoring short-wave, as used by [CVMix] KPP',
'K m/s')
514 cs%id_netS = register_diag_field(
'ocean_model',
'KPP_netSalt', diag%axesT1, time, &
515 'Effective net surface salt flux, as used by [CVMix] KPP',
'ppt m/s')
516 cs%id_Kt_KPP = register_diag_field(
'ocean_model',
'KPP_Kheat', diag%axesTi, time, &
517 'Heat diffusivity due to KPP, as calculated by [CVMix] KPP',
'm2/s')
518 cs%id_Kd_in = register_diag_field(
'ocean_model',
'KPP_Kd_in', diag%axesTi, time, &
519 'Diffusivity passed to KPP',
'm2/s', conversion=us%Z2_T_to_m2_s)
520 cs%id_Ks_KPP = register_diag_field(
'ocean_model',
'KPP_Ksalt', diag%axesTi, time, &
521 'Salt diffusivity due to KPP, as calculated by [CVMix] KPP',
'm2/s')
522 cs%id_Kv_KPP = register_diag_field(
'ocean_model',
'KPP_Kv', diag%axesTi, time, &
523 'Vertical viscosity due to KPP, as calculated by [CVMix] KPP',
'm2/s')
524 cs%id_NLTt = register_diag_field(
'ocean_model',
'KPP_NLtransport_heat', diag%axesTi, time, &
525 'Non-local transport (Cs*G(sigma)) for heat, as calculated by [CVMix] KPP',
'nondim')
526 cs%id_NLTs = register_diag_field(
'ocean_model',
'KPP_NLtransport_salt', diag%axesTi, time, &
527 'Non-local tranpsort (Cs*G(sigma)) for scalars, as calculated by [CVMix] KPP',
'nondim')
528 cs%id_NLT_dTdt = register_diag_field(
'ocean_model',
'KPP_NLT_dTdt', diag%axesTL, time, &
529 'Temperature tendency due to non-local transport of heat, as calculated by [CVMix] KPP',
'K/s')
530 cs%id_NLT_dSdt = register_diag_field(
'ocean_model',
'KPP_NLT_dSdt', diag%axesTL, time, &
531 'Salinity tendency due to non-local transport of salt, as calculated by [CVMix] KPP',
'ppt/s')
532 cs%id_NLT_temp_budget = register_diag_field(
'ocean_model',
'KPP_NLT_temp_budget', diag%axesTL, time, &
533 'Heat content change due to non-local transport, as calculated by [CVMix] KPP',
'W/m^2')
534 cs%id_NLT_saln_budget = register_diag_field(
'ocean_model',
'KPP_NLT_saln_budget', diag%axesTL, time, &
535 'Salt content change due to non-local transport, as calculated by [CVMix] KPP',
'kg/(sec*m^2)')
536 cs%id_Tsurf = register_diag_field(
'ocean_model',
'KPP_Tsurf', diag%axesT1, time, &
537 'Temperature of surface layer (10% of OBL depth) as passed to [CVMix] KPP',
'C')
538 cs%id_Ssurf = register_diag_field(
'ocean_model',
'KPP_Ssurf', diag%axesT1, time, &
539 'Salinity of surface layer (10% of OBL depth) as passed to [CVMix] KPP',
'ppt')
540 cs%id_Usurf = register_diag_field(
'ocean_model',
'KPP_Usurf', diag%axesCu1, time, &
541 'i-component flow of surface layer (10% of OBL depth) as passed to [CVMix] KPP',
'm/s')
542 cs%id_Vsurf = register_diag_field(
'ocean_model',
'KPP_Vsurf', diag%axesCv1, time, &
543 'j-component flow of surface layer (10% of OBL depth) as passed to [CVMix] KPP',
'm/s')
544 cs%id_EnhK = register_diag_field(
'ocean_model',
'EnhK', diag%axesTI, time, &
545 'Langmuir number enhancement to K as used by [CVMix] KPP',
'nondim')
546 cs%id_EnhVt2 = register_diag_field(
'ocean_model',
'EnhVt2', diag%axesTL, time, &
547 'Langmuir number enhancement to Vt2 as used by [CVMix] KPP',
'nondim')
548 cs%id_La_SL = register_diag_field(
'ocean_model',
'KPP_La_SL', diag%axesT1, time, &
549 'Surface-layer Langmuir number computed in [CVMix] KPP',
'nondim')
551 allocate( cs%N( szi_(g), szj_(g), szk_(g)+1 ) )
553 allocate( cs%OBLdepth( szi_(g), szj_(g) ) )
554 cs%OBLdepth(:,:) = 0.
555 allocate( cs%kOBL( szi_(g), szj_(g) ) )
557 allocate( cs%La_SL( szi_(g), szj_(g) ) )
559 allocate( cs%Vt2( szi_(g), szj_(g), szk_(g) ) )
561 if (cs%id_OBLdepth_original > 0)
allocate( cs%OBLdepth_original( szi_(g), szj_(g) ) )
563 allocate( cs%OBLdepthprev( szi_(g), szj_(g) ) ) ; cs%OBLdepthprev(:,:) = 0.0
564 if (cs%id_BulkDrho > 0)
allocate( cs%dRho( szi_(g), szj_(g), szk_(g) ) )
565 if (cs%id_BulkDrho > 0) cs%dRho(:,:,:) = 0.
566 if (cs%id_BulkUz2 > 0)
allocate( cs%Uz2( szi_(g), szj_(g), szk_(g) ) )
567 if (cs%id_BulkUz2 > 0) cs%Uz2(:,:,:) = 0.
568 if (cs%id_BulkRi > 0)
allocate( cs%BulkRi( szi_(g), szj_(g), szk_(g) ) )
569 if (cs%id_BulkRi > 0) cs%BulkRi(:,:,:) = 0.
570 if (cs%id_Sigma > 0)
allocate( cs%sigma( szi_(g), szj_(g), szk_(g)+1 ) )
571 if (cs%id_Sigma > 0) cs%sigma(:,:,:) = 0.
572 if (cs%id_Ws > 0)
allocate( cs%Ws( szi_(g), szj_(g), szk_(g) ) )
573 if (cs%id_Ws > 0) cs%Ws(:,:,:) = 0.
574 if (cs%id_N2 > 0)
allocate( cs%N2( szi_(g), szj_(g), szk_(g)+1 ) )
575 if (cs%id_N2 > 0) cs%N2(:,:,:) = 0.
576 if (cs%id_Kt_KPP > 0)
allocate( cs%Kt_KPP( szi_(g), szj_(g), szk_(g)+1 ) )
577 if (cs%id_Kt_KPP > 0) cs%Kt_KPP(:,:,:) = 0.
578 if (cs%id_Ks_KPP > 0)
allocate( cs%Ks_KPP( szi_(g), szj_(g), szk_(g)+1 ) )
579 if (cs%id_Ks_KPP > 0) cs%Ks_KPP(:,:,:) = 0.
580 if (cs%id_Kv_KPP > 0)
allocate( cs%Kv_KPP( szi_(g), szj_(g), szk_(g)+1 ) )
581 if (cs%id_Kv_KPP > 0) cs%Kv_KPP(:,:,:) = 0.
582 if (cs%id_Tsurf > 0)
allocate( cs%Tsurf( szi_(g), szj_(g)) )
583 if (cs%id_Tsurf > 0) cs%Tsurf(:,:) = 0.
584 if (cs%id_Ssurf > 0)
allocate( cs%Ssurf( szi_(g), szj_(g)) )
585 if (cs%id_Ssurf > 0) cs%Ssurf(:,:) = 0.
586 if (cs%id_Usurf > 0)
allocate( cs%Usurf( szib_(g), szj_(g)) )
587 if (cs%id_Usurf > 0) cs%Usurf(:,:) = 0.
588 if (cs%id_Vsurf > 0)
allocate( cs%Vsurf( szi_(g), szjb_(g)) )
589 if (cs%id_Vsurf > 0) cs%Vsurf(:,:) = 0.
590 if (cs%id_EnhVt2 > 0)
allocate( cs%EnhVt2( szi_(g), szj_(g), szk_(g)) )
591 if (cs%id_EnhVt2 > 0) cs%EnhVt2(:,:,:) = 0.
592 if (cs%id_EnhK > 0)
allocate( cs%EnhK( szi_(g), szj_(g), szk_(g)+1 ) )
593 if (cs%id_EnhK > 0) cs%EnhK(:,:,:) = 0.
595 id_clock_kpp_calc = cpu_clock_id(
'Ocean KPP calculate)', grain=clock_module)
596 id_clock_kpp_compute_bld = cpu_clock_id(
'(Ocean KPP comp BLD)', grain=clock_routine)
598 end function kpp_init
601 subroutine kpp_calculate(CS, G, GV, US, h, uStar, &
602 buoyFlux, Kt, Ks, Kv, nonLocalTransHeat,&
603 nonLocalTransScalar, waves)
606 type(
kpp_cs),
pointer :: cs
611 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: h
612 real,
dimension(SZI_(G),SZJ_(G)),
intent(in) :: ustar
613 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1),
intent(in) :: buoyflux
614 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1),
intent(inout) :: kt
617 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1),
intent(inout) :: ks
620 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1),
intent(inout) :: kv
623 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1),
intent(inout) :: nonlocaltransheat
624 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1),
intent(inout) :: nonlocaltransscalar
628 real,
dimension( G%ke ) :: cellheight
629 real,
dimension( G%ke+1 ) :: ifaceheight
630 real,
dimension( G%ke+1, 2) :: kdiffusivity
631 real,
dimension( G%ke+1 ) :: kviscosity
632 real,
dimension( G%ke+1, 2) :: nonlocaltrans
634 real :: surffricvel, surfbuoyflux
635 real :: sigma, sigmaratio
645 call hchksum(h,
"KPP in: h",g%HI,haloshift=0, scale=gv%H_to_m)
646 call hchksum(ustar,
"KPP in: uStar",g%HI,haloshift=0, scale=us%Z_to_m*us%s_to_T)
647 call hchksum(buoyflux,
"KPP in: buoyFlux",g%HI,haloshift=0)
648 call hchksum(kt,
"KPP in: Kt",g%HI,haloshift=0, scale=us%Z2_T_to_m2_s)
649 call hchksum(ks,
"KPP in: Ks",g%HI,haloshift=0, scale=us%Z2_T_to_m2_s)
652 nonlocaltrans(:,:) = 0.0
654 if (cs%id_Kd_in > 0)
call post_data(cs%id_Kd_in, kt, cs%diag)
656 call cpu_clock_begin(id_clock_kpp_calc)
657 buoy_scale = us%L_to_m**2*us%s_to_T**3
670 if (g%mask2dT(i,j)==0.) cycle
673 surffricvel = us%Z_to_m*us%s_to_T * ustar(i,j)
680 dh = h(i,j,k) * gv%H_to_m
682 hcorr = min( dh - cs%min_thickness, 0. )
683 dh = max( dh, cs%min_thickness )
684 cellheight(k) = ifaceheight(k) - 0.5 * dh
685 ifaceheight(k+1) = ifaceheight(k) - dh
689 surfbuoyflux = buoy_scale*buoyflux(i,j,1)
698 if (cs%SW_METHOD == sw_method_all_sw)
then 699 surfbuoyflux = buoy_scale * buoyflux(i,j,1)
700 elseif (cs%SW_METHOD == sw_method_mxl_sw)
then 702 surfbuoyflux = buoy_scale * (buoyflux(i,j,1) - buoyflux(i,j,int(cs%kOBL(i,j))+1))
703 elseif (cs%SW_METHOD == sw_method_lv1_sw)
then 704 surfbuoyflux = buoy_scale * (buoyflux(i,j,1) - buoyflux(i,j,2))
708 if (.not. (cs%MatchTechnique ==
'MatchBoth'))
then 709 kdiffusivity(:,:) = 0.
712 kdiffusivity(:,1) = us%Z2_T_to_m2_s * kt(i,j,:)
713 kdiffusivity(:,2) = us%Z2_T_to_m2_s * ks(i,j,:)
714 kviscosity(:) = us%Z2_T_to_m2_s * kv(i,j,:)
717 call cvmix_coeffs_kpp(kviscosity(:), &
733 cvmix_kpp_params_user=cs%KPP_params )
737 if (kviscosity(k) < 0. .or. kdiffusivity(k,1) < 0.)
then 738 call mom_error(fatal,
"KPP_calculate, after CVMix_coeffs_kpp: "// &
739 "Negative vertical viscosity or diffusivity has been detected. " // &
740 "This is likely related to the choice of MATCH_TECHNIQUE and INTERP_TYPE2." //&
741 "You might consider using the default options for these parameters." )
745 IF (cs%LT_K_ENHANCEMENT)
then 746 if (cs%LT_K_METHOD==lt_k_mode_constant)
then 747 langenhk = cs%KPP_K_ENH_FAC
748 elseif (cs%LT_K_METHOD==lt_k_mode_vr12)
then 750 langenhk = sqrt(1.+(1.5*cs%La_SL(i,j))**(-2) + &
751 (5.4*cs%La_SL(i,j))**(-4))
752 elseif (cs%LT_K_METHOD==lt_k_mode_rw16)
then 754 langenhk = min(2.25, 1. + 1./cs%La_SL(i,j))
761 if (cs%LT_K_SHAPE== lt_k_constant)
then 762 if (cs%id_EnhK > 0) cs%EnhK(i,j,:) = langenhk
763 kdiffusivity(k,1) = kdiffusivity(k,1) * langenhk
764 kdiffusivity(k,2) = kdiffusivity(k,2) * langenhk
765 kviscosity(k) = kviscosity(k) * langenhk
766 elseif (cs%LT_K_SHAPE == lt_k_scaled)
then 767 sigma = min(1.0,-ifaceheight(k)/cs%OBLdepth(i,j))
768 sigmaratio = sigma * (1. - sigma)**2 / 0.148148037
769 if (cs%id_EnhK > 0) cs%EnhK(i,j,k) = (1.0 + (langenhk - 1.)*sigmaratio)
770 kdiffusivity(k,1) = kdiffusivity(k,1) * ( 1. + &
771 ( langenhk - 1.)*sigmaratio)
772 kdiffusivity(k,2) = kdiffusivity(k,2) * ( 1. + &
773 ( langenhk - 1.)*sigmaratio)
774 kviscosity(k) = kviscosity(k) * ( 1. + &
775 ( langenhk - 1.)*sigmaratio)
789 if (surfbuoyflux < 0.0)
then 790 if (cs%NLT_shape == nlt_shape_cubic)
then 792 sigma = min(1.0,-ifaceheight(k)/cs%OBLdepth(i,j))
793 nonlocaltrans(k,1) = (1.0 - sigma)**2 * (1.0 + 2.0*sigma)
794 nonlocaltrans(k,2) = nonlocaltrans(k,1)
796 elseif (cs%NLT_shape == nlt_shape_parabolic)
then 798 sigma = min(1.0,-ifaceheight(k)/cs%OBLdepth(i,j))
799 nonlocaltrans(k,1) = (1.0 - sigma)**2
800 nonlocaltrans(k,2) = nonlocaltrans(k,1)
802 elseif (cs%NLT_shape == nlt_shape_linear)
then 804 sigma = min(1.0,-ifaceheight(k)/cs%OBLdepth(i,j))
805 nonlocaltrans(k,1) = (1.0 - sigma)
806 nonlocaltrans(k,2) = nonlocaltrans(k,1)
808 elseif (cs%NLT_shape == nlt_shape_cubic_lmd)
then 811 sigma = min(1.0,-ifaceheight(k)/cs%OBLdepth(i,j))
812 nonlocaltrans(k,1) = cs%CS2 * sigma*(1.0 -sigma)**2
813 nonlocaltrans(k,2) = nonlocaltrans(k,1)
820 nonlocaltransheat(i,j,:) = nonlocaltrans(:,1)
821 nonlocaltransscalar(i,j,:) = nonlocaltrans(:,2)
824 if (cs%KPPzeroDiffusivity)
then 825 kdiffusivity(:,1) = 0.0
826 kdiffusivity(:,2) = 0.0
832 if (cs%id_Vt2 > 0)
then 844 cs%OBLdepthprev(i,j)=cs%OBLdepth(i,j)
845 if (cs%id_sigma > 0)
then 847 if (cs%OBLdepth(i,j)>0.) cs%sigma(i,j,:) = -ifaceheight/cs%OBLdepth(i,j)
849 if (cs%id_Kt_KPP > 0) cs%Kt_KPP(i,j,:) = kdiffusivity(:,1)
850 if (cs%id_Ks_KPP > 0) cs%Ks_KPP(i,j,:) = kdiffusivity(:,2)
851 if (cs%id_Kv_KPP > 0) cs%Kv_KPP(i,j,:) = kviscosity(:)
854 if (.not. cs%passiveMode)
then 855 if (cs%KPPisAdditive)
then 857 kt(i,j,k) = kt(i,j,k) + us%m2_s_to_Z2_T * kdiffusivity(k,1)
858 ks(i,j,k) = ks(i,j,k) + us%m2_s_to_Z2_T * kdiffusivity(k,2)
859 kv(i,j,k) = kv(i,j,k) + us%m2_s_to_Z2_T * kviscosity(k)
860 if (cs%Stokes_Mixing) waves%KvS(i,j,k) = kv(i,j,k)
864 if (kdiffusivity(k,1) /= 0.) kt(i,j,k) = us%m2_s_to_Z2_T * kdiffusivity(k,1)
865 if (kdiffusivity(k,2) /= 0.) ks(i,j,k) = us%m2_s_to_Z2_T * kdiffusivity(k,2)
866 if (kviscosity(k) /= 0.) kv(i,j,k) = us%m2_s_to_Z2_T * kviscosity(k)
867 if (cs%Stokes_Mixing) waves%KvS(i,j,k) = kv(i,j,k)
877 call cpu_clock_end(id_clock_kpp_calc)
880 call hchksum(kt,
"KPP out: Kt", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
881 call hchksum(ks,
"KPP out: Ks", g%HI, haloshift=0, scale=us%Z2_T_to_m2_s)
885 if (cs%id_OBLdepth > 0)
call post_data(cs%id_OBLdepth, cs%OBLdepth, cs%diag)
886 if (cs%id_OBLdepth_original > 0)
call post_data(cs%id_OBLdepth_original,cs%OBLdepth_original,cs%diag)
887 if (cs%id_sigma > 0)
call post_data(cs%id_sigma, cs%sigma, cs%diag)
888 if (cs%id_Ws > 0)
call post_data(cs%id_Ws, cs%Ws, cs%diag)
889 if (cs%id_Vt2 > 0)
call post_data(cs%id_Vt2, cs%Vt2, cs%diag)
890 if (cs%id_uStar > 0)
call post_data(cs%id_uStar, ustar, cs%diag)
891 if (cs%id_buoyFlux > 0)
call post_data(cs%id_buoyFlux, buoyflux, cs%diag)
892 if (cs%id_Kt_KPP > 0)
call post_data(cs%id_Kt_KPP, cs%Kt_KPP, cs%diag)
893 if (cs%id_Ks_KPP > 0)
call post_data(cs%id_Ks_KPP, cs%Ks_KPP, cs%diag)
894 if (cs%id_Kv_KPP > 0)
call post_data(cs%id_Kv_KPP, cs%Kv_KPP, cs%diag)
895 if (cs%id_NLTt > 0)
call post_data(cs%id_NLTt, nonlocaltransheat, cs%diag)
896 if (cs%id_NLTs > 0)
call post_data(cs%id_NLTs, nonlocaltransscalar,cs%diag)
899 end subroutine kpp_calculate
903 subroutine kpp_compute_bld(CS, G, GV, US, h, Temp, Salt, u, v, tv, uStar, buoyFlux, Waves)
906 type(
kpp_cs),
pointer :: cs
910 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: h
911 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: temp
912 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: salt
913 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)),
intent(in) :: u
914 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)),
intent(in) :: v
916 real,
dimension(SZI_(G),SZJ_(G)),
intent(in) :: ustar
917 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1),
intent(in) :: buoyflux
921 integer :: i, j, k, km1
922 real,
dimension( G%ke ) :: cellheight
923 real,
dimension( G%ke+1 ) :: ifaceheight
924 real,
dimension( G%ke+1 ) :: n2_1d
925 real,
dimension( G%ke ) :: ws_1d
926 real,
dimension( G%ke ) :: deltarho
927 real,
dimension( G%ke ) :: deltau2
928 real,
dimension( G%ke ) :: surfbuoyflux2
929 real,
dimension( G%ke ) :: bulkri_1d
932 real,
dimension( 3*G%ke ) :: rho_1d
933 real,
dimension( 3*G%ke ) :: pres_1d
934 real,
dimension( 3*G%ke ) :: temp_1d
935 real,
dimension( 3*G%ke ) :: salt_1d
937 real :: surffricvel, surfbuoyflux, coriolis
940 real :: rho1, rhok, uk, vk, sigma, sigmaratio
942 real :: zbottomminusoffset
947 real :: surfhtemp, surftemp
948 real :: surfhsalt, surfsalt
949 real :: surfhu, surfu
950 real :: surfhv, surfv
953 integer :: kk, ksfc, ktmp
957 real,
dimension(G%ke) :: langenhvt2
958 real,
dimension(G%ke) :: u_h, v_h
959 real :: mld_guess, la
960 real :: surfhus, surfhvs, surfus, surfvs, wavedir, currentdir
961 real :: varup, vardn, m, varlo, varavg
962 real :: h10pct, h20pct,cmnfact, usx20pct, usy20pct, enhvt2
968 call hchksum(salt,
"KPP in: S",g%HI,haloshift=0)
969 call hchksum(temp,
"KPP in: T",g%HI,haloshift=0)
970 call hchksum(u,
"KPP in: u",g%HI,haloshift=0)
971 call hchksum(v,
"KPP in: v",g%HI,haloshift=0)
974 call cpu_clock_begin(id_clock_kpp_compute_bld)
977 gorho = us%L_T_to_m_s**2*us%m_to_Z * gv%g_Earth / gv%Rho0
978 buoy_scale = us%L_to_m**2*us%s_to_T**3
995 if (g%mask2dT(i,j)==0.) cycle
998 u_h(k) = 0.5 * us%L_T_to_m_s*(u(i,j,k)+u(i-1,j,k))
999 v_h(k) = 0.5 * us%L_T_to_m_s*(v(i,j,k)+v(i,j-1,k))
1003 coriolis = 0.25*us%s_to_T*( (g%CoriolisBu(i,j) + g%CoriolisBu(i-1,j-1)) + &
1004 (g%CoriolisBu(i-1,j) + g%CoriolisBu(i,j-1)) )
1005 surffricvel = us%Z_to_m*us%s_to_T * ustar(i,j)
1013 ifaceheight(1) = 0.0
1014 pref = 0. ;
if (
associated(tv%p_surf)) pref = tv%p_surf(i,j)
1019 dh = h(i,j,k) * gv%H_to_m
1021 hcorr = min( dh - cs%min_thickness, 0. )
1022 dh = max( dh, cs%min_thickness )
1023 cellheight(k) = ifaceheight(k) - 0.5 * dh
1024 ifaceheight(k+1) = ifaceheight(k) - dh
1027 sldepth_0d = cs%surf_layer_ext*max( max(-cellheight(k),-ifaceheight(2) ), cs%minOBLdepth )
1030 if (-1.0*ifaceheight(ktmp+1) >= sldepth_0d)
then 1048 delh = min( max(0.0, sldepth_0d - htot), h(i,j,ktmp)*gv%H_to_m )
1054 surfhtemp = surfhtemp + temp(i,j,ktmp) * delh
1055 surfhsalt = surfhsalt + salt(i,j,ktmp) * delh
1056 surfhu = surfhu + 0.5*us%L_T_to_m_s*(u(i,j,ktmp)+u(i-1,j,ktmp)) * delh
1057 surfhv = surfhv + 0.5*us%L_T_to_m_s*(v(i,j,ktmp)+v(i,j-1,ktmp)) * delh
1058 if (cs%Stokes_Mixing)
then 1059 surfhus = surfhus + 0.5*(waves%US_x(i,j,ktmp)+waves%US_x(i-1,j,ktmp)) * delh
1060 surfhvs = surfhvs + 0.5*(waves%US_y(i,j,ktmp)+waves%US_y(i,j-1,ktmp)) * delh
1064 surftemp = surfhtemp / htot
1065 surfsalt = surfhsalt / htot
1066 surfu = surfhu / htot
1067 surfv = surfhv / htot
1068 surfus = surfhus / htot
1069 surfvs = surfhvs / htot
1074 uk = 0.5*us%L_T_to_m_s*(u(i,j,k)+u(i-1,j,k)) - surfu
1075 vk = 0.5*us%L_T_to_m_s*(v(i,j,k)+v(i,j-1,k)) - surfv
1077 if (cs%Stokes_Mixing)
then 1081 uk = uk + (0.5*(waves%Us_x(i,j,k)+waves%US_x(i-1,j,k)) -surfus )
1082 vk = vk + (0.5*(waves%Us_y(i,j,k)+waves%Us_y(i,j-1,k)) -surfvs )
1085 deltau2(k) = uk**2 + vk**2
1093 pres_1d(kk+1) = pref
1094 pres_1d(kk+2) = pref
1095 pres_1d(kk+3) = pref
1096 temp_1d(kk+1) = surftemp
1097 temp_1d(kk+2) = temp(i,j,k)
1098 temp_1d(kk+3) = temp(i,j,km1)
1099 salt_1d(kk+1) = surfsalt
1100 salt_1d(kk+2) = salt(i,j,k)
1101 salt_1d(kk+3) = salt(i,j,km1)
1105 pref = pref + (gv%g_Earth * gv%H_to_RZ) * h(i,j,k)
1108 surfbuoyflux2(k) = buoy_scale * (buoyflux(i,j,1) - buoyflux(i,j,k+1))
1112 if (cs%LT_K_ENHANCEMENT .or. cs%LT_VT2_ENHANCEMENT)
then 1113 mld_guess = max( 1.*us%m_to_Z, abs(us%m_to_Z*cs%OBLdepthprev(i,j) ) )
1114 call get_langmuir_number(la, g, gv, us, mld_guess, ustar(i,j), i, j, &
1115 h=h(i,j,:), u_h=u_h, v_h=v_h, waves=waves)
1129 deltarho(k) = rho_1d(kk+2) - rho_1d(kk+1)
1130 n2_1d(k) = (gorho * (rho_1d(kk+2) - rho_1d(kk+3)) ) / &
1131 ((0.5*(h(i,j,km1) + h(i,j,k))+gv%H_subroundoff)*gv%H_to_m)
1132 cs%N(i,j,k) = sqrt( max( n2_1d(k), 0.) )
1134 n2_1d(g%ke+1 ) = 0.0
1135 cs%N(i,j,g%ke+1 ) = 0.0
1141 call cvmix_kpp_compute_turbulent_scales( &
1142 cs%surf_layer_ext, &
1147 cvmix_kpp_params_user=cs%KPP_params )
1150 cs%Vt2(i,j,:) = cvmix_kpp_compute_unresolved_shear( &
1151 zt_cntr=cellheight(1:g%ke), &
1153 n_iface=cs%N(i,j,:), &
1154 cvmix_kpp_params_user=cs%KPP_params )
1157 IF (cs%LT_VT2_ENHANCEMENT)
then 1158 IF (cs%LT_VT2_METHOD==lt_vt2_mode_constant)
then 1160 langenhvt2(k) = cs%KPP_VT2_ENH_FAC
1162 elseif (cs%LT_VT2_METHOD==lt_vt2_mode_vr12)
then 1164 enhvt2 = sqrt(1.+(1.5*cs%La_SL(i,j))**(-2) + &
1165 (5.4*cs%La_SL(i,j))**(-4))
1167 langenhvt2(k) = enhvt2
1169 elseif (cs%LT_VT2_METHOD==lt_vt2_mode_rw16)
then 1171 enhvt2 = 1. + 2.3*cs%La_SL(i,j)**(-0.5)
1173 langenhvt2(k) = enhvt2
1175 elseif (cs%LT_VT2_METHOD==lt_vt2_mode_lf17)
then 1176 cs%CS=cvmix_get_kpp_real(
'c_s',cs%KPP_params)
1178 wst = (max(0.,-buoy_scale*buoyflux(i,j,1))*(-cellheight(k)))**(1./3.)
1179 langenhvt2(k) = sqrt((0.15*wst**3. + 0.17*surffricvel**3.* &
1180 (1.+0.49*cs%La_SL(i,j)**(-2.))) / &
1181 (0.2*ws_1d(k)**3/(cs%cs*cs%surf_layer_ext*cs%vonKarman**4.)))
1193 cs%Vt2(i,j,k)=cs%Vt2(i,j,k)*langenhvt2(k)
1194 if (cs%id_EnhVt2 > 0) cs%EnhVt2(i,j,k)=langenhvt2(k)
1198 bulkri_1d = cvmix_kpp_compute_bulk_richardson( &
1199 zt_cntr = cellheight(1:g%ke), &
1200 delta_buoy_cntr=gorho*deltarho, &
1201 delta_vsqr_cntr=deltau2, &
1202 vt_sqr_cntr=cs%Vt2(i,j,:), &
1204 n_iface=cs%N(i,j,:))
1207 surfbuoyflux = buoy_scale * buoyflux(i,j,1)
1210 call cvmix_kpp_compute_obl_depth( &
1215 zt_cntr=cellheight, &
1216 surf_fric=surffricvel, &
1217 surf_buoy=surfbuoyflux, &
1218 coriolis=coriolis, &
1219 cvmix_kpp_params_user=cs%KPP_params )
1223 if (cs%deepOBLoffset>0.)
then 1224 zbottomminusoffset = ifaceheight(g%ke+1) + min(cs%deepOBLoffset,-0.1*ifaceheight(g%ke+1))
1225 cs%OBLdepth(i,j) = min( cs%OBLdepth(i,j), -zbottomminusoffset )
1229 if(cs%fixedOBLdepth) cs%OBLdepth(i,j) = cs%fixedOBLdepth_value
1230 cs%OBLdepth(i,j) = max( cs%OBLdepth(i,j), -ifaceheight(2) )
1231 cs%OBLdepth(i,j) = min( cs%OBLdepth(i,j), -ifaceheight(g%ke+1) )
1232 cs%kOBL(i,j) = cvmix_kpp_compute_kobl_depth( ifaceheight, cellheight, cs%OBLdepth(i,j) )
1237 if (cs%id_Ws > 0)
then 1238 call cvmix_kpp_compute_turbulent_scales( &
1239 -cellheight/cs%OBLdepth(i,j), &
1244 cvmix_kpp_params_user=cs%KPP_params)
1245 cs%Ws(i,j,:) = ws_1d(:)
1249 if (cs%id_N2 > 0) cs%N2(i,j,:) = n2_1d(:)
1250 if (cs%id_BulkDrho > 0) cs%dRho(i,j,:) = deltarho(:)
1251 if (cs%id_BulkRi > 0) cs%BulkRi(i,j,:) = bulkri_1d(:)
1252 if (cs%id_BulkUz2 > 0) cs%Uz2(i,j,:) = deltau2(:)
1253 if (cs%id_Tsurf > 0) cs%Tsurf(i,j) = surftemp
1254 if (cs%id_Ssurf > 0) cs%Ssurf(i,j) = surfsalt
1255 if (cs%id_Usurf > 0) cs%Usurf(i,j) = surfu
1256 if (cs%id_Vsurf > 0) cs%Vsurf(i,j) = surfv
1261 call cpu_clock_end(id_clock_kpp_compute_bld)
1264 if (cs%id_BulkRi > 0)
call post_data(cs%id_BulkRi, cs%BulkRi, cs%diag)
1265 if (cs%id_N > 0)
call post_data(cs%id_N, cs%N, cs%diag)
1266 if (cs%id_N2 > 0)
call post_data(cs%id_N2, cs%N2, cs%diag)
1267 if (cs%id_Tsurf > 0)
call post_data(cs%id_Tsurf, cs%Tsurf, cs%diag)
1268 if (cs%id_Ssurf > 0)
call post_data(cs%id_Ssurf, cs%Ssurf, cs%diag)
1269 if (cs%id_Usurf > 0)
call post_data(cs%id_Usurf, cs%Usurf, cs%diag)
1270 if (cs%id_Vsurf > 0)
call post_data(cs%id_Vsurf, cs%Vsurf, cs%diag)
1271 if (cs%id_BulkDrho > 0)
call post_data(cs%id_BulkDrho, cs%dRho, cs%diag)
1272 if (cs%id_BulkUz2 > 0)
call post_data(cs%id_BulkUz2, cs%Uz2, cs%diag)
1273 if (cs%id_EnhK > 0)
call post_data(cs%id_EnhK, cs%EnhK, cs%diag)
1274 if (cs%id_EnhVt2 > 0)
call post_data(cs%id_EnhVt2, cs%EnhVt2, cs%diag)
1275 if (cs%id_La_SL > 0)
call post_data(cs%id_La_SL, cs%La_SL, cs%diag)
1278 if (cs%n_smooth > 0)
call kpp_smooth_bld(cs,g,gv,h)
1280 end subroutine kpp_compute_bld
1284 subroutine kpp_smooth_bld(CS,G,GV,h)
1286 type(KPP_CS),
pointer :: CS
1287 type(ocean_grid_type),
intent(inout) :: G
1288 type(verticalGrid_type),
intent(in) :: GV
1289 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: h
1292 real,
dimension(SZI_(G),SZJ_(G)) :: OBLdepth_prev
1293 real,
dimension( G%ke ) :: cellHeight
1295 real,
dimension( G%ke+1 ) :: iFaceHeight
1297 real :: wc, ww, we, wn, ws
1300 integer :: i, j, k, s
1302 call cpu_clock_begin(id_clock_kpp_smoothing)
1305 call pass_var(cs%OBLdepth, g%Domain, halo=cs%n_smooth)
1307 if (cs%id_OBLdepth_original > 0) cs%OBLdepth_original = cs%OBLdepth
1311 obldepth_prev = cs%OBLdepth
1320 if (g%mask2dT(i,j)==0.) cycle
1322 ifaceheight(1) = 0.0
1327 dh = h(i,j,k) * gv%H_to_m
1329 hcorr = min( dh - cs%min_thickness, 0. )
1330 dh = max( dh, cs%min_thickness )
1331 cellheight(k) = ifaceheight(k) - 0.5 * dh
1332 ifaceheight(k+1) = ifaceheight(k) - dh
1336 ww = 0.125 * g%mask2dT(i-1,j)
1337 we = 0.125 * g%mask2dT(i+1,j)
1338 ws = 0.125 * g%mask2dT(i,j-1)
1339 wn = 0.125 * g%mask2dT(i,j+1)
1340 wc = 1.0 - (ww+we+wn+ws)
1342 cs%OBLdepth(i,j) = wc * obldepth_prev(i,j) &
1343 + ww * obldepth_prev(i-1,j) &
1344 + we * obldepth_prev(i+1,j) &
1345 + ws * obldepth_prev(i,j-1) &
1346 + wn * obldepth_prev(i,j+1)
1349 if (cs%deepen_only) cs%OBLdepth(i,j) = max(cs%OBLdepth(i,j), obldepth_prev(i,j))
1352 cs%OBLdepth(i,j) = min( cs%OBLdepth(i,j), -ifaceheight(g%ke+1) )
1353 cs%kOBL(i,j) = cvmix_kpp_compute_kobl_depth( ifaceheight, cellheight, cs%OBLdepth(i,j) )
1359 call cpu_clock_end(id_clock_kpp_smoothing)
1361 end subroutine kpp_smooth_bld
1366 subroutine kpp_get_bld(CS, BLD, G, US, m_to_BLD_units)
1371 real,
dimension(SZI_(G),SZJ_(G)),
intent(inout) :: bld
1372 real,
optional,
intent(in) :: m_to_bld_units
1378 scale = us%m_to_Z ;
if (
present(m_to_bld_units)) scale = m_to_bld_units
1381 do j = g%jsc, g%jec ;
do i = g%isc, g%iec
1382 bld(i,j) = scale * cs%OBLdepth(i,j)
1385 end subroutine kpp_get_bld
1388 subroutine kpp_nonlocaltransport_temp(CS, G, GV, h, nonLocalTrans, surfFlux, &
1391 type(
kpp_cs),
intent(in) :: cs
1394 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: h
1395 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1),
intent(in) :: nonlocaltrans
1396 real,
dimension(SZI_(G),SZJ_(G)),
intent(in) :: surfflux
1398 real,
intent(in) :: dt
1399 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(inout) :: scalar
1400 real,
intent(in) :: c_p
1403 real,
dimension( SZI_(G), SZJ_(G), SZK_(G) ) :: dtracer
1406 dtracer(:,:,:) = 0.0
1411 dtracer(i,j,k) = ( nonlocaltrans(i,j,k) - nonlocaltrans(i,j,k+1) ) / &
1412 ( h(i,j,k) + gv%H_subroundoff ) * surfflux(i,j)
1418 if (cs%applyNonLocalTrans)
then 1423 scalar(i,j,k) = scalar(i,j,k) + dt * dtracer(i,j,k)
1430 if (cs%id_QminusSW > 0)
call post_data(cs%id_QminusSW, surfflux, cs%diag)
1431 if (cs%id_NLT_dTdt > 0)
call post_data(cs%id_NLT_dTdt, dtracer, cs%diag)
1432 if (cs%id_NLT_temp_budget > 0)
then 1433 dtracer(:,:,:) = 0.0
1438 dtracer(i,j,k) = (nonlocaltrans(i,j,k) - nonlocaltrans(i,j,k+1)) * &
1439 surfflux(i,j) * c_p * gv%H_to_kg_m2
1443 call post_data(cs%id_NLT_temp_budget, dtracer, cs%diag)
1446 end subroutine kpp_nonlocaltransport_temp
1451 subroutine kpp_nonlocaltransport_saln(CS, G, GV, h, nonLocalTrans, surfFlux, dt, scalar)
1453 type(
kpp_cs),
intent(in) :: cs
1456 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(in) :: h
1457 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1),
intent(in) :: nonlocaltrans
1458 real,
dimension(SZI_(G),SZJ_(G)),
intent(in) :: surfflux
1460 real,
intent(in) :: dt
1461 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)),
intent(inout) :: scalar
1464 real,
dimension( SZI_(G), SZJ_(G), SZK_(G) ) :: dtracer
1467 dtracer(:,:,:) = 0.0
1472 dtracer(i,j,k) = ( nonlocaltrans(i,j,k) - nonlocaltrans(i,j,k+1) ) / &
1473 ( h(i,j,k) + gv%H_subroundoff ) * surfflux(i,j)
1479 if (cs%applyNonLocalTrans)
then 1484 scalar(i,j,k) = scalar(i,j,k) + dt * dtracer(i,j,k)
1491 if (cs%id_netS > 0)
call post_data(cs%id_netS, surfflux, cs%diag)
1492 if (cs%id_NLT_dSdt > 0)
call post_data(cs%id_NLT_dSdt, dtracer, cs%diag)
1493 if (cs%id_NLT_saln_budget > 0)
then 1494 dtracer(:,:,:) = 0.0
1499 dtracer(i,j,k) = (nonlocaltrans(i,j,k) - nonlocaltrans(i,j,k+1)) * &
1500 surfflux(i,j) * gv%H_to_kg_m2
1504 call post_data(cs%id_NLT_saln_budget, dtracer, cs%diag)
1507 end subroutine kpp_nonlocaltransport_saln
1513 subroutine kpp_end(CS)
1516 if (.not.
associated(cs))
return 1520 end subroutine kpp_end
Ocean grid type. See mom_grid for details.
Calculates density of sea water from T, S and P.
A structure that can be parsed to read and document run-time parameters.
Provides the ocean grid type.
Wraps the MPP cpu clock functions.
The MOM6 facility to parse input files for runtime parameters.
An overloaded interface to log the values of various types of parameters.
Interfaces to non-domain-oriented communication subroutines, including the MOM6 reproducing sums faci...
Interface for surface waves.
Make a diagnostic available for averaging or output.
Describes various unit conversion factors.
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Describes the decomposed MOM domain and has routines for communications across PEs.
Routines for error handling and I/O management.
Container for all surface wave related parameters.
Control structure for containing KPP parameters/data.
Provides subroutines for quantities specific to the equation of state.
An overloaded interface to log version information about modules.
Describes the vertical ocean grid, including unit conversion factors.
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
Provides a transparent vertical ocean grid type and supporting routines.
Provides transparent structures with groups of MOM6 variables and supporting routines.
Provides the K-Profile Parameterization (KPP) of Large et al., 1994, via CVMix.
Do a halo update on an array.
An overloaded interface to read and log the values of various types of parameters.
Provides checksumming functions for debugging.
A control structure for the equation of state.