The dynamic k-equation (Source File: tkeeq.F90)


INTERFACE:

    subroutine tkeeq(nlev,dt,u_taus,u_taub,z0s,z0b,h,NN,SS)
DESCRIPTION:

The transport equation for the turbulent kinetic energy, $ k$, follows immediately from the contraction of the Reynolds-stress tensor. In the case of a Boussinesq-fluid, this equation can be written as

$\displaystyle \dot{k} = {\cal D}_k + P + G - \epsilon \; , \quad$ (150)

where $ \dot{k}$ denotes the material derivative of $ k$. $ P$ and $ G$ are the production of $ k$ by mean shear and buoyancy, respectively, and $ \epsilon$ the rate of dissipation. $ {\cal D}_k$ represents the sum of the viscous and turbulent transport terms. For horizontally homogeneous flows, the transport term $ {\cal D}_k$ appearing in (150) is presently expressed by a simple gradient formulation,

$\displaystyle {\cal D}_k = \dfrac{\partial}{\partial {z}} \left( \dfrac{\nu_t}{\sigma_k} \dfrac{\partial {k}}{\partial {z}} \right) \; , \quad$ (151)

where $ \sigma_k$ is the constant Schmidt-number for $ k$.

In horizontally homogeneous flows, the shear and the buoyancy production, $ P$ and $ G$, can be written as

\begin{displaymath}\begin{array}{rcl} P &=& - {\langle u'w' \rangle} \dfrac{\par...
...ad \\ [3mm] G &=& {\langle w'b' \rangle} \; , \quad \end{array}\end{displaymath} (152)

see (48). Their computation is discussed in section 4.7.21.

The rate of dissipation, $ \epsilon$, can be either obtained directly from its parameterised transport equation as discussed in section 4.7.28, or from any other model yielding an appropriate description of the dissipative length-scale, $ l$. Then, $ \epsilon$ follows from the well-known cascading relation of turbulence,

$\displaystyle \epsilon = (c_\mu^0)^3 \frac{k^{\frac{3}{2}}}{l} \; , \quad$ (153)

where $ c_\mu^0$ is a constant of the model.


USES:

    use turbulence,   only: P,B,num
    use turbulence,   only: tke,tkeo,k_min,eps
    use turbulence,   only: k_bc, k_ubc, k_lbc, ubc_type, lbc_type
    use turbulence,   only: sig_k
    use util,         only: Dirichlet,Neumann
 
    IMPLICIT NONE
INPUT PARAMETERS:
 
    number of vertical layers
    integer,  intent(in)                :: nlev
 
    time step (s)
    REALTYPE, intent(in)                :: dt
 
    surface and bottom
    friction velocity (m/s)
    REALTYPE, intent(in)                :: u_taus,u_taub
 
    surface and bottom
    roughness length (m)
    REALTYPE, intent(in)                :: z0s,z0b
 
    layer thickness (m)
    REALTYPE, intent(in)                :: h(0:nlev)
 
    square of shear and buoyancy
    frequency (1/s^2)
    REALTYPE, intent(in)                :: NN(0:nlev),SS(0:nlev)
REVISION HISTORY:
    Original author(s): Lars Umlauf
                       (re-write after first version of
                        H. Burchard and K. Bolding)
    $Log: tkeeq.F90,v $
    Revision 1.9  2005-11-15 11:35:02  lars
    documentation finish for print
    Revision 1.8  2005/11/03 20:53:37  hb
    Patankar trick reverted to older versions for
    stabilising 3D computations
    Revision 1.7  2005/08/11 13:11:50  lars
    Added explicit loops for diffusivities for 3-D z-level support.
    Thanks to Vicente Fernandez.
    Revision 1.6  2005/06/27 13:44:07  kbk
    modified + removed traling blanks
    Revision 1.5  2003/03/28 09:20:35  kbk
    added new copyright to files
    Revision 1.4  2003/03/10 09:02:06  gotm
    Added new Generic Turbulence Model + improved documentation and cleaned up code

Karsten Bolding 2012-01-24