Skip to content

Commit

Permalink
Add ability to read external file distribution for KH.
Browse files Browse the repository at this point in the history
  • Loading branch information
wfcooke committed Nov 30, 2016
1 parent 2d66f99 commit 1f1dd99
Showing 1 changed file with 33 additions and 0 deletions.
33 changes: 33 additions & 0 deletions src/parameterizations/lateral/MOM_hor_visc.F90
Original file line number Diff line number Diff line change
Expand Up @@ -88,13 +88,15 @@ module MOM_hor_visc

use MOM_diag_mediator, only : post_data, register_diag_field, safe_alloc_ptr
use MOM_diag_mediator, only : diag_ctrl, time_type
use MOM_domains, only : pass_var
use MOM_error_handler, only : MOM_error, FATAL, WARNING
use MOM_file_parser, only : get_param, log_version, param_file_type
use MOM_grid, only : ocean_grid_type
use MOM_lateral_mixing_coeffs, only : VarMix_CS
use MOM_MEKE_types, only : MEKE_type
use MOM_open_boundary, only : ocean_OBC_type, OBC_FLATHER
use MOM_verticalGrid, only : verticalGrid_type
use fms_mod, only : read_data

implicit none ; private

Expand Down Expand Up @@ -131,11 +133,17 @@ module MOM_hor_visc
logical :: bound_Coriolis ! If true & SMAGORINSKY_AH is used, the biharmonic
! viscosity is modified to include a term that
! scales quadratically with the velocity shears.
logical :: use_Kh_bg_2d ! Read 2d background viscosity from a file.
real :: Kh_bg_min ! The minimum value allowed for Laplacian horizontal
! viscosity. The default is 0.0

real ALLOCABLE_, dimension(NIMEM_,NJMEM_) :: &
Kh_bg_xx, &! The background Laplacian viscosity at h points, in units
! of m2 s-1. The actual viscosity may be the larger of this
! viscosity and the Smagorinsky viscosity.
Kh_bg_2d, &! The background Laplacian viscosity at h points, in units
! of m2 s-1. The actual viscosity may be the larger of this
! viscosity and the Smagorinsky viscosity.
Ah_bg_xx, &! The background biharmonic viscosity at h points, in units
! of m4 s-1. The actual viscosity may be the larger of this
! viscosity and the Smagorinsky viscosity.
Expand Down Expand Up @@ -518,6 +526,9 @@ subroutine horizontal_viscosity(u, v, h, diffu, diffv, MEKE, VarMix, G, GV, CS,
Kh = Kh + 0.25*( (MEKE%Ku(I,J)+MEKE%Ku(I+1,J+1)) &
+(MEKE%Ku(I+1,J)+MEKE%Ku(I,J+1)) )
endif
! Place a floor on the viscosity, if desired.
Kh = MAX(Kh,CS%Kh_bg_min)

if (CS%better_bound_Kh) then
if (Kh >= hrat_min*CS%Kh_Max_xy(I,J)) then
visc_bound_rem = 0.0
Expand Down Expand Up @@ -785,6 +796,9 @@ subroutine hor_visc_init(Time, G, param_file, diag, CS)
call get_param(param_file, mod, "KH", Kh, &
"The background Laplacian horizontal viscosity.", &
units = "m2 s-1", default=0.0)
call get_param(param_file, mod, "KH_BG_MIN", CS%Kh_bg_min, &
"The minimum value allowed for Laplacian horizontal viscosity, KH.", &
units = "m2 s-1", default=0.0)
call get_param(param_file, mod, "KH_VEL_SCALE", Kh_vel_scale, &
"The velocity scale which is multiplied by the grid \n"//&
"spacing to calculate the Laplacian viscosity. \n"//&
Expand Down Expand Up @@ -877,6 +891,12 @@ subroutine hor_visc_init(Time, G, param_file, diag, CS)
"is strongly encouraged, and no slip BCs are not used with \n"//&
"the biharmonic viscosity.", default=.false.)

call get_param(param_file, mod, "USE_KH_BG_2D", CS%use_Kh_bg_2d, &
"If true, read a file containing 2-d background harmonic \n"//&
"viscosities. The final viscosity is the maximum of the other "//&
"terms and this background value.", default=.false.)


if (CS%bound_Kh .or. CS%bound_Ah .or. CS%better_bound_Kh .or. CS%better_bound_Ah) &
call get_param(param_file, mod, "DT", dt, &
"The (baroclinic) dynamics time step.", units = "s", &
Expand Down Expand Up @@ -914,6 +934,13 @@ subroutine hor_visc_init(Time, G, param_file, diag, CS)
ALLOC_(CS%reduction_xx(isd:ied,jsd:jed)) ; CS%reduction_xx(:,:) = 0.0
ALLOC_(CS%reduction_xy(IsdB:IedB,JsdB:JedB)) ; CS%reduction_xy(:,:) = 0.0

if (CS%use_Kh_bg_2d) then
ALLOC_(CS%Kh_bg_2d(isd:ied,jsd:jed)) ; CS%Kh_bg_2d(:,:) = 0.0
call read_data('INPUT/KH_background_2d.nc', 'Kh', CS%Kh_bg_2d, &
domain=G%domain%mpp_domain, timelevel=1)
call pass_var(CS%Kh_bg_2d, G%domain)
endif

if (CS%biharmonic) then
ALLOC_(CS%Idx2dyCu(IsdB:IedB,jsd:jed)) ; CS%Idx2dyCu(:,:) = 0.0
ALLOC_(CS%Idx2dyCv(isd:ied,JsdB:JedB)) ; CS%Idx2dyCv(:,:) = 0.0
Expand Down Expand Up @@ -986,6 +1013,9 @@ subroutine hor_visc_init(Time, G, param_file, diag, CS)
if (CS%Smagorinsky_Kh) CS%LAPLAC_CONST_xx(i,j) = Smag_Lap_const * grid_sp_h2

CS%Kh_bg_xx(i,j) = MAX(Kh, Kh_vel_scale * sqrt(grid_sp_h2))

if (CS%use_Kh_bg_2d) CS%Kh_bg_xx(i,j) = MAX(CS%Kh_bg_2d(i,j), CS%Kh_bg_xx(i,j))

if (CS%bound_Kh .and. .not.CS%better_bound_Kh) then
CS%Kh_Max_xx(i,j) = Kh_Limit * grid_sp_h2
CS%Kh_bg_xx(i,j) = MIN(CS%Kh_bg_xx(i,j), CS%Kh_Max_xx(i,j))
Expand All @@ -996,6 +1026,9 @@ subroutine hor_visc_init(Time, G, param_file, diag, CS)
if (CS%Smagorinsky_Kh) CS%LAPLAC_CONST_xy(I,J) = Smag_Lap_const * grid_sp_q2

CS%Kh_bg_xy(I,J) = MAX(Kh, Kh_vel_scale * sqrt(grid_sp_q2))

if (CS%use_Kh_bg_2d) CS%Kh_bg_xy(i,j) = MAX(CS%Kh_bg_2d(i,j), CS%Kh_bg_xy(i,j))

if (CS%bound_Kh .and. .not.CS%better_bound_Kh) then
CS%Kh_Max_xy(I,J) = Kh_Limit * grid_sp_q2
CS%Kh_bg_xy(I,J) = MIN(CS%Kh_bg_xy(I,J), CS%Kh_Max_xy(I,J))
Expand Down

0 comments on commit 1f1dd99

Please sign in to comment.