From c1c5c5711f3c78b756a023c422a6aeed001c6b0f Mon Sep 17 00:00:00 2001 From: Robert Hallberg Date: Wed, 24 Jan 2024 09:26:50 -0500 Subject: [PATCH] (*)Use thickness_to_dz with BRYAN_LEWIS_DIFFUSIVITY Use thickness_to_dz to convert the layer thicknesses into the geometric depths used to set the Bryan-Lewis background diffusivities instead of multiplication by GV%H_to_m, thereby avoiding division by the Boussinesq reference density in non-Boussinesq mode. No answers are changed in any Boussinsesq configurations, but answers will change in non-Boussinesq configurations that have BRYAN_LEWIS_DIFFUSIVITY = True. --- src/parameterizations/vertical/MOM_bkgnd_mixing.F90 | 6 +++++- 1 file changed, 5 insertions(+), 1 deletion(-) diff --git a/src/parameterizations/vertical/MOM_bkgnd_mixing.F90 b/src/parameterizations/vertical/MOM_bkgnd_mixing.F90 index 693b9395bd..ee04f3f195 100644 --- a/src/parameterizations/vertical/MOM_bkgnd_mixing.F90 +++ b/src/parameterizations/vertical/MOM_bkgnd_mixing.F90 @@ -13,6 +13,7 @@ module MOM_bkgnd_mixing use MOM_file_parser, only : openParameterBlock, closeParameterBlock use MOM_forcing_type, only : forcing use MOM_grid, only : ocean_grid_type +use MOM_interface_heights, only : thickness_to_dz use MOM_unit_scaling, only : unit_scale_type use MOM_verticalGrid, only : verticalGrid_type use MOM_variables, only : thermo_var_ptrs, vertvisc_type @@ -338,6 +339,7 @@ subroutine calculate_bkgnd_mixing(h, tv, N2_lay, Kd_lay, Kd_int, Kv_bkgnd, j, G, real, dimension(SZK_(GV)+1) :: Kv_col !< Viscosities at the interfaces [m2 s-1] real, dimension(SZI_(G)) :: Kd_sfc !< Surface value of the diffusivity [H Z T-1 ~> m2 s-1 or kg m-1 s-1] real, dimension(SZI_(G)) :: depth !< Distance from surface of an interface [H ~> m or kg m-2] + real, dimension(SZI_(G),SZK_(GV)) :: dz !< Height change across layers [Z ~> m] real :: depth_c !< depth of the center of a layer [H ~> m or kg m-2] real :: I_Hmix !< inverse of fixed mixed layer thickness [H-1 ~> m-1 or m2 kg-1] real :: I_2Omega !< 1/(2 Omega) [T ~> s] @@ -364,10 +366,12 @@ subroutine calculate_bkgnd_mixing(h, tv, N2_lay, Kd_lay, Kd_int, Kv_bkgnd, j, G, ! Set up the background diffusivity. if (CS%Bryan_Lewis_diffusivity) then + call thickness_to_dz(h, tv, dz, j, G, GV) + do i=is,ie depth_int(1) = 0.0 do k=2,nz+1 - depth_int(k) = depth_int(k-1) + GV%H_to_m*h(i,j,k-1) + depth_int(k) = depth_int(k-1) + US%Z_to_m*dz(i,k-1) enddo call CVMix_init_bkgnd(max_nlev=nz, &