The dynamic psi-equation


INTERFACE:

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

This model has been formulated by Umlauf and Burchard (2003), who introduced a `generic' variable,

$\displaystyle \psi = (c_\mu^0)^p k^m l^n \; , \quad$ (167)

where $ k$ is the turbulent kinetic energy computed from (152) and $ l$ is the dissipative length-scale defined in (155). For appropriate choices of the exponents $ p$, $ m$, and $ n$, the variable $ \psi$ can be directly identified with the classic length-scale determining variables like the rate of dissipation, $ \epsilon$, or the product $ kl$ used by Mellor and Yamada (1982) (see section 4.14 and section 4.15). Some examples are compiled in table 8.

Table 8: Exponents $ p$, $ n$, $ m$ defined in (167), and their relation to the variable of the second equation in some well-known two-equation models.


The transport equation for $ \psi$ can written as

$\displaystyle \dot{\psi} = {\cal D}_\psi + \frac{\psi}{k} ( c_{\psi_1} P + c_{\psi_3} G - c_{\psi 2} \epsilon ) \; , \quad$ (168)

where $ \dot{\psi}$ denotes the material derivative of $ \psi$, see Umlauf and Burchard (2003). The production terms $ P$ and $ G$ follow from (154). $ {\cal D}_\psi$ represents the sum of the viscous and turbulent transport terms. The rate of dissipation can computed by solving (167) for $ l$ and inserting the result into (155).

For horizontally homogeneous flows, the transport terms $ {\cal D}_\psi$ appearing in (168) are expressed by a simple gradient formulation,

$\displaystyle {\cal D}_\psi = \dfrac{\partial}{\partial {z}} \left( \dfrac{\nu_t}{\sigma_\psi} \dfrac{\partial {\psi}}{\partial {z}} \right) \quad .$ (169)

For appropriate choices of the parameters, most of the classic transport equations can be directly recovered from the generic equation (168). An example is the transport equation for the inverse turbulent time scale, $ \omega \propto \epsilon / k$, which has been formulated by Wilcox (1988) and extended to buoyancy affected flows by Umlauf et al. (2003). The precise definition of $ \omega$ follows from table 8, and its transport equation can be written as

$\displaystyle \dot{\omega} = {\cal D}_\omega + \frac{\omega}{k} ( c_{\omega_1} P + c_{\omega_3} G - c_{\omega 2} \epsilon ) \; , \quad$ (170)

which is clearly a special case of (168). Model constants for this and other traditional models are given in table 9.

Table 9: Model constants of some standard models, converted to the notation used here. The Schmidt-numbers for the model of Mellor and Yamada (1982) are valid only in the logarithmic boundary-layer, because the diffusion models (157) and (164) are slightly different from (153) and (169). There is no indication that one class of diffusion models is superior.
  $ c_\mu^0$ $ \sigma_k^\psi$ $ \sigma_\psi$ $ c_{\psi 1}$ $ c_{\psi 2}$ $ c_{\psi 3}$  
$ k$-$ \epsilon$, Rodi (1987) : $ 0.5477$ $ 1.0$ $ 1.3$ $ 1.44$ $ 1.92$ (see eq. (116))  
$ k$-$ kl$, Mellor and Yamada (1982) : $ 0.5544$ $ 1.96$ $ 1.96$ $ 0.9$ $ 0.5$ $ 0.9$  
$ k$-$ \omega$, Wilcox (1988) : $ 0.5477$ $ 2$ $ 2$ $ 0.555$ $ 0.833$ (see eq. (116))  
$ k$-$ \tau$ Zeierman and Wolfshtein (1986): $ 0.5477$ $ 1.46$ $ 10.8$ $ 0.173$ $ 0.225$ (--)  


Apart from having to code only one equation to recover all of the traditional models, the main advantage of the generic equation is its flexibility. After choosing meaningful values for physically relevant parameters like the von Kármán constant, $ \kappa$, the temporal decay rate for homogeneous turbulence, $ d$, some parameters related to breaking surface waves, etc, a two-equation model can be generated, which has exactly the required properties. This is discussed in great detail in Umlauf and Burchard (2003). All algorithms have been implemented in GOTM and are described in section 4.7.3.


USES:

   use turbulence, only: P,B,num
   use turbulence, only: tke,tkeo,k_min,eps,eps_min,L
   use turbulence, only: cpsi1,cpsi2,cpsi3plus,cpsi3minus,sig_psi
   use turbulence, only: gen_m,gen_n,gen_p
   use turbulence, only: cm0,cde,galp,length_lim
   use turbulence, only: psi_bc, psi_ubc, psi_lbc, ubc_type, lbc_type
   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 and Hans Burchard

Karsten Bolding 2012-12-28