Back to home page

darwin3

 
 

    


Warning, /doc/phys_pkgs/gmredi.rst is written in an unsupported language. File is not indexed.

view on githubraw file Latest commit c54b782c on 2024-04-24 14:55:28 UTC
8679f9097b Jeff*0001 .. _sub_phys_pkg_gmredi:
                0002 
7337b40cfd Navi*0003 GMREDI: Gent-McWilliams/Redi Eddy Parameterization
                0004 **************************************************
8679f9097b Jeff*0005 
7337b40cfd Navi*0006 Introduction
                0007 ============
8679f9097b Jeff*0008 
7337b40cfd Navi*0009 Package :filelink:`gmredi <pkg/gmredi>` parameterizes
                0010 the effects of unresolved mesoscale eddies on tracer distribution,
                0011 i.e., temperature, salinity and other tracers.
                0012 
                0013 There are two parts to the Redi/GM subgrid-scale parameterization of geostrophic
                0014 eddies. The first, the Redi scheme (Redi 1982 :cite:`redi1982`), aims to mix tracer properties along
8679f9097b Jeff*0015 isentropes (neutral surfaces) by means of a diffusion operator oriented
a4576c7cde Juli*0016 along the local isentropic surface. The second part, GM
7337b40cfd Navi*0017 (Gent and McWiliams 1990 :cite:`gen-mcw:90`, Gent et al. 1995 :cite:`gen-eta:95`), adiabatically
                0018 rearranges tracers through an advective flux where the advecting flow
8679f9097b Jeff*0019 is a function of slope of the isentropic surfaces.
                0020 
7337b40cfd Navi*0021 A description of key package equations used is
                0022 given in :numref:`ssub_phys_pkg_gmredi_descr`.
                0023 CPP options enable or disable different aspects of the package
                0024 (:numref:`ssub_phys_pkg_gmredi_config`). Run-time options, flags, and filenames
                0025 are set in ``data.gmredi``
                0026 (:numref:`ssub_phys_pkg_gmredi_runtime`). Available diagnostics
                0027 output is listed in :numref:`ssub_phys_pkg_gmredi_diagnostics`.
                0028 
                0029 .. _ssub_phys_pkg_gmredi_descr:
                0030 
                0031 Description
                0032 ===========
                0033 
                0034 The first GCM implementation of the Redi scheme was by Cox (1987) :cite:`cox87` in the GFDL ocean
8679f9097b Jeff*0035 circulation model. The original approach failed to distinguish between
                0036 isopycnals and surfaces of locally referenced potential density (now
7337b40cfd Navi*0037 called neutral surfaces), which are proper isentropes for the ocean. As
8679f9097b Jeff*0038 will be discussed later, it also appears that the Cox implementation is
                0039 susceptible to a computational mode. Due to this mode, the Cox scheme
                0040 requires a background lateral diffusion to be present to conserve the
                0041 integrity of the model fields.
                0042 
                0043 The GM parameterization was then added to the GFDL code in the form of a
7337b40cfd Navi*0044 non-divergent bolus velocity. The method defines two streamfunctions
8679f9097b Jeff*0045 expressed in terms of the isoneutral slopes subject to the boundary
                0046 condition of zero value on upper and lower boundaries. The horizontal
                0047 bolus velocities are then the vertical derivative of these functions.
7337b40cfd Navi*0048 Here in lies a problem highlighted by Griffies et al. (1998) :cite:`gretal:98`: the
                0049 bolus velocities involve multiple derivatives on the potential density field,
                0050 which can consequently give rise to noise. Griffies et al. point out that the GM
8679f9097b Jeff*0051 bolus fluxes can be identically written as a skew flux which involves
                0052 fewer differential operators. Further, combining the skew flux
                0053 formulation and Redi scheme, substantial cancellations take place to the
                0054 point that the horizontal fluxes are unmodified from the lateral
                0055 diffusion parameterization.
                0056 
b938a3c63b antn*0057 .. _GM_redi_desc:
                0058 
8679f9097b Jeff*0059 Redi scheme: Isopycnal diffusion
7337b40cfd Navi*0060 --------------------------------
8679f9097b Jeff*0061 
                0062 The Redi scheme diffuses tracers along isopycnals and introduces a term
                0063 in the tendency (rhs) of such a tracer (here :math:`\tau`) of the form:
                0064 
7337b40cfd Navi*0065 .. math:: \nabla \cdot ( \kappa_\rho {\bf K}_{\rm Redi} \nabla \tau )
8679f9097b Jeff*0066 
                0067 where :math:`\kappa_\rho` is the along isopycnal diffusivity and
7337b40cfd Navi*0068 :math:`{\bf K}_{\rm Redi}` is a rank 2 tensor that projects the gradient of
8679f9097b Jeff*0069 :math:`\tau` onto the isopycnal surface. The unapproximated projection
                0070 tensor is:
                0071 
                0072 .. math::
                0073 
a4576c7cde Juli*0074    {\bf K}_{\rm Redi} = \frac{1}{1 + |{\bf S}|^2}
7337b40cfd Navi*0075    \begin{pmatrix}
8679f9097b Jeff*0076    1 + S_y^2& -S_x S_y & S_x \\
                0077    -S_x S_y  & 1 + S_x^2 & S_y \\
7337b40cfd Navi*0078    S_x & S_y & |{\bf S}|^2 \\
                0079    \end{pmatrix}
8679f9097b Jeff*0080 
5b172de0d2 Jean*0081 Here, :math:`S_x = \partial_x \sigma / (- \partial_z \sigma)`,
                0082 :math:`S_y = \partial_y \sigma / (- \partial_z \sigma)` are the components
                0083 of the isoneutral slope, and :math:`|{\bf S}|^2 = S_x^2 + S_y^2`.
8679f9097b Jeff*0084 
                0085 The first point to note is that a typical slope in the ocean interior is
                0086 small, say of the order :math:`10^{-4}`. A maximum slope might be of
                0087 order :math:`10^{-2}` and only exceeds such in unstratified regions
7337b40cfd Navi*0088 where the slope is ill-defined. It is, therefore, justifiable, and
                0089 customary, to make the small-slope approximation, i.e., :math:`|{\bf S}| \ll 1`. Then
                0090 Redi projection tensor then simplifies to:
8679f9097b Jeff*0091 
                0092 .. math::
7337b40cfd Navi*0093    {\bf K}_{\rm Redi} =
                0094    \begin{pmatrix}
8679f9097b Jeff*0095    1 & 0 & S_x \\
                0096    0 & 1 & S_y \\
7337b40cfd Navi*0097    S_x & S_y & |{\bf S}|^2 \\
                0098    \end{pmatrix}
8679f9097b Jeff*0099 
a4576c7cde Juli*0100 .. _GM_bolus_desc:
9c8516d9da Jeff*0101 
8679f9097b Jeff*0102 GM parameterization
7337b40cfd Navi*0103 -------------------
8679f9097b Jeff*0104 
7337b40cfd Navi*

** Warning **

Wide character in print at /usr/local/share/lxr/source line 1030, <$git> line 106.

0105 The GM parameterization aims to represent the advective or “transport” 8679f9097b Jeff*

** Warning **

Wide character in print at /usr/local/share/lxr/source line 1030, <$git> line 107.

0106 effect of geostrophic eddies by means of a “bolus” velocity, 7337b40cfd Navi*0107 :math:`{\bf u}^\star`. The divergence of this advective flux is added to 8679f9097b Jeff*0108 the tracer tendency equation (on the rhs): 0109 7337b40cfd Navi*0110 .. math:: - \nabla \cdot ( \tau {\bf u}^\star ) 8679f9097b Jeff*0111 7337b40cfd Navi*0112 The bolus velocity :math:`{\bf u}^\star` is defined as the rotational part 0113 of a streamfunction 0114 :math:`{\bf F}^\star = (F_x^\star, F_y^\star, 0)`: 8679f9097b Jeff*0115 0116 .. math:: 0117 7337b40cfd Navi*0118 {\bf u}^\star = \nabla \times {\bf F}^\star = 0119 \begin{pmatrix} 8679f9097b Jeff*0120 - \partial_z F_y^\star \\ 0121 + \partial_z F_x^\star \\ 0122 \partial_x F_y^\star - \partial_y F_x^\star 7337b40cfd Navi*0123 \end{pmatrix} 8679f9097b Jeff*0124 0125 and thus is automatically non-divergent. In the GM parameterization, the 0126 streamfunction is specified in terms of the isoneutral slopes 0127 :math:`S_x` and :math:`S_y`: 0128 0129 .. math:: 0130 0131 \begin{aligned} 7337b40cfd Navi*0132 F_x^\star & = -\kappa_{\rm GM} S_y\\ 0133 F_y^\star & = \kappa_{\rm GM} S_x 0134 \end{aligned} 8679f9097b Jeff*0135 0136 with boundary conditions :math:`F_x^\star=F_y^\star=0` on upper and 7337b40cfd Navi*0137 lower boundaries. :math:`\kappa_{\rm GM}` is colloquially called the isopycnal "thickness diffusivity" 0138 or the "GM diffusivity". The bolus transport in the GM 0139 parameterization is thus given by: 8679f9097b Jeff*0140 0141 .. math:: 0142 7337b40cfd Navi*0143 {\bf u}^\star = 0144 \begin{pmatrix} 8679f9097b Jeff*0145 u^\star \\ 0146 v^\star \\ 0147 w^\star 7337b40cfd Navi*0148 \end{pmatrix} = 0149 \begin{pmatrix} 0bad585a21 Navi*0150 - \partial_z (\kappa_{\rm GM} S_x) \\ 0151 - \partial_z (\kappa_{\rm GM} S_y) \\ 7337b40cfd Navi*0152 \partial_x (\kappa_{\rm GM} S_x) + \partial_y (\kappa_{\rm GM} S_y) 0153 \end{pmatrix} 8679f9097b Jeff*0154 7337b40cfd Navi*0155 This is the "advective form" of the GM parameterization as applied by Danabasoglu and McWilliams (1995) :cite:`danabasoglu:95`, 0156 employed in the GFDL Modular Ocean Model (MOM) versions 1 and 2. To use the advective form in MITgcm, set 0157 :varlink:`GM_AdvForm` ``=.TRUE.`` in ``data.gmredi`` 0158 (also requires ``#define`` :varlink:`GM_BOLUS_ADVEC` and :varlink:`GM_EXTRA_DIAGONAL`). 0159 As implemented in the MITgcm code, :math:`{\bf u}^\star` is simply added to Eulerian :math:`\vec{\bf u}` 0160 (i.e., MITgcm variables :varlink:`uVel`, :varlink:`vVel`, :varlink:`wVel`) 0161 and passed to tracer advection subroutines (:numref:`advection_schemes`) 0162 unless :varlink:`GM_AdvSeparate` ``=.TRUE.`` in ``data.gmredi``, in which case the bolus transport is computed separately. 8679f9097b Jeff*0163 7337b40cfd Navi*0164 Note that in MITgcm, the variables for the GM bolus 0165 streamfunction :varlink:`GM_PsiX` and :varlink:`GM_PsiY` are defined: 8679f9097b Jeff*0166 0167 .. math:: 5b172de0d2 Jean*0168 \begin{aligned} 7337b40cfd Navi*0169 \begin{pmatrix} 0170 \sf{GM\_PsiX} \\ 0171 \sf{GM\_PsiY} 0172 \end{pmatrix} = 0173 \begin{pmatrix} 0bad585a21 Navi*0174 \kappa_{\rm GM} S_x \\ 0175 \kappa_{\rm GM} S_y 7337b40cfd Navi*0176 \end{pmatrix} = 0177 \begin{pmatrix} 8679f9097b Jeff*0178 F_y^\star \\ 0179 -F_x^\star 7337b40cfd Navi*0180 \end{pmatrix} 5b172de0d2 Jean*0181 \end{aligned} 0182 :label: GM_bolus_psi 8679f9097b Jeff*0183 9c8516d9da Jeff*0184 .. _sub_gmredi_skewflux: 0185 8679f9097b Jeff*0186 Griffies Skew Flux 7337b40cfd Navi*0187 ------------------ 8679f9097b Jeff*0188 7337b40cfd Navi*0189 Griffies (1998) :cite:`gr:98` notes that the discretization of bolus velocities involves multiple 8679f9097b Jeff*0190 layers of differencing and interpolation that potentially lead to noisy 0191 fields and computational modes. He pointed out that the bolus flux can 0192 be re-written in terms of a non-divergent flux and a skew-flux: 0193 0194 .. math:: 0195 0196 \begin{aligned} 7337b40cfd Navi*0197 {\bf u}^\star \tau a4576c7cde Juli*0198 & = 7337b40cfd Navi*0199 \begin{pmatrix} 0bad585a21 Navi*0200 - \partial_z ( \kappa_{\rm GM} S_x ) \tau \\ 0201 - \partial_z ( \kappa_{\rm GM} S_y ) \tau \\ 7337b40cfd Navi*0202 \Big[ \partial_x (\kappa_{\rm GM} S_x) + \partial_y (\kappa_{\rm GM} S_y) \Big] \tau 0203 \end{pmatrix} 8679f9097b Jeff*0204 \\ a4576c7cde Juli*0205 & = 7337b40cfd Navi*0206 \begin{pmatrix} 0bad585a21 Navi*0207 - \partial_z ( \kappa_{\rm GM} S_x \tau) \\ 0208 - \partial_z ( \kappa_{\rm GM} S_y \tau) \\ 0209 \partial_x ( \kappa_{\rm GM} S_x \tau) + \partial_y ( \kappa_{\rm GM} S_y \tau) 7337b40cfd Navi*0210 \end{pmatrix} 0211 + \kappa_{\rm GM} \begin{pmatrix} 0212 S_x \partial_z \tau \\ 0213 S_y \partial_z \tau \\ 0214 - S_x \partial_x \tau - S_y \partial_y \tau 0215 \end{pmatrix} 0216 \end{aligned} 8679f9097b Jeff*0217 0218 The first vector is non-divergent and thus has no effect on the tracer 0219 field and can be dropped. The remaining flux can be written: 0220 0bad585a21 Navi*0221 .. math:: \bf{u}^\star \tau = - \kappa_{\rm GM} \bf{K}_{\rm GM} \bf{\nabla} \tau 8679f9097b Jeff*0222 0223 where 0224 0225 .. math:: 0226 7337b40cfd Navi*0227 {\bf K}_{\rm GM} = 0228 \begin{pmatrix} 0229 0 & 0 & -S_x \\ 0230 0 & 0 & -S_y \\ 0231 S_x & S_y & 0 0232 \end{pmatrix} 8679f9097b Jeff*0233 0234 is an anti-symmetric tensor. 0235 0236 This formulation of the GM parameterization involves fewer derivatives 0237 than the original and also involves only terms that already appear in 0238 the Redi mixing scheme. Indeed, a somewhat fortunate cancellation 0239 becomes apparent when we use the GM parameterization in conjunction with 0240 the Redi isoneutral mixing scheme: 0241 0242 .. math:: 0243 7337b40cfd Navi*0244 \kappa_\rho {\bf K}_{\rm Redi} \nabla \tau a4576c7cde Juli*0245 - {\bf u}^\star \tau = 7337b40cfd Navi*0246 ( \kappa_\rho {\bf K}_{\rm Redi} + \kappa_{\rm GM} {\bf K}_{\rm GM} ) \nabla \tau 8679f9097b Jeff*0247 7337b40cfd Navi*0248 If the Redi and GM diffusivities are equal, :math:`\kappa_{\rm GM} = \kappa_{\rho}`, then 8679f9097b Jeff*0249 0250 .. math:: 7337b40cfd Navi*0251 \kappa_\rho {\bf K}_{\rm Redi} + \kappa_{\rm GM} {\bf K}_{\rm GM} = 8679f9097b Jeff*0252 \kappa_\rho 7337b40cfd Navi*0253 \begin{pmatrix} 8679f9097b Jeff*0254 1 & 0 & 0 \\ 0255 0 & 1 & 0 \\ a4576c7cde Juli*0256 2 S_x & 2 S_y & |{\bf S}|^2 7337b40cfd Navi*0257 \end{pmatrix} 8679f9097b Jeff*0258 7337b40cfd Navi*0259 which only differs from the variable Laplacian diffusion tensor by the two 8679f9097b Jeff*0260 non-zero elements in the :math:`z`-row. 0261 0262 .. admonition:: Subroutine 0263 :class: note 0264 7337b40cfd Navi*0265 S/R GMREDI_CALC_TENSOR (:filelink:`pkg/gmredi/gmredi_calc_tensor.F`) 8679f9097b Jeff*0266 5b172de0d2 Jean*0267 :math:`\sigma_x`: **sigmaX** (argument on entry) 0268 0269 :math:`\sigma_y`: **sigmaY** (argument on entry) 0270 0271 :math:`\sigma_z`: **sigmaR** (argument on entry) 0272 0273 .. _sub_gmredi_in_P: 8679f9097b Jeff*0274 5b172de0d2 Jean*0275 Redi and GM schemes in pressure coordinate 0276 ------------------------------------------ 8679f9097b Jeff*0277 5b172de0d2 Jean*0278 When using pressure as a vertical coordinate (see :numref:`isomorphic-equations`), 0279 the Redi scheme can be applied in the same way as in :math:`z-`\coordinates, 0280 considering the slope of isoneutral surface relative to the model isobaric 0281 surface, to rotate the diffusion operator along isoneutral surfaces. 0282 0283 The two components of the slope relative to :math:`p-`\ coordinates are 0284 :math:`S_x^p = \partial_x \sigma / \partial_p \sigma`, 0285 :math:`S_y^p = \partial_y \sigma / \partial_p \sigma`. 0286 Note that for convienience and consistency with the current implementation, 0287 the sign of the slope in :math:`z-` or :math:`p-` coordinates is kept unchanged, i.e., 0288 identical to the sign of horizontal density gradient. 0289 The negative sign is added back in the Redi tensor expression: 0290 0291 .. math:: 0292 {\bf K}_{\rm Redi}^p = 0293 \begin{pmatrix} 0294 1 & 0 & -S_x^p \\ 0295 0 & 1 & -S_y^p \\ 0296 -S_x^p & -S_y^p & |{\bf S^p}|^2 \\ 0297 \end{pmatrix} 8679f9097b Jeff*0298 5b172de0d2 Jean*0299 In contrast, the GM scheme should instead consider the slope of isoneutral 0300 surface relative to geopotential surface (constant :math:`z-`\ surface), so that its 0301 effect will decrease the available potential energy and flatten the isopycnal. 0302 Since :math:`dp = -(\rho g) dz`, the slope to consider would be, in two dimensions: 8679f9097b Jeff*0303 5b172de0d2 Jean*0304 .. math:: S_x^z = \partial_x \sigma / (- \partial_z \sigma) = \frac{1}{\rho g} S_x^p 0305 0306 The effect on tracer :math:`\tau` from the bolus transport (:math:`{\bf u}^\star`) 0307 advection would be: 0308 0309 .. math:: [ {\bf u}^\star \cdot \nabla \tau ]^z 0310 & = u^\star \partial_x \tau + w^\star \partial_z \tau 0311 \\ 0312 & = \rho g \partial_p (\kappa_{\rm GM} \frac{1}{\rho g} S_x^p) \partial_x \tau 0313 - \partial_x (\kappa_{\rm GM} \frac{1}{\rho g} S_x^p) (\rho g) \partial_p \tau 0314 \\ 0315 & = {\bf u}^{\star p} \cdot \nabla^p \tau 0316 0317 .. math:: 0318 {\rm with:}~~ 0319 {\bf u}^{\star p} = 0320 \begin{pmatrix} 0321 u^{\star p} \\ 0322 v^{\star p} \\ 0323 \omega^{\star p} 0324 \end{pmatrix} = 0325 \begin{pmatrix} 0326 \rho g \partial_p (\kappa_{\rm GM} \frac{1}{\rho g} S_x^p) \\ 0327 \rho g \partial_p (\kappa_{\rm GM} \frac{1}{\rho g} S_y^p) \\ 0328 - \rho g \partial_x (\kappa_{\rm GM} \frac{1}{\rho g} S_x^p) 0329 - \rho g \partial_y (\kappa_{\rm GM} \frac{1}{\rho g} S_y^p) 0330 \end{pmatrix} 0331 .. :label: GM_bolus_in_p (note to me: this is commented out) 0332 0333 This formulation above has not been implemented yet and instead only a simplified 0334 version is available that ignores the difference between isobaric surfaces and 0335 horizontal surfaces, as if in the above expression, the density :math:`\rho` were 0336 uniform. This approximation seems valid for the ocean where the 0337 isopycnal slope is generally much larger than the isobaric slope relative to 0338 horizontal surface. 0339 0340 With this approxmation, the expression of the bolus transport simplifies 0341 and becomes "isomorphic" to the :math:`z-`\ coordinate form :eq:`GM_bolus_psi`, 0342 with the sign reversal of all three components of the bolus transport, 0343 due to the expression of the curl in a left-handed coordinate system. 8679f9097b Jeff*0344 a4576c7cde Juli*0345 .. _sub_gmredi_visbeck: 0346 7337b40cfd Navi*0347 Visbeck et al. 1997 GM diffusivity :math:`\kappa_{GM}(x,y)` 0348 ----------------------------------------------------------- 8679f9097b Jeff*0349 7337b40cfd Navi*0350 Visbeck et al. (1997) :cite:`visbeck:97` suggest making the eddy coefficient, 0351 :math:`\kappa_{\rm GM}`, a function of 0352 the Eady growth rate, :math:`|f|/\sqrt{\rm Ri}`. The formula involves a 8679f9097b Jeff*0353 non-dimensional constant, :math:`\alpha`, and a length-scale :math:`L`: 0354 7337b40cfd Navi*0355 .. math:: \kappa_{\rm GM} = \alpha L^2 \overline{ \frac{|f|}{\sqrt{\rm Ri}} }^z 8679f9097b Jeff*0356 0357 where the Eady growth rate has been depth averaged (indicated by the 0358 over-line). A local Richardson number is defined 7337b40cfd Navi*0359 :math:`{\rm Ri} = N^2 / (\partial_z u)^2` which, when combined with thermal wind gives: 8679f9097b Jeff*0360 0361 .. math:: 0362 7337b40cfd Navi*0363 \frac{1}{\rm Ri} = \frac{(\partial u/\partial z)^2}{N^2} = 0364 \frac{ \left ( \dfrac{g}{f \rho_0} | \nabla \sigma | \right )^2 }{N^2} = 8679f9097b Jeff*0365 \frac{ M^4 }{ |f|^2 N^2 } 0366 7337b40cfd Navi*0367 where :math:`M^2 = g | \nabla \sigma| / \rho_0`. Substituting into 0bad585a21 Navi*0368 the formula for :math:`\kappa_{\rm GM}` gives: 8679f9097b Jeff*0369 0370 .. math:: 0371 0bad585a21 Navi*0372 \kappa_{\rm GM} = \alpha L^2 \overline{ \frac{M^2}{N} }^z = 8679f9097b Jeff*0373 \alpha L^2 \overline{ \frac{M^2}{N^2} N }^z = 7337b40cfd Navi*0374 \alpha L^2 \overline{ |{\bf S}| N }^z 8679f9097b Jeff*0375 a4576c7cde Juli*0376 Marshall et al. 2012 GM diffusivity :math:`\kappa_{GM}(x,y)` 0377 ------------------------------------------------------------ 0378 0379 Marshall et al. (2012) :cite:`marshall:12b` via the GEOMETRIC framework suggest 0380 that the eddy coefficient :math:`\kappa_{\rm GM}` be a linear function of a 0381 parameterized total eddy energy :math:`E`, a non-dimensional constant 0382 :math:`\alpha` (namelist parameter :varlink:`GEOM_alpha`; note 0383 bounded in magnitude by 1 in the QG setting) and the 0384 model stratification. A key impact of GEOMETRIC is to make the sensitivity of 0385 the Antarctic Circumpolar Current (ACC) and AMOC to changes in the Southern 0386 Ocean wind forcing closer to those in analogous high resolution models (e.g., 0387 Mak et al. 2018 :cite:`mak:18`). Following Mak et al. (2022) :cite:`mak:22`, the 0388 depth-integrated coefficient is given by 0389 0390 .. math:: 0391 \kappa_{\rm GM} = \alpha \frac{\int E\; \mathrm{d}z}{\int \Gamma (M^2 / N)\; \mathrm{d}z} \Gamma(z) = \alpha \frac{\int E\; \mathrm{d}z}{\int \Gamma |{\bf S}| N\; \mathrm{d}z} \Gamma(z), 0392 0393 where notation is as :numref:`sub_gmredi_visbeck`, noting that :math:`|{\bf S}| 0394 = M^2 / N^2` is the magnitude of the isopycnal slopes. :math:`\Gamma(z)` is some 0395 vertical structure function that can be prescribed if required; default is 0396 :math:`\Gamma(z)\equiv1`, with an option for :math:`\Gamma(z)\equiv N^2 / 0397 N^2_{\rm ref}` following Ferreria et al. (2005) :cite:`ferriera:05`, 0398 activated by setting parameter :varlink:`GEOM_vert_struc` ``= .TRUE.``. The 0399 depth-integrated eddy energy follows its own prognostic equation as 0400 0401 .. math:: 0402 \frac{\mathrm{d}}{\mathrm{d} t} \int E\; \mathrm{d}z + 0403 \nabla \cdot \left( (\tilde{\boldsymbol{u}} - |c|\boldsymbol{e}_x ) \int E\; \mathrm{d}z \right) = 0404 \int \kappa_{\rm gm} |{\bf S}|^2 N^2\; \mathrm{d}z - 0405 \lambda \int E\; \mathrm{d}z + 0406 \nu_E \nabla^2 \int E\; \mathrm{d}z 0407 0408 Terms above, respectively, are the time-evolution, advection, source, dissipation and diffusion of 0409 depth-integrated parameterized total eddy energy, with 0410 :math:`\tilde{\boldsymbol{u}}` the depth-averaged velocity, :math:`c` a 0411 long Rossby phase speed (computed from a WKB approximation, rotated accordingly if ``usingCurvilinearGrid = .TRUE.``), :math:`\lambda` a 0412 dissipation rate (namelist parameter :varlink:`GEOM_lmbda`), and :math:`\nu_E` an energy 0413 diffusion coefficient (namelist parameter :varlink:`GEOM_diffKh_EKE`). Note that bit-for-bit 0414 restart reproducibility GEOM requires an additional pickup file (``pickup_gmredi``). 0415 0416 .. admonition:: note 0417 :class: note 0418 0419 The present GEOMETRIC implementation is strictly for the GM coefficient and 0420 overrides any other specifications in ``data.gmredi`` for the GM coefficient. 0421 The choice of Redi coefficient, clipping/tapering options and others however 0422 are still controlled by the analogous entries in ``data.gmredi``. 0423 9c8516d9da Jeff*0424 .. _sub_gmredi_tapering_stability: 0425 8679f9097b Jeff*0426 Tapering and stability 7337b40cfd Navi*0427 ---------------------- 8679f9097b Jeff*0428 0429 Experience with the GFDL model showed that the GM scheme has to be 0430 matched to the convective parameterization. This was originally 0431 expressed in connection with the introduction of the KPP boundary layer 7337b40cfd Navi*0432 scheme (Large et al. 1994 :cite:`lar-eta:94`) but in fact, as subsequent experience with the MIT model has 8679f9097b Jeff*0433 found, is necessary for any convective parameterization. 0434 0435 Slope clipping 0436 ++++++++++++++ 0437 0438 Deep convection sites and the mixed layer are indicated by homogenized, 0439 unstable or nearly unstable stratification. The slopes in such regions 0440 can be either infinite, very large with a sign reversal or simply very 0441 large. From a numerical point of view, large slopes lead to large 0442 variations in the tensor elements (implying large bolus flow) and can be 7337b40cfd Navi*0443 numerically unstable. This was first recognized by Cox (1987) :cite:`cox87` who implemented 8679f9097b Jeff*

** Warning **

Wide character in print at /usr/local/share/lxr/source line 1030, <$git> line 445.

0444 “slope clipping” in the isopycnal mixing tensor. Here, the slope 0445 magnitude is simply restricted by an upper limit: 0446 0447 .. math:: 0448 0449 \begin{aligned} 7337b40cfd Navi*0450 |\nabla_h \sigma| & = \sqrt{ \sigma_x^2 + \sigma_y^2 }\\ a4576c7cde Juli*0451 S_{\rm lim} & = - \frac{|\nabla_h \sigma|}{ S_{\max} }, 7337b40cfd Navi*0452 \quad \mbox{where $S_{\max}>0$ is a parameter} \\ 0453 \sigma_z^\star & = \min( \sigma_z, S_{\rm lim} ) \\ 0454 {[s_x, s_y]} & = - \frac{ [\sigma_x, \sigma_y] }{\sigma_z^\star} 0455 \end{aligned} 8679f9097b Jeff*0456 0457 Notice that this algorithm assumes stable stratification through the 7337b40cfd Navi*

** Warning **

Wide character in print at /usr/local/share/lxr/source line 1030, <$git> line 459.

0458 “min” function. In the case where the fluid is well stratified 0bad585a21 Navi*0459 (:math:`\sigma_z < S_{\rm lim}`) then the slopes evaluate to: 8679f9097b Jeff*0460 7337b40cfd Navi*0461 .. math:: {[s_x, s_y]} = - \frac{ [\sigma_x, \sigma_y] }{\sigma_z} 8679f9097b Jeff*0462 0bad585a21 Navi*0463 while in the limited regions (:math:`\sigma_z > S_{\rm lim}`) the slopes 8679f9097b Jeff*0464 become: 0465 7337b40cfd Navi*0466 .. math:: {[s_x, s_y]} = \frac{ [\sigma_x, \sigma_y] }{|\nabla_h \sigma| / S_{\max}} 8679f9097b Jeff*0467 0468 so that the slope magnitude is limited :math:`\sqrt{s_x^2 + s_y^2} = 7337b40cfd Navi*0469 S_{\max}`. 8679f9097b Jeff*0470 0471 The slope clipping scheme is activated in the model by setting 7337b40cfd Navi*

** Warning **

Wide character in print at /usr/local/share/lxr/source line 1030, <$git> line 473.

0472 :varlink:`GM_taper_scheme` ``= ’clipping’`` in ``data.gmredi``. 8679f9097b Jeff*0473 0474 Even using slope clipping, it is normally the case that the vertical 0475 diffusion term (with coefficient :math:`\kappa_\rho{\bf K}_{33} = 7337b40cfd Navi*0476 \kappa_\rho S_{\max}^2`) is large and must be time-stepped using an 0477 implicit procedure (see :numref:`implicit-backward-stepping`). Fig. 8679f9097b Jeff*0478 [fig-mixedlayer] shows the mixed layer depth resulting from a) using the 0479 GM scheme with clipping and b) no GM scheme (horizontal diffusion). The 0480 classic result of dramatically reduced mixed layers is evident. Indeed, 0481 the deep convection sites to just one or two points each and are much 0482 shallower than we might prefer. This, it turns out, is due to the over 0483 zealous re-stratification due to the bolus transport parameterization. 0484 Limiting the slopes also breaks the adiabatic nature of the GM/Redi 0485 parameterization, re-introducing diabatic fluxes in regions where the 0486 limiting is in effect. 0487 7337b40cfd Navi*0488 .. admonition:: Subroutine 0489 :class: note 0490 0491 S/R GMREDI_SLOPE_LIMIT (:filelink:`pkg/gmredi/gmredi_slope_limit.F`) 8679f9097b Jeff*0492 7337b40cfd Navi*0493 :math:`\sigma_x, s_x`: **SlopeX** (argument) 0494 0495 :math:`\sigma_y, s_y`: **SlopeY** (argument) 0496 0497 :math:`\sigma_z`: **dSigmadRReal** (argument) 8679f9097b Jeff*0498 7337b40cfd Navi*0499 :math:`z_\sigma^{*}`: **dRdSigmaLtd** (argument) 0500 0501 Tapering: Gerdes, Koberle and Willebrand, 1991 (GKW91) 0502 ++++++++++++++++++++++++++++++++++++++++++++++++++++++ 0503 0504 The tapering scheme used in Gerdes et al. (1991) :cite:`gkw:91` (GKW91) 0505 addressed two issues with the clipping 8679f9097b Jeff*0506 method: the introduction of large vertical fluxes in addition to 0507 convective adjustment fluxes is avoided by tapering the GM/Redi slopes 0508 back to zero in low-stratification regions; the adjustment of slopes is 0509 replaced by a tapering of the entire GM/Redi tensor. This means the 0510 direction of fluxes is unaffected as the amplitude is scaled. 0511 0512 The scheme inserts a tapering function, :math:`f_1(S)`, in front of the 0513 GM/Redi tensor: 0514 7337b40cfd Navi*0515 .. math:: f_1(S) = \min \left[ 1, \left( \frac{S_{\max}}{|{\bf S}|}\right)^2 \right] 8679f9097b Jeff*0516 7337b40cfd Navi*0517 where :math:`S_{\max}` is the maximum slope you want allowed. Where the 0518 slopes, :math:`|{\bf S}|<S_{\max}` then :math:`f_1(S) = 1` and the tensor is 0519 un-tapered but where :math:`|{\bf S}| \ge S_{\max}` then :math:`f_1(S)` scales 8679f9097b Jeff*0520 down the tensor so that the effective vertical diffusivity term 7337b40cfd Navi*0521 :math:`\kappa f_1(S) |{\bf S}|^2 = \kappa S_{\max}^2`. 8679f9097b Jeff*0522 0523 The GKW91 tapering scheme is activated in the model by setting 7337b40cfd Navi*

** Warning **

Wide character in print at /usr/local/share/lxr/source line 1030, <$git> line 525.

0524 :varlink:`GM_taper_scheme` ``= ’gkw91’`` in ``data.gmredi``. 0525 0526 .. figure:: figs/tapers.* 0527 :width: 70% 0528 :align: center 0529 :alt: Tapering for GM scheme 0530 :name: tapers 0531 a4576c7cde Juli*0532 Taper functions used in GKW91 and DM95. 7337b40cfd Navi*0533 0534 .. figure:: figs/effective_slopes.* 0535 :width: 70% 0536 :align: center 0537 :alt: Tapering for GM scheme 0538 :name: effective_slopes 8679f9097b Jeff*0539 7337b40cfd Navi*0540 Effective slope as a function of 'true' slope using Cox slope clipping, GKW91 limiting and DM95 limiting. 8679f9097b Jeff*0541 7337b40cfd Navi*0542 Tapering: Danabasoglu and McWilliams, 1995 (DM95) 0543 +++++++++++++++++++++++++++++++++++++++++++++++++ 8679f9097b Jeff*0544 7337b40cfd Navi*0545 The tapering scheme used by Danabasoglu and McWilliams (1995) :cite:`danabasoglu:95` (DM95) 0546 followed a similar procedure but used a different tapering function, :math:`f_1(S)`: 8679f9097b Jeff*0547 7337b40cfd Navi*0548 .. math:: f_1(S) = \frac{1}{2} \left[ 1+\tanh \left( \frac{S_c - |{\bf S}|}{S_d} \right) \right] 8679f9097b Jeff*0549 0550 where :math:`S_c = 0.004` is a cut-off slope and :math:`S_d=0.001` is a 0551 scale over which the slopes are smoothly tapered. Functionally, the 0552 operates in the same way as the GKW91 scheme but has a substantially 7337b40cfd Navi*0553 lower cut-off, turning off the GM/Redi parameterization for weaker 8679f9097b Jeff*0554 slopes. 0555 0556 The DM95 tapering scheme is activated in the model by setting 7337b40cfd Navi*

** Warning **

Wide character in print at /usr/local/share/lxr/source line 1030, <$git> line 558.

0557 :varlink:`GM_taper_scheme` ``= ’dm95’`` in ``data.gmredi``. 8679f9097b Jeff*0558 7337b40cfd Navi*0559 Tapering: Large, Danabasoglu and Doney, 1997 (LDD97) 0560 ++++++++++++++++++++++++++++++++++++++++++++++++++++ 8679f9097b Jeff*0561 7337b40cfd Navi*0562 The tapering used in Large et al. (1997) :cite:`lar-eta:97` (LDD97) 0563 is based on the DM95 tapering scheme, but also 8679f9097b Jeff*0564 tapers the scheme with an additional function of height, :math:`f_2(z)`, 7337b40cfd Navi*0565 so that the GM/Redi subgrid-scale fluxes are reduced near the surface: 8679f9097b Jeff*0566 7337b40cfd Navi*0567 .. math:: f_2(z) = \frac{1}{2} \left[ 1 + \sin \left(\pi \frac{z}{D} - \frac{\pi}{2} \right) \right] 8679f9097b Jeff*0568 7337b40cfd Navi*0569 where :math:`D = (c / f) |{\bf S}|` is a depth scale, with :math:`f` the 0570 Coriolis parameter and :math:`c=2`m/s (corresponding to the first baroclinic wave speed, so that :math:`c/f` is the Rossby radius). 0571 This tapering that varies with depth 0572 was introduced to fix some spurious interaction with the mixed-layer KPP 8679f9097b Jeff*0573 parameterization. 0574 0575 The LDD97 tapering scheme is activated in the model by setting 7337b40cfd Navi*

** Warning **

Wide character in print at /usr/local/share/lxr/source line 1030, <$git> line 577.

0576 :varlink:`GM_taper_scheme` ``= ’ldd97’`` in ``data.gmredi``. 0577 0578 .. _ssub_phys_pkg_gmredi_config: 0579 0580 GMREDI configuration and compiling 0581 ================================== 0582 0583 Compile-time options 0584 -------------------- 0585 0586 As with all MITgcm packages, GMREDI can be turned on or off at compile time 0587 (see :numref:`building_code`) 0588 0589 - using the ``packages.conf`` file by adding ``gmredi`` to it 0590 0591 - or using :filelink:`genmake2 <tools/genmake2>` adding ``-enable=gmredi`` or 0592 ``-disable=gmredi`` switches 0593 0594 - **required packages and CPP options**: 0595 :filelink:`gmredi <pkg/gmredi>` requires 0596 0597 Parts of the :filelink:`gmredi <pkg/gmredi>` code can be enabled or disabled at 0598 compile time via CPP preprocessor flags. These options are set in 0599 :filelink:`GMREDI_OPTIONS.h <pkg/gmredi/GMREDI_OPTIONS.h>`. 0600 :numref:`tab_phys_pkg_gmredi_cpp` summarizes the most important ones. For additional 0601 options see :filelink:`GMREDI_OPTIONS.h <pkg/gmredi/GMREDI_OPTIONS.h>`. 0602 0603 .. tabularcolumns:: |\Y{.375}|\Y{.1}|\Y{.55}| 0604 0605 .. csv-table:: Some of the most relevant CPP preprocessor flags in the :filelink:`gmredi <pkg/gmredi>` package. 0606 :header: "CPP option", "Default", "Description" 0607 :widths: 30, 10, 60 0608 :name: tab_phys_pkg_gmredi_cpp 0609 0610 :varlink:`GM_NON_UNITY_DIAGONAL`, #define, allows the leading diagonal (top two rows) to be non-unity 0611 :varlink:`GM_EXTRA_DIAGONAL`, #define, allows different values of :math:`\kappa_{\rm GM}` and :math:`\kappa_{\rho}`; also required for advective form 0612 :varlink:`GM_BOLUS_ADVEC`, #define, allows use of the advective form (bolus velocity) 0613 :varlink:`GM_BOLUS_BVP`, #define, allows use of Boundary-Value-Problem method to evaluate bolus transport 0614 :varlink:`ALLOW_GM_LEITH_QG`, #undef, allow QG Leith variable viscosity to be added to GMRedi coefficient 0615 :varlink:`GM_VISBECK_VARIABLE_K`, #undef, allows Visbeck et al. formulation to compute :math:`\kappa_{\rm GM}` a4576c7cde Juli*0616 :varlink:`GM_GEOM_VARIABLE_K`, #undef, allows Marshall et al. formulation to compute :math:`\kappa_{\rm GM}` 7337b40cfd Navi*0617 0618 .. _ssub_phys_pkg_gmredi_runtime: 0619 0620 Run-time parameters 0621 =================== 0622 0623 Run-time parameters (see :numref:`tab_phys_pkg_gmredi_runtimeparms`) are set in 0624 ``data.gmredi`` (read in :filelink:`pkg/gmredi/gmredi_readparms.F`). 0625 0626 Enabling the package 0627 -------------------- 0628 0629 :filelink:`gmredi <pkg/gmredi>` package is switched on/off at run-time by 0630 setting :varlink:`useGMREDI` ``= .TRUE.,`` in ``data.pkg``. 0631 0632 General flags and parameters 0633 ---------------------------- 0634 0635 :numref:`tab_phys_pkg_gmredi_runtimeparms` lists most run-time parameters. 0636 0637 .. tabularcolumns:: |\Y{.275}|\Y{.20}|\Y{.525}| 0638 0639 .. table:: Run-time parameters and default values 0640 :class: longtable 0641 :name: tab_phys_pkg_gmredi_runtimeparms 0642 0643 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0644 | Name | Default value | Description | 0645 +====================================+==============================+=========================================================================+ 0646 | :varlink:`GM_AdvForm` | FALSE | use advective form (bolus velocity); FALSE uses skewflux form | 0647 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0648 | :varlink:`GM_AdvSeparate` | FALSE | do advection by Eulerian and bolus velocity separately | 0649 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0650 | :varlink:`GM_background_K` | 0.0 | thickness diffusivity :math:`\kappa_{\rm GM}` (m\ :sup:`2`\ /s) | 0651 | | | (GM bolus transport) | 0652 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0653 | :varlink:`GM_isopycK` | :varlink:`GM_background_K` | isopycnal diffusivity :math:`\kappa_{\rho}` (m\ :sup:`2`\ /s) | 0654 | | | (Redi tensor) | 0655 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0656 | :varlink:`GM_maxSlope` | 1.0E-02 | maximum slope (tapering/clipping) | 0657 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0658 | :varlink:`GM_Kmin_horiz` | 0.0 | minimum horizontal diffusivity (m\ :sup:`2`\ /s) | 0659 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0660 | :varlink:`GM_Small_Number` | 1.0E-20 | :math:`\epsilon` used in computing the slope | 0661 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0662 | :varlink:`GM_slopeSqCutoff` | 1.0E+48 | :math:`|{\bf S}|^2` cut-off value for zero taper function | 0663 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0664 | :varlink:`GM_taper_scheme` | ' ' | taper scheme option ('orig', 'clipping', 'fm07', 'stableGmAdjTap', | 0665 | | | 'linear', 'ac02', 'gkw91', 'dm95', 'ldd97') | 0666 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0667 | :varlink:`GM_maxTransLay` | 500.0 | maximum transition layer thickness (m) | 0668 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0669 | :varlink:`GM_facTrL2ML` | 5.0 | maximum trans. layer thick. as a factor of local mixed-layer depth | 0670 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0671 | :varlink:`GM_facTrL2dz` | 1.0 | minimum trans. layer thick. as a factor of local dr | 0672 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0673 | :varlink:`GM_Scrit` | 0.004 | :math:`S_c` parameter for 'dm95' and 'ldd97 ' tapering function | 0674 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0675 | :varlink:`GM_Sd` | 0.001 | :math:`S_d` parameter for 'dm95' and 'ldd97' tapering function | 0676 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0677 | :varlink:`GM_UseBVP` | FALSE | use Boundary-Value-Problem method for bolus transport | 0678 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0679 | :varlink:`GM_BVP_ModeNumber` | 1 | vertical mode number used for speed :math:`c` in BVP transport | 0680 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0681 | :varlink:`GM_BVP_cMin` | 1.0E-01 | minimum value for wave speed parameter :math:`c` in BVP (m/s) | 0682 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0683 | :varlink:`GM_UseSubMeso` | FALSE | use sub-mesoscale eddy parameterization for bolus transport | 0684 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0685 | :varlink:`subMeso_Ceff` | 7.0E-02 | efficiency coefficient of mixed-layer eddies | 0686 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0687 | :varlink:`subMeso_invTau` | 2.0E-06 | inverse of mixing timescale in sub-meso parameterization (s\ :sup:`-1`) | 0688 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0689 | :varlink:`subMeso_LfMin` | 1.0E+03 | minimum value for length-scale :math:`L_f` (m) | 0690 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0691 | :varlink:`subMeso_Lmax` | 110.0E+03 | maximum horizontal grid-scale length (m) | 0692 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0693 | :varlink:`GM_Visbeck_alpha` | 0.0 | :math:`\alpha` parameter for Visbeck et al. scheme (non-dim.) | 0694 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0695 | :varlink:`GM_Visbeck_length` | 200.0E+03 | :math:`L` length scale parameter for Visbeck et al. scheme (m) | 0696 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0697 | :varlink:`GM_Visbeck_depth` | 1000.0 | depth (m) over which to average in computing Visbeck | 0698 | | | :math:`\kappa_{\rm GM}` | 0699 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0700 | :varlink:`GM_Visbeck_maxSlope` | :varlink:`GM_maxSlope` | maximum slope used in computing Visbeck et al. :math:`\kappa_{\rm GM}` | 0701 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0702 | :varlink:`GM_Visbeck_minVal_K` | 0.0 | minimum :math:`\kappa_{\rm GM}` (m\ :sup:`2`\ /s) using Visbeck et al. | 0703 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0704 | :varlink:`GM_Visbeck_maxVal_K` | 2500.0 | maximum :math:`\kappa_{\rm GM}` (m\ :sup:`2`\ /s) using Visbeck et al. | 0705 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ c54b782cac rosi*0706 | :varlink:`GM_useGEOM` | FALSE | use GEOM scheme if #define :varlink:`GM_GEOM_VARIABLE_K` | a4576c7cde Juli*0707 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0708 | :varlink:`GEOM_alpha` | 0.06 | :math:`\alpha` parameter for GEOM scheme (non-dim.) | 0709 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0710 | :varlink:`GEOM_lmbda` | 1.16E-07 | eddy energy dissipation rate for GEOM scheme (s\ :sup:`-1`) | 0711 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0712 | :varlink:`GEOM_diffKh_EKE` | 5.0E+02 | eddy energy diffusion coefficient for GEOM scheme (m\ :sup:`2`\ /s) | 0713 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0714 | :varlink:`GEOM_ini_EKE` | 1.0E-03 | initial value for parameterized total depth eddy energy for GEOM scheme | 0715 | | | (m\ :sup:`3`\ / s\ :sup:`2`) | 0716 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0717 | :varlink:`GEOM_vert_struc` | FALSE | use imposed vertical structure function in GEOM scheme | 0718 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0719 | :varlink:`GEOM_vert_struc_min` | 0.1 | minimum value of vertical structure function in GEOM scheme (non-dim.) | 0720 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0721 | :varlink:`GEOM_vert_struc_max` | 1.0 | maximum value of vertical structure function in GEOM scheme (non-dim.) | 0722 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0723 | :varlink:`GEOM_minVal_K` | 0.0 | minimum :math:`\kappa_{\rm GM}` (m\ :sup:`2`\ /s) using GEOM scheme | 0724 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0725 | :varlink:`GEOM_maxVal_K` | 2500.0 | maximum :math:`\kappa_{\rm GM}` (m\ :sup:`2`\ /s) using GEOM scheme | 0726 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 7337b40cfd Navi*0727 | :varlink:`GM_useLeithQG` | FALSE | add Leith QG viscosity to GMRedi tensor | 0728 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0729 | :varlink:`GM_iso2dFile` | ' ' | input file for 2D (:math:`x,y`) scaling of isopycnal diffusivity | 0730 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0731 | :varlink:`GM_iso1dFile` | ' ' | input file for 1D vert. scaling of isopycnal diffusivity | 0732 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0733 | :varlink:`GM_bol2dFile` | ' ' | input file for 2D (:math:`x,y`) scaling of thickness diffusivity | 0734 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0735 | :varlink:`GM_bol1dFile` | ' ' | input file for 1D vert. scaling of thickness diffusivity | 0736 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0737 | :varlink:`GM_background_K3dFile` | ' ' | input file for 3D (:math:`x,y,r`) :varlink:`GM_background_K` | 0738 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0739 | :varlink:`GM_isopycK3dFile` | ' ' | input file for 3D (:math:`x,y,r`) :varlink:`GM_isopycK` | 0740 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 0741 | :varlink:`GM_MNC` | :varlink:`useMNC` | write GMREDI snapshot output using :filelink:`/pkg/mnc` | 0742 +------------------------------------+------------------------------+-------------------------------------------------------------------------+ 8679f9097b Jeff*0743 0744 .. _ssub_phys_pkg_gmredi_diagnostics: 0745 7337b40cfd Navi*0746 GMREDI Diagnostics 0747 ================== 8679f9097b Jeff*0748 0749 :: 0750 7337b40cfd Navi*0751 ---------------------------------------------------------------------- 0752 <-Name->|Levs|<- code ->|<-- Units -->|<- Description 0753 ---------------------------------------------------------------------- 0754 GM_VisbK| 1 |SM P M1|m^2/s |Mixing coefficient from Visbeck etal parameterization 0755 GM_hTrsL| 1 |SM P M1|m |Base depth (>0) of the Transition Layer 0756 GM_baseS| 1 |SM P M1|1 |Slope at the base of the Transition Layer 0757 GM_rLamb| 1 |SM P M1|1/m |Slope vertical gradient at Trans. Layer Base (=recip.Lambda) 0758 SubMesLf| 1 |SM P M1|m |Sub-Meso horiz. Length Scale (Lf) 0759 SubMpsiX| 1 |UU M1|m^2/s |Sub-Meso transp.stream-funct. magnitude (Psi0): U component 0760 SubMpsiY| 1 |VV M1|m^2/s |Sub-Meso transp.stream-funct. magnitude (Psi0): V component 0761 GM_Kux | 18 |UU P MR|m^2/s |K_11 element (U.point, X.dir) of GM-Redi tensor 0762 GM_Kvy | 18 |VV P MR|m^2/s |K_22 element (V.point, Y.dir) of GM-Redi tensor 0763 GM_Kuz | 18 |UU MR|m^2/s |K_13 element (U.point, Z.dir) of GM-Redi tensor 0764 GM_Kvz | 18 |VV MR|m^2/s |K_23 element (V.point, Z.dir) of GM-Redi tensor 0765 GM_Kwx | 18 |UM LR|m^2/s |K_31 element (W.point, X.dir) of GM-Redi tensor 0766 GM_Kwy | 18 |VM LR|m^2/s |K_32 element (W.point, Y.dir) of GM-Redi tensor 0767 GM_Kwz | 18 |WM P LR|m^2/s |K_33 element (W.point, Z.dir) of GM-Redi tensor 0768 GM_PsiX | 18 |UU LR|m^2/s |GM Bolus transport stream-function : U component 0769 GM_PsiY | 18 |VV LR|m^2/s |GM Bolus transport stream-function : V component 0770 GM_KuzTz| 18 |UU MR|degC.m^3/s |Redi Off-diagonal Temperature flux: X component 0771 GM_KvzTz| 18 |VV MR|degC.m^3/s |Redi Off-diagonal Temperature flux: Y component 0772 GM_KwzTz| 18 |WM LR|degC.m^3/s |Redi main-diagonal vertical Temperature flux 0773 GM_ubT | 18 |UUr MR|degC.m^3/s |Zonal Mass-Weight Bolus Transp of Pot Temp 0774 GM_vbT | 18 |VVr MR|degC.m^3/s |Meridional Mass-Weight Bolus Transp of Pot Temp 0775 GM_BVPcW| 1 |SU P M1|m/s |WKB wave speed (at Western edge location) 0776 GM_BVPcS| 1 |SV P M1|m/s |WKB wave speed (at Southern edge location) a4576c7cde Juli*0777 GM_GEOMK| 15 |SMRP MR|m/s^2 |GEOM 3d kgm field 0778 GEOMeE | 1 |SM P M1|m^3/s^2 |GEOM parameterized depth-int eddy energy 0779 GEOMstru| 15 |SMRP MR| |spatial structure function 0780 GEOMEgen| 1 |SM P M1|m^3/s^3 |GEOM eddy energy generation tendency 0781 GEOMEdis| 1 |SM P M1|m^3/s^3 |GEOM eddy energy dissipation tendency 0782 GEOMEadv| 1 |SM P M1|m^3/s^3 |GEOM eddy energy advective tendency 0783 GEOMEwav| 1 |SM P M1|m^3/s^3 |GEOM eddy energy wave advection tendency 0784 GEOMElap| 1 |SM P M1|m^3/s^3 |GEOM eddy energy diffusion tendency 0785 GEOM_c1 | 1 |SM P M1|m/s |first baroclinic wave phase speed 0786 GEOMcros| 1 |SM P M1|m/s |signed long Rossby wave phase speed 7337b40cfd Navi*0787 0788 Experiments and tutorials that use GMREDI 0789 ========================================= 8679f9097b Jeff*0790 0bad585a21 Navi*0791 - Southern Ocean Reentrant Channel Example, in :filelink:`verification/tutorial_reentrant_channel`, 0792 described in :numref:`sec_eg_reentrant_channel` 0793 0794 - Global Ocean Simulation, in :filelink:`verification/tutorial_global_oce_latlon`, 0795 described in :numref:`sec_global_oce_latlon` 8679f9097b Jeff*0796 7337b40cfd Navi*0797 - Front Relax experiment, in :filelink:`verification/front_relax` 8679f9097b Jeff*0798 7337b40cfd Navi*0799 - Ideal 2D Ocean experiment, in :filelink:`verification/ideal_2D_oce`.