|
|
|||
Warning, /doc/algorithm/adv-schemes.rst is written in an unsupported language. File is not indexed.
view on githubraw file Latest commit c6adb04b on 2022-04-22 17:45:23 UTC68aadaa6bd Phob*0001 Linear advection schemes 0002 ======================== 0003 0004 The advection schemes known as centered second order, centered fourth 0005 order, first order upwind and upwind biased third order are known as 0006 linear advection schemes because the coefficient for interpolation of 0007 the advected tracer are linear and a function only of the flow, not the 0008 tracer field it self. We discuss these first since they are most 0009 commonly used in the field and most familiar. 0010 9c8516d9da Jeff*0011 .. _adv_cent_2ord: 0012 68aadaa6bd Phob*0013 Centered second order advection-diffusion 0014 ----------------------------------------- 0015 0016 The basic discretization, centered second order, is the default. It is 0017 designed to be consistent with the continuity equation to facilitate 0018 conservation properties analogous to the continuum. However, centered 0019 second order advection is notoriously noisy and must be used in 0020 conjunction with some finite amount of diffusion to produce a sensible 0021 solution. 0022 0023 The advection operator is discretized: 0024 0025 .. math:: c6adb04bff jpga*0026 {\cal A}_c \Delta r_f h_c G_{\rm adv}^\tau = 68aadaa6bd Phob*0027 \delta_i F_x + \delta_j F_y + \delta_k F_r 0028 :label: cent_2nd_ord 0029 0030 where the area integrated fluxes are given by: 0031 0032 .. math:: 0033 0034 \begin{aligned} 0035 F_x & = & U \overline{ \tau }^i \\ 0036 F_y & = & V \overline{ \tau }^j \\ 0037 F_r & = & W \overline{ \tau }^k\end{aligned} 0038 0039 The quantities :math:`U`, :math:`V` and :math:`W` are volume fluxes. 0040 defined as: 0041 0042 .. math:: 0043 0044 \begin{aligned} 0045 U & = & \Delta y_g \Delta r_f h_w u \\ 0046 V & = & \Delta x_g \Delta r_f h_s v \\ 0047 W & = & {\cal A}_c w\end{aligned} 0048 0049 For non-divergent flow, this discretization can be shown to conserve the 0050 tracer both locally and globally and to globally conserve tracer 0051 variance, :math:`\tau^2`. The proof is given in 0052 Adcroft (1995) :cite:`adcroft:95` and Adcroft et al. (1997) :cite:`adcroft:97` . 0053 0054 .. admonition:: S/R :filelink:`GAD_C2_ADV_X <pkg/generic_advdiff/gad_c2_adv_x.F>` 0055 :class: note 0056 0057 | :math:`F_x` : :varlink:`uT` ( argument ) 0058 | :math:`U` : :varlink:`uTrans` ( argument ) 0059 | :math:`\tau` : :varlink:`tracer` ( argument ) 0060 0061 .. admonition:: S/R :filelink:`GAD_C2_ADV_Y <pkg/generic_advdiff/gad_c2_adv_y.F>` 0062 :class: note 0063 0064 | :math:`F_y` : :varlink:`vT` ( argument ) 0065 | :math:`V` : :varlink:`vTrans` ( argument ) 0066 | :math:`\tau` : :varlink:`tracer` ( argument ) 0067 0068 .. admonition:: S/R :filelink:`GAD_C2_ADV_R <pkg/generic_advdiff/gad_c2_adv_r.F>` 0069 :class: note 0070 0071 | :math:`F_r` : :varlink:`wT` ( argument ) 0072 | :math:`W` : :varlink:`rTrans` ( argument ) 0073 | :math:`\tau` : :varlink:`tracer` ( argument ) 0074 0075 Third order upwind bias advection 0076 --------------------------------- 0077 0078 Upwind biased third order advection offers a relatively good compromise** Warning **
Wide character in print at /usr/local/share/lxr/source line 1030, <$git> line 80.
0079 between accuracy and smoothness. It is not a “positive” scheme meaning 0080 false extrema are permitted but the amplitude of such are significantly 0081 reduced over the centered second order method. 0082 0083 The third order upwind fluxes are discretized: 0084 0085 .. math:: 0086 0087 \begin{aligned} 0088 F_x & = & U \overline{\tau - \frac{1}{6} \delta_{ii} \tau}^i 0089 + \frac{1}{2} |U| \delta_i \frac{1}{6} \delta_{ii} \tau \\ 0090 F_y & = & V \overline{\tau - \frac{1}{6} \delta_{ii} \tau}^j 0091 + \frac{1}{2} |V| \delta_j \frac{1}{6} \delta_{jj} \tau \\ 0092 F_r & = & W \overline{\tau - \frac{1}{6} \delta_{ii} \tau}^k 0093 + \frac{1}{2} |W| \delta_k \frac{1}{6} \delta_{kk} \tau \end{aligned} 0094 0095 At boundaries, :math:`\delta_{\hat{n}} \tau` is set to zero allowing 0096 :math:`\delta_{nn}` to be evaluated. We are currently examine the 0097 accuracy of this boundary condition and the effect on the solution. 0098 0099 .. admonition:: S/R :filelink:`GAD_U3_ADV_X <pkg/generic_advdiff/gad_u3_adv_x.F>` 0100 :class: note 0101 0102 | :math:`F_x` : :varlink:`uT` ( argument ) 0103 | :math:`U` : :varlink:`uTrans` ( argument ) 0104 | :math:`\tau` : :varlink:`tracer` ( argument ) 0105 0106 .. admonition:: S/R :filelink:`GAD_U3_ADV_Y <pkg/generic_advdiff/gad_u3_adv_y.F>` 0107 :class: note 0108 0109 | :math:`F_y` : :varlink:`vT` ( argument ) 0110 | :math:`V` : :varlink:`vTrans` ( argument ) 0111 | :math:`\tau` : :varlink:`tracer` ( argument ) 0112 0113 .. admonition:: S/R :filelink:`GAD_U3_ADV_R <pkg/generic_advdiff/gad_u3_adv_r.F>` 0114 :class: note 0115 0116 | :math:`F_r` : :varlink:`wT` ( argument ) 0117 | :math:`W` : :varlink:`rTrans` ( argument ) 0118 | :math:`\tau` : :varlink:`tracer` ( argument ) 0119 0120 Centered fourth order advection 0121 ------------------------------- 0122 0123 Centered fourth order advection is formally the most accurate scheme we 0124 have implemented and can be used to great effect in high resolution 0125 simulations where dynamical scales are well resolved. However, the scheme 0126 is noisy, like the centered second order method, and so must be used with 0127 some finite amount of diffusion. Bi-harmonic is recommended since it is 0128 more scale selective and less likely to diffuse away the well resolved 0129 gradient the fourth order scheme worked so hard to create. 0130 0131 The centered fourth order fluxes are discretized: 0132 0133 .. math:: 0134 0135 \begin{aligned} 0136 F_x & = & U \overline{\tau - \frac{1}{6} \delta_{ii} \tau}^i \\ 0137 F_y & = & V \overline{\tau - \frac{1}{6} \delta_{ii} \tau}^j \\ 0138 F_r & = & W \overline{\tau - \frac{1}{6} \delta_{ii} \tau}^k\end{aligned} 0139 0140 As for the third order scheme, the best discretization near boundaries 0141 is under investigation but currently :math:`\delta_i \tau=0` on a 0142 boundary. 0143 0144 .. admonition:: S/R :filelink:`GAD_C4_ADV_X <pkg/generic_advdiff/gad_c4_adv_x.F>` 0145 :class: note 0146 0147 | :math:`F_x` : :varlink:`uT` ( argument ) 0148 | :math:`U` : :varlink:`uTrans` ( argument ) 0149 | :math:`\tau` : :varlink:`tracer` ( argument ) 0150 0151 .. admonition:: S/R :filelink:`GAD_C4_ADV_Y <pkg/generic_advdiff/gad_c4_adv_y.F>` 0152 :class: note 0153 0154 | :math:`F_y` : :varlink:`vT` ( argument ) 0155 | :math:`V` : :varlink:`vTrans` ( argument ) 0156 | :math:`\tau` : :varlink:`tracer` ( argument ) 0157 0158 .. admonition:: S/R :filelink:`GAD_C4_ADV_R <pkg/generic_advdiff/gad_c4_adv_r.F>` 0159 :class: note 0160 0161 | :math:`F_r` : :varlink:`wT` ( argument ) 0162 | :math:`W` : :varlink:`rTrans` ( argument ) 0163 | :math:`\tau` : :varlink:`tracer` ( argument ) 0164 0165 First order upwind advection 0166 ---------------------------- 0167 0168 Although the upwind scheme is the underlying scheme for the robust or** Warning **
Wide character in print at /usr/local/share/lxr/source line 1030, <$git> line 170.
0169 non-linear methods given in :numref:`nonlinear_adv`, we haven’t actually implemented this method 0170 for general use. It would be very diffusive and it is unlikely that it 0171 could ever produce more useful results than the positive higher order 0172 schemes. 0173 0174 Upwind bias is introduced into many schemes using the *abs* function and 0175 it allows the first order upwind flux to be written: 0176 0177 .. math:: 0178 0179 \begin{aligned} 0180 F_x & = & U \overline{ \tau }^i - \frac{1}{2} |U| \delta_i \tau \\ 0181 F_y & = & V \overline{ \tau }^j - \frac{1}{2} |V| \delta_j \tau \\ 0182 F_r & = & W \overline{ \tau }^k - \frac{1}{2} |W| \delta_k \tau\end{aligned} 0183 0184 If for some reason the above method is desired, the second order 0185 flux limiter scheme described in :numref:`secondord_FL` reduces to the above scheme if the 0186 limiter is set to zero. 0187 0188 .. _nonlinear_adv: 0189 0190 Non-linear advection schemes 0191 ============================ 0192 0193 Non-linear advection schemes invoke non-linear interpolation and are 0194 widely used in computational fluid dynamics (non-linear does not refer 0195 to the non-linearity of the advection operator). The flux limited 0196 advection schemes belong to the class of finite volume methods which 0197 neatly ties into the spatial discretization of the model. 0198 0199 When employing the flux limited schemes, first order upwind or 0200 direct-space-time method, the time-stepping is switched to forward in 0201 time. 0202 0203 .. _secondord_FL: 0204 0205 Second order flux limiters 0206 -------------------------- 0207 0208 The second order flux limiter method can be cast in several ways but is 0209 generally expressed in terms of other flux approximations. For example, 0210 in terms of a first order upwind flux and second order Lax-Wendroff 0211 flux, the limited flux is given as: 0212 c6adb04bff jpga*0213 .. math:: F = (1 - \psi(r)) F_1 + \psi(r) F_{\rm LW} 68aadaa6bd Phob*0214 :label: limited_flux 0215 0216 where :math:`\psi(r)` is the limiter function, 0217 0218 .. math:: F_1 = u \overline{\tau}^i - \frac{1}{2} |u| \delta_i \tau 0219 0220 is the upwind flux, 0221 c6adb04bff jpga*0222 .. math:: F_{\rm LW} = u \overline{\tau}^i - \frac{1}{2} c |u| \delta_i \tau 68aadaa6bd Phob*0223 c6adb04bff jpga*0224 is the Lax-Wendroff flux and :math:`c = \frac{|u| \Delta t}{\Delta x}` is 68aadaa6bd Phob*0225 the Courant (CFL) number. 0226 0227 The limiter function, :math:`\psi(r)`, takes the slope ratio 0228 0229 .. math:: 0230 0231 \begin{aligned} 0232 r = \frac{ \tau_{i-1} - \tau_{i-2} }{ \tau_{i} - \tau_{i-1} } & \forall & u > 0 0233 \\ 0234 r = \frac{ \tau_{i+1} - \tau_{i} }{ \tau_{i} - \tau_{i-1} } & \forall & u < 0\end{aligned} 0235 0236 as its argument. There are many choices of limiter function but we 0237 only provide the Superbee limiter (Roe 1995 :cite:`roe:85`): 0238 0239 .. math:: \psi(r) = \max[0,\min[1,2r],\min[2,r]] 0240 0241 .. admonition:: S/R :filelink:`GAD_FLUXLIMIT_ADV_X <pkg/generic_advdiff/gad_fluxlimit_adv_x.F>` 0242 :class: note 0243 0244 | :math:`F_x` : :varlink:`uT` ( argument ) 0245 | :math:`U` : :varlink:`uTrans` ( argument ) 0246 | :math:`\tau` : :varlink:`tracer` ( argument ) 0247 0248 .. admonition:: S/R :filelink:`GAD_FLUXLIMIT_ADV_Y <pkg/generic_advdiff/gad_fluxlimit_adv_y.F>` 0249 :class: note 0250 0251 | :math:`F_y` : :varlink:`vT` ( argument ) 0252 | :math:`V` : :varlink:`vTrans` ( argument ) 0253 | :math:`\tau` : :varlink:`tracer` ( argument ) 0254 0255 .. admonition:: S/R :filelink:`GAD_FLUXLIMIT_ADV_R <pkg/generic_advdiff/gad_fluxlimit_adv_r.F>` 0256 :class: note 0257 0258 | :math:`F_r` : :varlink:`wT` ( argument ) 0259 | :math:`W` : :varlink:`rTrans` ( argument ) 0260 | :math:`\tau` : :varlink:`tracer` ( argument ) 0261 0262 Third order direct space-time 0263 ----------------------------- 0264 0265 The direct space-time method deals with space and time discretization 0266 together (other methods that treat space and time separately are known** Warning **
Wide character in print at /usr/local/share/lxr/source line 1030, <$git> line 268.
0267 collectively as the “Method of Lines”). The Lax-Wendroff scheme falls 0268 into this category; it adds sufficient diffusion to a second order flux 0269 that the forward-in-time method is stable. The upwind biased third order 0270 DST scheme is: 0271 0272 .. math:: 0273 \begin{aligned}F = u \left( \tau_{i-1} 0274 + d_0 (\tau_{i}-\tau_{i-1}) + d_1 (\tau_{i-1}-\tau_{i-2}) \right) 0275 \phantom{W} & \forall & u > 0 \\ 0276 F = u \left( \tau_{i} 0277 - d_0 (\tau_{i}-\tau_{i-1}) - d_1 (\tau_{i+1}-\tau_{i}) \right) 0278 \phantom{W} & \forall & u < 0\end{aligned} 0279 :label: F_posneg-u 0280 0281 where 0282 0283 .. math:: 0284 0285 \begin{aligned} 0286 d_0 & = & \frac{1}{6} ( 2 - |c| ) ( 1 - |c| ) \\ 0287 d_1 & = & \frac{1}{6} ( 1 - |c| ) ( 1 + |c| )\end{aligned} 0288 0289 The coefficients :math:`d_0` and :math:`d_1` approach :math:`1/3` and 0290 :math:`1/6` respectively as the Courant number, :math:`c`, vanishes. In 0291 this limit, the conventional third order upwind method is recovered. For 0292 finite Courant number, the deviations from the linear method are 0293 analogous to the diffusion added to centered second order advection in 0294 the Lax-Wendroff scheme. 0295 0296 The DST3 method described above must be used in a forward-in-time manner 0297 and is stable for :math:`0 \le |c| \le 1`. Although the scheme appears 0298 to be forward-in-time, it is in fact third order in time and the 0299 accuracy increases with the Courant number! For low Courant number, DST3 0300 produces very similar results (indistinguishable in 0301 :numref:`advect-1d-lo`) to the linear third order method but for large 0302 Courant number, where the linear upwind third order method is unstable, 0303 the scheme is extremely accurate (:numref:`advect-1d-hi`) with only 0304 minor overshoots. 0305 0306 .. admonition:: S/R :filelink:`GAD_DST3_ADV_X <pkg/generic_advdiff/gad_dst3_adv_x.F>` 0307 :class: note 0308 0309 | :math:`F_x` : :varlink:`uT` ( argument ) 0310 | :math:`U` : :varlink:`uTrans` ( argument ) 0311 | :math:`\tau` : :varlink:`tracer` ( argument ) 0312 0313 .. admonition:: S/R :filelink:`GAD_DST3_ADV_Y <pkg/generic_advdiff/gad_dst3_adv_y.F>` 0314 :class: note 0315 0316 | :math:`F_y` : :varlink:`vT` ( argument ) 0317 | :math:`V` : :varlink:`vTrans` ( argument ) 0318 | :math:`\tau` : :varlink:`tracer` ( argument ) 0319 0320 .. admonition:: S/R :filelink:`GAD_DST3_ADV_R <pkg/generic_advdiff/gad_dst3_adv_r.F>` 0321 :class: note 0322 0323 | :math:`F_r` : :varlink:`wT` ( argument ) 0324 | :math:`W` : :varlink:`rTrans` ( argument ) 0325 | :math:`\tau` : :varlink:`tracer` ( argument ) 0326 0327 Third order direct space-time with flux limiting 0328 ------------------------------------------------ 0329 0330 The overshoots in the DST3 method can be controlled with a flux limiter. 0331 The limited flux is written: 0332 0333 .. math:: 0334 F = \frac{1}{2}(u+|u|)\left( \tau_{i-1} + \psi(r^+)(\tau_{i} - \tau_{i-1} )\right) 0335 + \frac{1}{2}(u-|u|)\left( \tau_{i-1} + \psi(r^-)(\tau_{i} - \tau_{i-1} )\right) 0336 :label: dst3_limiter 0337 0338 where 0339 0340 .. math:: 0341 0342 \begin{aligned} 0343 r^+ & = & \frac{\tau_{i-1} - \tau_{i-2}}{\tau_{i} - \tau_{i-1}} \\ 0344 r^- & = & \frac{\tau_{i+1} - \tau_{i}}{\tau_{i} - \tau_{i-1}}\end{aligned} 0345 0346 and the limiter is the Sweby limiter: 0347 0348 .. math:: \psi(r) = \max[0, \min[\min(1,d_0+d_1r],\frac{1-c}{c}r ]] 0349 0350 .. admonition:: S/R :filelink:`GAD_DST3FL_ADV_X <pkg/generic_advdiff/gad_dst3fl_adv_x.F>` 0351 :class: note 0352 0353 | :math:`F_x` : :varlink:`uT` ( argument ) 0354 | :math:`U` : :varlink:`uTrans` ( argument ) 0355 | :math:`\tau` : :varlink:`tracer` ( argument ) 0356 0357 .. admonition:: S/R :filelink:`GAD_DST3FL_ADV_Y <pkg/generic_advdiff/gad_dst3fl_adv_y.F>` 0358 :class: note 0359 0360 | :math:`F_y` : :varlink:`vT` ( argument ) 0361 | :math:`V` : :varlink:`vTrans` ( argument ) 0362 | :math:`\tau` : :varlink:`tracer` ( argument ) 0363 0364 .. admonition:: S/R :filelink:`GAD_DST3FL_ADV_R <pkg/generic_advdiff/gad_dst3fl_adv_r.F>` 0365 :class: note 0366 0367 | :math:`F_r` : :varlink:`wT` ( argument ) 0368 | :math:`W` : :varlink:`rTrans` ( argument ) 0369 | :math:`\tau` : :varlink:`tracer` ( argument ) 0370 0371 Multi-dimensional advection 0372 --------------------------- 0373 0374 In many of the aforementioned advection schemes the behavior in multiple 0375 dimensions is not necessarily as good as the one dimensional behavior. 0376 For instance, a shape preserving monotonic scheme in one dimension can 0377 have severe shape distortion in two dimensions if the two components of 0378 horizontal fluxes are treated independently. There is a large body of 0379 literature on the subject dealing with this problem and among the fixes 0380 are operator and flux splitting methods, corner flux methods, and more. 0381 We have adopted a variant on the standard splitting methods that allows 0382 the flux calculations to be implemented as if in one dimension: 0383 0384 .. math:: 0385 \begin{aligned} 0bad585a21 Navi*0386 \tau^{n+1/3} & = \tau^{n} 68aadaa6bd Phob*0387 - \Delta t \left( \frac{1}{\Delta x} \delta_i F^x(\tau^{n}) 0388 - \tau^{n} \frac{1}{\Delta x} \delta_i u \right) \\ 0bad585a21 Navi*0389 \tau^{n+2/3} & = \tau^{n+1/3} 68aadaa6bd Phob*0390 - \Delta t \left( \frac{1}{\Delta y} \delta_j F^y(\tau^{n+1/3}) 0391 - \tau^{n} \frac{1}{\Delta y} \delta_i v \right) \\ 0bad585a21 Navi*0392 \tau^{n+3/3} & = \tau^{n+2/3} 68aadaa6bd Phob*0393 - \Delta t \left( \frac{1}{\Delta r} \delta_k F^x(\tau^{n+2/3}) 0394 - \tau^{n} \frac{1}{\Delta r} \delta_i w \right)\end{aligned} 0395 :label: tau_multiD 0396 0397 In order to incorporate this method into the general model algorithm, we 0398 compute the effective tendency rather than update the tracer so that 0399 other terms such as diffusion are using the :math:`n` time-level and not 0400 the updated :math:`n+3/3` quantities: 0401 0bad585a21 Navi*0402 .. math:: G^{n+1/2}_{\rm adv} = \frac{1}{\Delta t} ( \tau^{n+3/3} - \tau^{n} ) 68aadaa6bd Phob*0403 0404 So that the over all time-stepping looks likes: 0405 0bad585a21 Navi*0406 .. math:: \tau^{n+1} = \tau^{n} + \Delta t \left( G^{n+1/2}_{\rm adv} + G_{\rm diff}(\tau^{n}) + G^{n}_{\rm forcing} \right) 68aadaa6bd Phob*0407 0408 .. admonition:: S/R :filelink:`GAD_ADVECTION <pkg/generic_advdiff/gad_advection.F>` 0409 :class: note 0410 0411 | :math:`\tau` : :varlink:`tracer` ( argument ) 0412 | :math:`G^{n+1/2}_{adv}` : :varlink:`gTracer` ( argument ) 0413 | :math:`F_x, F_y, F_r` : :varlink:`aF` ( local ) 0414 | :math:`U` : :varlink:`uTrans` ( local ) 0415 | :math:`V` : :varlink:`vTrans` ( local ) 0416 | :math:`W` : :varlink:`rTrans` ( local ) 0417 0418 A schematic of multi-dimension time stepping for the cube sphere configuration is show in :numref:`multiDim_CS` . 0419 0420 .. figure:: figs/multiDim_CS.* 0421 :width: 60% 0422 :align: center 0423 :alt: multiDim_CS 0424 :name: multiDim_CS 0425 0426 Multi-dimensional advection time-stepping with cubed-sphere topology. 0427 0428 Comparison of advection schemes 0429 =============================== 0430 c6adb04bff jpga*0431 :numref:`adv_scheme_summary` shows a summary of the different advection schemes available in MITgcm. 68aadaa6bd Phob*** Warning **
Wide character in print at /usr/local/share/lxr/source line 1030, <$git> line 433.
0432 “AB” stands for Adams-Bashforth and “DST” for direct space-time. The code corresponds to the number used 0433 to select the corresponding advection scheme in the parameter file (e.g., ``tempAdvScheme=3`` in file 0434 ``data`` selects the 3rd order upwind advection scheme for temperature advection). 0435 0436 .. table:: MITgcm Advection Schemes 0437 :name: adv_scheme_summary 0438 0439 +-------------------------------------------------------------+------+-----+---------------+---------+---------------------------------------------------+ 0440 | Advection Scheme | Code | Use | Use | Stencil | Comments | 0441 | | | AB? | multi-dim? | (1-D) | | 0442 +=============================================================+======+=====+===============+=========+===================================================+ 0443 | 1st order upwind | 1 | no | yes\ :sup:`*` | 3 | linear :math:`\tau`, non-linear :math:`\vec{v}` | 0444 +-------------------------------------------------------------+------+-----+---------------+---------+---------------------------------------------------+ 0445 | centered 2nd order | 2 | yes | no | 3 | linear | 0446 +-------------------------------------------------------------+------+-----+---------------+---------+---------------------------------------------------+ 0447 | 3rd order upwind | 3 | yes | no | 5 | linear :math:`\tau` | 0448 +-------------------------------------------------------------+------+-----+---------------+---------+---------------------------------------------------+ 0449 | centered 4th order | 4 | yes | no | 5 | linear | 0450 +-------------------------------------------------------------+------+-----+---------------+---------+---------------------------------------------------+ 0451 | 2nd order DST (Lax-Wendroff) | 20 | no | yes\ :sup:`*` | 3 | linear :math:`\tau`, non-linear :math:`\vec{v}` | 0452 +-------------------------------------------------------------+------+-----+---------------+---------+---------------------------------------------------+ 0453 | 3rd order DST | 30 | no | yes\ :sup:`*` | 5 | linear :math:`\tau`, non-linear :math:`\vec{v}` | 0454 +-------------------------------------------------------------+------+-----+---------------+---------+---------------------------------------------------+ 0455 | 2nd order flux limiters | 77 | no | yes\ :sup:`*` | 5 | non-linear | 0456 +-------------------------------------------------------------+------+-----+---------------+---------+---------------------------------------------------+ 0457 | 3rd order DST flux limiter | 33 | no | yes\ :sup:`*` | 5 | non-linear | 0458 +-------------------------------------------------------------+------+-----+---------------+---------+---------------------------------------------------+** Warning **
Wide character in print at /usr/local/share/lxr/source line 1030, <$git> line 460.
0459 | piecewise parabolic w/“null” limiter | 40 | no | yes | 7 | non-linear | 0460 +-------------------------------------------------------------+------+-----+---------------+---------+---------------------------------------------------+** Warning **
Wide character in print at /usr/local/share/lxr/source line 1030, <$git> line 462.
0461 | piecewise parabolic w/“mono” limiter | 41 | no | yes | 7 | non-linear | 0462 +-------------------------------------------------------------+------+-----+---------------+---------+---------------------------------------------------+** Warning **
Wide character in print at /usr/local/share/lxr/source line 1030, <$git> line 464.
0463 | piecewise parabolic w/“weno” limiter | 42 | no | yes | 7 | non-linear | 0464 +-------------------------------------------------------------+------+-----+---------------+---------+---------------------------------------------------+** Warning **
Wide character in print at /usr/local/share/lxr/source line 1030, <$git> line 466.
0465 | piecewise quartic w/“null” limiter | 50 | no | yes | 9 | non-linear | 0466 +-------------------------------------------------------------+------+-----+---------------+---------+---------------------------------------------------+** Warning **
Wide character in print at /usr/local/share/lxr/source line 1030, <$git> line 468.
0467 | piecewise quartic w/“mono” limiter | 51 | no | yes | 9 | non-linear | 0468 +-------------------------------------------------------------+------+-----+---------------+---------+---------------------------------------------------+** Warning **
Wide character in print at /usr/local/share/lxr/source line 1030, <$git> line 470.
0469 | piecewise quartic w/“weno” limiter | 52 | no | yes | 9 | non-linear | 0470 +-------------------------------------------------------------+------+-----+---------------+---------+---------------------------------------------------+ 0471 | 7th order one-step method w/monotonicity preserving limiter | 7 | no | yes | 9 | non-linear | 0472 +-------------------------------------------------------------+------+-----+---------------+---------+---------------------------------------------------+ 0473 | second order-moment Prather | 80 | no | yes | 3 | non-linear | 0474 +-------------------------------------------------------------+------+-----+---------------+---------+---------------------------------------------------+ 0475 | second order-moment Prather w/limiter | 81 | no | yes | 3 | non-linear | 0476 +-------------------------------------------------------------+------+-----+---------------+---------+---------------------------------------------------+ 0477 0478 yes\ :sup:`*` indicates that either the multi-dim advection algorithm or standard approach can be utilized, controlled by a namelist parameter :varlink:`multiDimAdvection` 0479 (in these cases, given that these schemes was designed to use multi-dim advection, using the standard approach is not recommended). 0480 The minimum size of the required tile overlap region (:varlink:`OLx`, :varlink:`OLx`) 0481 is (stencil size -1)/2. The minimum overlap required by the model in general is 2, 0482 so for some of the above choices the advection scheme will not cost anything in terms of an additional overlap requirement, 0483 but especially given a small tile size, using scheme 7 for example would require costly additional overlap points** Warning **
Wide character in print at /usr/local/share/lxr/source line 1030, <$git> line 485.
0484 (note a cube sphere grid with a “wet-corner” requires doubling this overlap!)** Warning **
Wide character in print at /usr/local/share/lxr/source line 1030, <$git> line 486.
0485 In the ‘Comments’ column, :math:`\tau` refers to tracer advection, :math:`\vec{v}` momentum advection. 0486 0487 Shown in :numref:`advect-1d-lo` and :numref:`advect-1d-hi` is a 1-D comparison of advection schemes. Here we advect both a smooth hill and a hill with a more abrupt shock. 0488 :numref:`advect-1d-lo` shown the result for a weak flow (low Courant number) whereas :numref:`advect-1d-hi` shows the result for a stronger flow (high Courant number). 0489 0490 .. figure:: figs/advect-1d-lo.* 0491 :width: 100% 0492 :align: center 0493 :alt: advect-1d-lo 0494 :name: advect-1d-lo 0495 0496 Comparison of 1-D advection schemes: Courant number is 0.05 with 60 points and solutions are shown for T=1 (one complete period). a) Shows the upwind biased schemes; first order upwind, DST3, third order upwind and second order upwind. b) Shows the centered schemes; Lax-Wendroff, DST4, centered second order, centered fourth order and finite volume fourth order. c) Shows the second order flux limiters: minmod, Superbee, MC limiter and the van Leer limiter. d) Shows the DST3 method with flux limiters due to Sweby with :math:`\mu =1` , :math:`\mu =c/(1-c)` and a fourth order DST method with Sweby limiter, :math:`\mu =c/(1-c)` . 0497 0498 .. figure:: figs/advect-1d-hi.* 0499 :width: 100% 0500 :align: center 0501 :alt: advect-1d-hi 0502 :name: advect-1d-hi 0503 0504 Comparison of 1-D advection schemes: Courant number is 0.89 with 60 points and solutions are shown for T=1 (one complete period). a) Shows the upwind biased schemes; first order upwind and DST3. Third order upwind and second order upwind are unstable at this Courant number. b) Shows the centered schemes; Lax-Wendroff, DST4. Centered second order, centered fourth order and finite volume fourth order are unstable at this Courant number. c) Shows the second order flux limiters: minmod, Superbee, MC limiter and the van Leer limiter. d) Shows the DST3 method with flux limiters due to Sweby with :math:`\mu =1` , :math:`\mu =c/(1-c)` and a fourth order DST method with Sweby limiter, :math:`\mu =c/(1-c)` . 0505 0506 :numref:`advect-2d-lo-diag`, :numref:`advect-2d-mid-diag` and 0507 :numref:`advect-2d-hi-diag` show solutions to a simple diagonal advection 0508 problem using a selection of schemes for low, moderate and high Courant 0509 numbers, respectively. The top row shows the linear schemes, integrated 0510 with the Adams-Bashforth method. Theses schemes are clearly unstable for 0511 the high Courant number and weakly unstable for the moderate Courant 0512 number. The presence of false extrema is very apparent for all Courant 0513 numbers. The middle row shows solutions obtained with the unlimited but 0514 multi-dimensional schemes. These solutions also exhibit false extrema 0515 though the pattern now shows symmetry due to the multi-dimensional 0516 scheme. Also, the schemes are stable at high Courant number where the** Warning **
Wide character in print at /usr/local/share/lxr/source line 1030, <$git> line 518.
0517 linear schemes weren’t. The bottom row (left and middle) shows the 0518 limited schemes and most obvious is the absence of false extrema. The 0519 accuracy and stability of the unlimited non-linear schemes is retained 0520 at high Courant number but at low Courant number the tendency is to 0521 lose amplitude in sharp peaks due to diffusion. The one dimensional 0522 tests shown in :numref:`advect-1d-lo` and :numref:`advect-1d-hi` show 0523 this phenomenon. 0524 0525 Finally, the bottom left and right panels use the same advection scheme 0526 but the right does not use the multi-dimensional method. At low Courant 0527 number this appears to not matter but for moderate Courant number severe 0528 distortion of the feature is apparent. Moreover, the stability of the 0529 multi-dimensional scheme is determined by the maximum Courant number 0530 applied of each dimension while the stability of the method of lines is 0531 determined by the sum. Hence, in the high Courant number plot, the 0532 scheme is unstable. 0533 0534 .. figure:: figs/advect-2d-lo-diag.* 0535 :width: 100% 0536 :align: center 0537 :alt: advect-2d-lo-diag 0538 :name: advect-2d-lo-diag 0539 0540 Comparison of advection schemes in two dimensions; diagonal advection of a resolved Gaussian feature. Courant number is 0.01 with 30 :math:`\times` 30 points and solutions are shown for T=1/2. White lines indicate zero crossing (ie. the presence of false minima). The left column shows the second order schemes; top) centered second order with Adams-Bashforth, middle) Lax-Wendroff and bottom) Superbee flux limited. The middle column shows the third order schemes; top) upwind biased third order with Adams-Bashforth, middle) third order direct space-time method and bottom) the same with flux limiting. The top right panel shows the centered fourth order scheme with Adams-Bashforth and right middle panel shows a fourth order variant on the DST method. Bottom right panel shows the Superbee flux limiter (second order) applied independently in each direction (method of lines). 0541 0542 .. figure:: figs/advect-2d-mid-diag.* 0543 :width: 100% 0544 :align: center 0545 :alt: advect-2d-mid-diag 0546 :name: advect-2d-mid-diag 0547 0548 Comparison of advection schemes in two dimensions; diagonal advection of a resolved Gaussian feature. Courant number is 0.27 with 30 :math:`\times` 30 points and solutions are shown for T=1/2. White lines indicate zero crossing (ie. the presence of false minima). The left column shows the second order schemes; top) centered second order with Adams-Bashforth, middle) Lax-Wendroff and bottom) Superbee flux limited. The middle column shows the third order schemes; top) upwind biased third order with Adams-Bashforth, middle) third order direct space-time method and bottom) the same with flux limiting. The top right panel shows the centered fourth order scheme with Adams-Bashforth and right middle panel shows a fourth order variant on the DST method. Bottom right panel shows the Superbee flux limiter (second order) applied independently in each direction (method of lines). 0549 0550 .. figure:: figs/advect-2d-hi-diag.* 0551 :width: 100% 0552 :align: center 0553 :alt: advect-2d-hi-diag 0554 :name: advect-2d-hi-diag 0555 0556 Comparison of advection schemes in two dimensions; diagonal advection of a resolved Gaussian feature. Courant number is 0.47 with 30 :math:`\times` 30 points and solutions are shown for T=1/2. White lines indicate zero crossings and initial maximum values (ie. the presence of false extrema). The left column shows the second order schemes; top) centered second order with Adams-Bashforth, middle) Lax-Wendroff and bottom) Superbee flux limited. The middle column shows the third order schemes; top) upwind biased third order with Adams-Bashforth, middle) third order direct space-time method and bottom) the same with flux limiting. The top right panel shows the centered fourth order scheme with Adams-Bashforth and right middle panel shows a fourth order variant on the DST method. Bottom right panel shows the Superbee flux limiter (second order) applied independently in each direction (method of lines). 0557 0558 With many advection schemes implemented in the code two questions arise:** Warning **
Wide character in print at /usr/local/share/lxr/source line 1030, <$git> line 560.
0559 “Which scheme is best?” and “Why don’t you just offer the best advection** Warning **
Wide character in print at /usr/local/share/lxr/source line 1030, <$git> line 561.
0560 scheme?”. Unfortunately, no one advection scheme is “the best” for all 0561 particular applications and for new applications it is often a matter of 0562 trial to determine which is most suitable. Here are some guidelines but 0563 these are not the rule; 0564 0565 - If you have a coarsely resolved model, using a positive or upwind 0566 biased scheme will introduce significant diffusion to the solution 0567 and using a centered higher order scheme will introduce more noise. 0568 In this case, simplest may be best. 0569 0570 - If you have a high resolution model, using a higher order scheme will 0571 give a more accurate solution but scale-selective diffusion might 0572 need to be employed. The flux limited methods offer similar accuracy 0573 in this regime. 0574 0575 - If your solution has shocks or propagating fronts then a flux limited 0576 scheme is almost essential. 0577 0578 - If your time-step is limited by advection, the multi-dimensional 0579 non-linear schemes have the most stability (up to Courant number 1). 0580 0581 - If you need to know how much diffusion/dissipation has occurred you 0582 will have a lot of trouble figuring it out with a non-linear method. 0583 0584 - The presence of false extrema is non-physical and this alone is the 0585 strongest argument for using a positive scheme. 0586
| [ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
|
This page was automatically generated from https://github.com/darwinproject/darwin3 by the 2.3.7-MITgcm-0.1 LXR engine. The LXR team |
|