Skip to content
This repository has been archived by the owner on Apr 18, 2018. It is now read-only.

Latest commit

 

History

History
490 lines (359 loc) · 16.9 KB

gmredi.rst

File metadata and controls

490 lines (359 loc) · 16.9 KB

GMREDI: Gent-McWilliams/Redi SGS Eddy Parameterization

There are two parts to the Redi/GM parameterization of geostrophic eddies. The first, the Redi scheme redi1982, aims to mix tracer properties along isentropes (neutral surfaces) by means of a diffusion operator oriented along the local isentropic surface. The second part, GM gen-mcw:90,gen-eta:95 , adiabatically re-arranges tracers through an advective flux where the advecting flow is a function of slope of the isentropic surfaces.

The first GCM implementation of the Redi scheme was by cox87 in the GFDL ocean circulation model. The original approach failed to distinguish between isopycnals and surfaces of locally referenced potential density (now called neutral surfaces) which are proper isentropes for the ocean. As will be discussed later, it also appears that the Cox implementation is susceptible to a computational mode. Due to this mode, the Cox scheme requires a background lateral diffusion to be present to conserve the integrity of the model fields.

The GM parameterization was then added to the GFDL code in the form of a non-divergent bolus velocity. The method defines two stream-functions expressed in terms of the isoneutral slopes subject to the boundary condition of zero value on upper and lower boundaries. The horizontal bolus velocities are then the vertical derivative of these functions. Here in lies a problem highlighted by gretal:98: the bolus velocities involve multiple derivatives on the potential density field, which can consequently give rise to noise. Griffies et al. point out that the GM bolus fluxes can be identically written as a skew flux which involves fewer differential operators. Further, combining the skew flux formulation and Redi scheme, substantial cancellations take place to the point that the horizontal fluxes are unmodified from the lateral diffusion parameterization.

Redi scheme: Isopycnal diffusion

The Redi scheme diffuses tracers along isopycnals and introduces a term in the tendency (rhs) of such a tracer (here τ) of the form:

$$\bf{\nabla} \cdot \kappa_\rho \bf{K}_{Redi} \bf{\nabla} \tau$$

where κρ is the along isopycnal diffusivity and $\bf{K}_{Redi}$ is a rank 2 tensor that projects the gradient of τ onto the isopycnal surface. The unapproximated projection tensor is:

$$\begin{aligned} \bf{K}_{Redi} = \frac{1}{1 + |S|^2} \left( \begin{array}{ccc} 1 + S_y^2& -S_x S_y & S_x \\\ -S_x S_y & 1 + S_x^2 & S_y \\\ S_x & S_y & |S|^2 \\\ \end{array} \right) \end{aligned}$$

Here, Sx =  − ∂xσ/∂zσ and Sy =  − ∂yσ/∂zσ are the components of the isoneutral slope.

The first point to note is that a typical slope in the ocean interior is small, say of the order 10 − 4. A maximum slope might be of order 10 − 2 and only exceeds such in unstratified regions where the slope is ill defined. It is therefore justifiable, and customary, to make the small slope approximation, |S| <  < 1. The Redi projection tensor then becomes:

$$\begin{aligned} \bf{K}_{Redi} = \left( \begin{array}{ccc} 1 & 0 & S_x \\\ 0 & 1 & S_y \\\ S_x & S_y & |S|^2 \\\ \end{array} \right) \end{aligned}$$

GM parameterization

The GM parameterization aims to represent the “advective” or “transport” effect of geostrophic eddies by means of a “bolus” velocity, $\bf{u}^\star$. The divergence of this advective flux is added to the tracer tendency equation (on the rhs):

$$- \bf{\nabla} \cdot \tau \bf{u}^\star$$

The bolus velocity $\bf{u}^\star$ is defined as the rotational of a streamfunction $\bf{F}^\star$=(Fx, Fy, 0):

$$\begin{aligned} \bf{u}^\star = \nabla \times \bf{F}^\star = \left( \begin{array}{c} - \partial_z F_y^\star \\\ + \partial_z F_x^\star \\\ \partial_x F_y^\star - \partial_y F_x^\star \end{array} \right), \end{aligned}$$

and thus is automatically non-divergent. In the GM parameterization, the streamfunction is specified in terms of the isoneutral slopes Sx and Sy:

$$\begin{aligned} \begin{aligned} F_x^\star & = & -\kappa_{GM} S_y \\\ F_y^\star & = & \kappa_{GM} S_x\end{aligned} \end{aligned}$$

with boundary conditions Fx = Fy = 0 on upper and lower boundaries. In the end, the bolus transport in the GM parameterization is given by:

$$\begin{aligned} \bf{u}^\star = \left( \begin{array}{c} u^\star \\\ v^\star \\\ w^\star \end{array} \right) = \left( \begin{array}{c} - \partial_z (\kappa_{GM} S_x) \\\ - \partial_z (\kappa_{GM} S_y) \\\ \partial_x (\kappa_{GM} S_x) + \partial_y (\kappa_{GM} S_y) \end{array} \right) \end{aligned}$$

This is the form of the GM parameterization as applied by Donabasaglu, 1997, in MOM versions 1 and 2.

Note that in the MITgcm, the variables containing the GM bolus streamfunction are:

$$\begin{aligned} \left( \begin{array}{c} GM\_PsiX \\\ GM\_PsiY \end{array} \right) = \left( \begin{array}{c} \kappa_{GM} S_x \\\ \kappa_{GM} S_y \end{array} \right)= \left( \begin{array}{c} F_y^\star \\\ -F_x^\star \end{array} \right). \end{aligned}$$

Griffies Skew Flux

gr:98 notes that the discretisation of bolus velocities involves multiple layers of differencing and interpolation that potentially lead to noisy fields and computational modes. He pointed out that the bolus flux can be re-written in terms of a non-divergent flux and a skew-flux:

$$\begin{aligned} \begin{aligned} \bf{u}^\star \tau & = & \left( \begin{array}{c} - \partial_z ( \kappa_{GM} S_x ) \tau \\\ - \partial_z ( \kappa_{GM} S_y ) \tau \\\ (\partial_x \kappa_{GM} S_x + \partial_y \kappa_{GM} S_y)\tau \end{array} \right) \\\ & = & \left( \begin{array}{c} - \partial_z ( \kappa_{GM} S_x \tau) \\\ - \partial_z ( \kappa_{GM} S_y \tau) \\\ \partial_x ( \kappa_{GM} S_x \tau) + \partial_y ( \kappa_{GM} S_y \tau) \end{array} \right) + \left( \begin{array}{c} \kappa_{GM} S_x \partial_z \tau \\\ \kappa_{GM} S_y \partial_z \tau \\\ - \kappa_{GM} S_x \partial_x \tau - \kappa_{GM} S_y \partial_y \tau \end{array} \right)\end{aligned} \end{aligned}$$

The first vector is non-divergent and thus has no effect on the tracer field and can be dropped. The remaining flux can be written:

$$\bf{u}^\star \tau = - \kappa_{GM} \bf{K}_{GM} \bf{\nabla} \tau$$

where

$$\begin{aligned} \bf{K}_{GM} = \left( \begin{array}{ccc} 0 & 0 & -S_x \\\ 0 & 0 & -S_y \\\ S_x & S_y & 0 \end{array} \right) \end{aligned}$$

is an anti-symmetric tensor.

This formulation of the GM parameterization involves fewer derivatives than the original and also involves only terms that already appear in the Redi mixing scheme. Indeed, a somewhat fortunate cancellation becomes apparent when we use the GM parameterization in conjunction with the Redi isoneutral mixing scheme:

$$\kappa_\rho \bf{K}_{Redi} \bf{\nabla} \tau - u^\star \tau = ( \kappa_\rho \bf{K}_{Redi} + \kappa_{GM} \bf{K}_{GM} ) \bf{\nabla} \tau$$

In the instance that κGM = κρ then

$$\begin{aligned} \kappa_\rho \bf{K}_{Redi} + \kappa_{GM} \bf{K}_{GM} = \kappa_\rho \left( \begin{array}{ccc} 1 & 0 & 0 \\\ 0 & 1 & 0 \\\ 2 S_x & 2 S_y & |S|^2 \end{array} \right) \end{aligned}$$

which differs from the variable Laplacian diffusion tensor by only two non-zero elements in the z-row.

Subroutine

S/R GMREDI_CALC_TENSOR (pkg/gmredi/gmredi_calc_tensor.F)

σx: SlopeX (argument on entry)

σy: SlopeY (argument on entry)

σz: SlopeY (argument)

Sx: SlopeX (argument on exit)

Sy: SlopeY (argument on exit)

Variable κGM

visbeck:97 suggest making the eddy coefficient, κGM, a function of the Eady growth rate, $|f|/\sqrt{Ri}$. The formula involves a non-dimensional constant, α, and a length-scale L:

$$\kappa_{GM} = \alpha L^2 \overline{ \frac{|f|}{\sqrt{Ri}} }^z$$

where the Eady growth rate has been depth averaged (indicated by the over-line). A local Richardson number is defined Ri = N2/(∂u/∂z)2 which, when combined with thermal wind gives:

$$\frac{1}{Ri} = \frac{(\frac{\partial u}{\partial z})^2}{N^2} = \frac{ ( \frac{g}{f \rho_o} | {\bf \nabla} \sigma | )^2 }{N^2} = \frac{ M^4 }{ |f|^2 N^2 }$$

where M2 is defined $M^2 = \frac{g}{\rho_o} |{\bf \nabla} \sigma|$. Substituting into the formula for κGM gives:

$$\kappa_{GM} = \alpha L^2 \overline{ \frac{M^2}{N} }^z = \alpha L^2 \overline{ \frac{M^2}{N^2} N }^z = \alpha L^2 \overline{ |S| N }^z$$

Tapering and stability

Experience with the GFDL model showed that the GM scheme has to be matched to the convective parameterization. This was originally expressed in connection with the introduction of the KPP boundary layer scheme lar-eta:94 but in fact, as subsequent experience with the MIT model has found, is necessary for any convective parameterization.

Subroutine

S/R GMREDI_SLOPE_LIMIT (pkg/gmredi/gmredi_slope_limit.F)

σx, sx: SlopeX (argument)

σy, sy: SlopeY (argument)

σz: dSigmadRReal (argument)

zσ*: dRdSigmaLtd (argument)

Taper functions used in GKW91 and DM95. Effective slope as a function of 'true' slope using Cox slope clipping, GKW91 limiting and DM95 limiting.

Slope clipping

Deep convection sites and the mixed layer are indicated by homogenized, unstable or nearly unstable stratification. The slopes in such regions can be either infinite, very large with a sign reversal or simply very large. From a numerical point of view, large slopes lead to large variations in the tensor elements (implying large bolus flow) and can be numerically unstable. This was first recognized by cox87 who implemented “slope clipping” in the isopycnal mixing tensor. Here, the slope magnitude is simply restricted by an upper limit:

$$\begin{aligned} \begin{aligned} |\nabla \sigma| & = & \sqrt{ \sigma_x^2 + \sigma_y^2 } \\\ S_{lim} & = & - \frac{|\nabla \sigma|}{ S_{max} } \;\;\;\;\;\;\;\; \mbox{where $S_{max}$ is a parameter} \\\ \sigma_z^\star & = & \min( \sigma_z , S_{lim} ) \\\ {[s_x,s_y]} & = & - \frac{ [\sigma_x,\sigma_y] }{\sigma_z^\star}\end{aligned} \end{aligned}$$

Notice that this algorithm assumes stable stratification through the “min” function. In the case where the fluid is well stratified (σz < Slim) then the slopes evaluate to:

$${[s_x,s_y]} = - \frac{ [\sigma_x,\sigma_y] }{\sigma_z}$$

while in the limited regions (σz > Slim) the slopes become:

$${[s_x,s_y]} = \frac{ [\sigma_x,\sigma_y] }{|\nabla \sigma|/S_{max}}$$

so that the slope magnitude is limited $\sqrt{s_x^2 + s_y^2} = S_{max}$.

The slope clipping scheme is activated in the model by setting GM_taper_scheme = ’clipping’ in data.gmredi.

Even using slope clipping, it is normally the case that the vertical diffusion term (with coefficient $\kappa_\rho{\bf K}_{33} = \kappa_\rho S_{max}^2$) is large and must be time-stepped using an implicit procedure (see section on discretisation and code later). Fig. [fig-mixedlayer] shows the mixed layer depth resulting from a) using the GM scheme with clipping and b) no GM scheme (horizontal diffusion). The classic result of dramatically reduced mixed layers is evident. Indeed, the deep convection sites to just one or two points each and are much shallower than we might prefer. This, it turns out, is due to the over zealous re-stratification due to the bolus transport parameterization. Limiting the slopes also breaks the adiabatic nature of the GM/Redi parameterization, re-introducing diabatic fluxes in regions where the limiting is in effect.

Tapering: Gerdes, Koberle and Willebrand, Clim. Dyn. 1991

The tapering scheme used in gkw:91 addressed two issues with the clipping method: the introduction of large vertical fluxes in addition to convective adjustment fluxes is avoided by tapering the GM/Redi slopes back to zero in low-stratification regions; the adjustment of slopes is replaced by a tapering of the entire GM/Redi tensor. This means the direction of fluxes is unaffected as the amplitude is scaled.

The scheme inserts a tapering function, f1(S), in front of the GM/Redi tensor:

$$f_1(S) = \min \left[ 1, \left( \frac{S_{max}}{|S|}\right)^2 \right]$$

where Smax is the maximum slope you want allowed. Where the slopes, |S| < Smax then f1(S) = 1 and the tensor is un-tapered but where |S| ≥ Smax then f1(S) scales down the tensor so that the effective vertical diffusivity term κf1(S)|S|2 = κSmax2.

The GKW91 tapering scheme is activated in the model by setting GM_taper_scheme = ’gkw91’ in data.gmredi.

Tapering: Danabasoglu and McWilliams, J. Clim. 1995

The tapering scheme used by followed a similar procedure but used a different tapering function, f1(S):

$$f_1(S) = \frac{1}{2} \left( 1+\tanh \left[ \frac{S_c - |S|}{S_d} \right] \right)$$

where Sc = 0.004 is a cut-off slope and Sd = 0.001 is a scale over which the slopes are smoothly tapered. Functionally, the operates in the same way as the GKW91 scheme but has a substantially lower cut-off, turning off the GM/Redi SGS parameterization for weaker slopes.

The DM95 tapering scheme is activated in the model by setting GM_taper_scheme = ’dm95’ in data.gmredi.

Tapering: Large, Danabasoglu and Doney, JPO 1997

The tapering used in lar-eta:97 is based on the DM95 tapering scheme, but also tapers the scheme with an additional function of height, f2(z), so that the GM/Redi SGS fluxes are reduced near the surface:

$$f_2(z) = \frac{1}{2} \left( 1 + \sin(\pi \frac{z}{D} - \frac{\pi}{2})\right)$$

where D = Lρ|S| is a depth-scale and Lρ = c/f with c = 2 m s − 1. This tapering with height was introduced to fix some spurious interaction with the mixed-layer KPP parameterization.

The LDD97 tapering scheme is activated in the model by setting GM_taper_scheme = ’ldd97’ in data.gmredi.

Figure missing Mixed layer depth using GM parameterization with a) Cox slope clipping and for comparison b) using horizontal constant diffusion.

Figure missing Mixed layer depth using GM parameterization with a) Cox slope clipping and for comparison b) using horizontal constant diffusion.

Package Reference

------------------------------------------------------------------------
<-Name->|Levs|<-parsing code->|<--  Units   -->|<- Tile (max=80c) 
------------------------------------------------------------------------
GM_VisbK|  1 |SM P    M1      |m^2/s           |Mixing coefficient from Visbeck etal parameterization
GM_Kux  | 15 |UU P 177MR      |m^2/s           |K_11 element (U.point, X.dir) of GM-Redi tensor
GM_Kvy  | 15 |VV P 176MR      |m^2/s           |K_22 element (V.point, Y.dir) of GM-Redi tensor
GM_Kuz  | 15 |UU   179MR      |m^2/s           |K_13 element (U.point, Z.dir) of GM-Redi tensor
GM_Kvz  | 15 |VV   178MR      |m^2/s           |K_23 element (V.point, Z.dir) of GM-Redi tensor
GM_Kwx  | 15 |UM   181LR      |m^2/s           |K_31 element (W.point, X.dir) of GM-Redi tensor
GM_Kwy  | 15 |VM   180LR      |m^2/s           |K_32 element (W.point, Y.dir) of GM-Redi tensor
GM_Kwz  | 15 |WM P    LR      |m^2/s           |K_33 element (W.point, Z.dir) of GM-Redi tensor
GM_PsiX | 15 |UU   184LR      |m^2/s           |GM Bolus transport stream-function : X component
GM_PsiY | 15 |VV   183LR      |m^2/s           |GM Bolus transport stream-function : Y component
GM_KuzTz| 15 |UU   186MR      |degC.m^3/s      |Redi Off-diagonal Tempetature flux: X component
GM_KvzTz| 15 |VV   185MR      |degC.m^3/s      |Redi Off-diagonal Tempetature flux: Y component

Experiments and tutorials that use gmredi

  • Global Ocean tutorial, in tutorial_global_oce_latlon verification directory, described in section [sec:eg-global]
  • Front Relax experiment, in front_relax verification directory.
  • Ideal 2D Ocean experiment, in ideal_2D_oce verification directory.