vv_momentum_3d - $y$-momentum eq. (Source File: vv_momentum_3d.F90)

INTERFACE:

    subroutine vv_momentum_3d(n,bdy3d)
DESCRIPTION:

Here, the budget equation for layer-averaged momentum in northern direction, $q_k$, is calculated. The physical equation is given as equation (2), the layer-integrated equation as (27), and after curvilinear transformation as (39). In this routine, first the Coriolis rotation term, $fp_k$ is calculated, either as direct transport averaging, or following Espelid et al. (2000) by using velocity averages (in case the compiler option NEW_CORI is set).

As a next step, explicit forcing terms (advection, diffusion, internal pressure gradient, surface stresses) are added up (into the variable ex(k)), the eddy viscosity is horizontally interpolated to the V-point, and the barotropic pressure gradient is calculated (the latter includes the pressure gradient correction for drying points, see section 5.5). Afterwards, the matrix is set up for each water column, and it is solved by means of a tri-diagonal matrix solver.

In case that the compiler option STRUCTURE_FRICTION is switched on, the frictional effect of structures in the water column is calculated by adding the quadratic frictional term $C v \sqrt{u^2+v^2}$ (with a minus sign on the right hand side) numerically implicitly to the $v$-equation, with the friction coefficient $C$. The explicit part of this term, $C\sqrt{u^2+v^2}$, is calculated in the routine structure_friction_3d.F90.

Finally, the new velocity profile is shifted such that its vertical integral is identical to the time integral of the vertically integrated transport. If the compiler option MUDFLAT is defined, this fitting of profiles is made with respect to the new surface elevation, otherwise to the old surface elevation.

When GETM is run as a slice model (compiler option SLICE_MODEL is activated), the result for $j=2$ is copied to $j=1$ and $j=3$. USES:

    use exceptions
    use parameters, only: g,avmmol,rho_0
    use domain, only: imin,imax,jmin,jmax,kmax,H,HV,min_depth
    use domain, only: dry_v,corv,au,av,az
 #if defined CURVILINEAR || defined SPHERICAL
    use domain, only: dyv,arvd1,dxc,dyx,dyc,dxx
 #else
    use domain, only: dx,dy
 #endif
    use variables_2d, only: Vint,D
    use bdy_3d, only: do_bdy_3d
    use variables_3d, only: dt,cnpar,kvmin,uu,vv,huo,hvo,hvn,vvEx,ww,hun
    use variables_3d, only: num,nuh,sseo,ssvn,rrv
    use variables_3d, only: ssvo
 #ifdef _MOMENTUM_TERMS_
    use variables_3d, only: tdv_v,cor_v,ipg_v,epg_v,vsd_v,hsd_v,adv_v
 #endif
 #ifdef STRUCTURE_FRICTION
    use variables_3d, only: sf
 #endif
 #ifndef NO_BAROCLINIC
    use variables_3d, only: idpdy
 #endif
    use halo_zones, only: update_3d_halo,wait_halo,V_TAG
    use meteo, only: tausy,airp
    use m3d, only: ip_fac
    use m3d, only: vel_check,min_vel,max_vel
    use getm_timers, only: tic, toc, TIM_VVMOMENTUM, TIM_VVMOMENTUMH
  $ use omp_lib
    IMPLICIT NONE
INPUT PARAMETERS:
    integer, intent(in)                 :: n
    logical, intent(in)                 :: bdy3d
REVISION HISTORY:
    Original author(s): Hans Burchard & Karsten Bolding
LOCAL VARIABLES:
    integer                   :: i,j,k,rc
 #ifdef NEW_CORI
    REALTYPE,dimension(I3DFIELD) :: work3d
 #endif
    REALTYPE, POINTER         :: dif(:)
    REALTYPE, POINTER         :: auxn(:),auxo(:)
    REALTYPE, POINTER         :: a1(:),a2(:)
    REALTYPE, POINTER         :: a3(:),a4(:)
    REALTYPE, POINTER         :: Res(:),ex(:)
    REALTYPE                  :: zp,zm,zy,ResInt,Diff,Uloc
    REALTYPE                  :: gamma=g*rho_0
    REALTYPE                  :: cord_curv=_ZERO_
    REALTYPE                  :: gammai,rho_0i
    integer                   :: status