#include "cppdefs.h" MODULE sed_surface_mod #if defined NONLINEAR && defined SEDIMENT ! !git $Id$ !svn $Id: sed_surface.F 1151 2023-02-09 03:08:53Z arango $ !==================================================== John C. Warner === ! Copyright (c) 2002-2023 The ROMS/TOMS Group Hernan G. Arango ! ! Licensed under a MIT/X style license ! ! See License_ROMS.md ! !======================================================================= ! ! ! This routine computed sediment surface layer (sediment-water ! ! interface) properties. ! ! ! ! References: ! ! ! ! Warner, J.C., C.R. Sherwood, R.P. Signell, C.K. Harris, and H.G. ! ! Arango, 2008: Development of a three-dimensional, regional, ! ! coupled wave, current, and sediment-transport model, Computers ! ! & Geosciences, 34, 1284-1306. ! ! ! !======================================================================= ! implicit none PRIVATE PUBLIC :: sed_surface CONTAINS ! !*********************************************************************** SUBROUTINE sed_surface (ng, tile) !*********************************************************************** ! USE mod_param USE mod_ocean USE mod_sedbed USE mod_stepping ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile ! ! Local variable declarations. ! character (len=*), parameter :: MyFile = & & __FILE__ ! # include "tile.h" ! # ifdef PROFILE CALL wclock_on (ng, iNLM, 16, __LINE__, MyFile) # endif CALL sed_surface_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & nstp(ng), nnew(ng), & & SEDBED(ng) % bed_frac, & & SEDBED(ng) % bottom) # ifdef PROFILE CALL wclock_off (ng, iNLM, 16, __LINE__, MyFile) # endif ! RETURN END SUBROUTINE sed_surface ! !*********************************************************************** SUBROUTINE sed_surface_tile (ng, tile, & & LBi, UBi, LBj, UBj, & & IminS, ImaxS, JminS, JmaxS, & & nstp, nnew, & & bed_frac, bottom) !*********************************************************************** ! USE mod_param USE mod_scalars USE mod_sediment ! USE bc_3d_mod, ONLY : bc_r3d_tile # ifdef DISTRIBUTE USE mp_exchange_mod, ONLY : mp_exchange3d, mp_exchange4d # endif ! ! Imported variable declarations. ! integer, intent(in) :: ng, tile integer, intent(in) :: LBi, UBi, LBj, UBj integer, intent(in) :: IminS, ImaxS, JminS, JmaxS integer, intent(in) :: nstp, nnew ! # ifdef ASSUMED_SHAPE real(r8), intent(inout) :: bed_frac(LBi:,LBj:,:,:) real(r8), intent(inout) :: bottom(LBi:,LBj:,:) # else real(r8), intent(inout) :: bed_frac(LBi:UBi,LBj:UBj,Nbed,NST) real(r8), intent(inout) :: bottom(LBi:UBi,LBj:UBj,MBOTP) # endif ! ! Local variable declarations. ! integer :: i, ised, j real(r8) :: cff1, cff2, cff3, cff4 real(r8), parameter :: eps = 1.0E-14_r8 # include "set_bounds.h" ! ! Update mean surface properties. ! Sd50 must be positive definite, due to BBL routines. ! Srho must be >1000, due to (s-1) in BBL routines ! J_LOOP : DO j=Jstr,Jend DO i=Istr,Iend cff1=1.0_r8 cff2=1.0_r8 cff3=1.0_r8 cff4=1.0_r8 DO ised=1,NST cff1=cff1*tau_ce(ised,ng)**bed_frac(i,j,1,ised) cff2=cff2*Sd50(ised,ng)**bed_frac(i,j,1,ised) cff3=cff3*(wsed(ised,ng)+eps)**bed_frac(i,j,1,ised) cff4=cff4*Srho(ised,ng)**bed_frac(i,j,1,ised) END DO bottom(i,j,itauc)=cff1 bottom(i,j,isd50)=MIN(cff2,Zob(ng)) bottom(i,j,iwsed)=cff3 bottom(i,j,idens)=MAX(cff4,1050.0_r8) END DO END DO J_LOOP !! !! Move bottom roughness computations from BBL routines to here. !! ! !----------------------------------------------------------------------- ! Apply periodic or gradient boundary conditions to property arrays. !----------------------------------------------------------------------- ! CALL bc_r3d_tile (ng, tile, & & LBi, UBi, LBj, UBj, 1, MBOTP, & & bottom) # ifdef DISTRIBUTE CALL mp_exchange3d (ng, tile, iNLM, 1, & & LBi, UBi, LBj, UBj, 1, MBOTP, & & NghostPoints, & & EWperiodic(ng), NSperiodic(ng), & & bottom) # endif ! RETURN END SUBROUTINE sed_surface_tile #endif END MODULE sed_surface_mod