PDE-Specialized Explicit Runge-Kutta Methods

Low Storage Explicit Runge-Kutta Methods

OrdinaryDiffEq.CarpenterKennedy2N54Type
CarpenterKennedy2N54(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
                       step_limiter! = OrdinaryDiffEq.trivial_limiter!,
                       thread = OrdinaryDiffEq.False(),
                       williamson_condition = true)

A fourth-order, five-stage explicit low-storage method of Carpenter and Kennedy (free 3rd order Hermite interpolant). Fixed timestep only. Designed for hyperbolic PDEs (stability properties).

Like SSPRK methods, this method also takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

@article{carpenter1994fourth, title={Fourth-order 2N-storage Runge-Kutta schemes}, author={Carpenter, Mark H and Kennedy, Christopher A}, year={1994} }

source
OrdinaryDiffEq.SHLDDRK64Type
SHLDDRK64(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
            step_limiter! = OrdinaryDiffEq.trivial_limiter!,
            thread = OrdinaryDiffEq.False(),
            williamson_condition = true)

A fourth-order, six-stage explicit low-storage method. Fixed timestep only.

Like SSPRK methods, this method also takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

  • D. Stanescu, W. G. Habashi. 2N-Storage Low Dissipation and Dispersion Runge-Kutta Schemes for Computational Acoustics. Journal of Computational Physics, 143(2), pp 674-681, 1998. doi: https://doi.org/10.1006/jcph.1998.5986
source
Missing docstring.

Missing docstring for SHLDDRK52. Check Documenter's build log for details.

Missing docstring.

Missing docstring for SHLDDRK_2N. Check Documenter's build log for details.

OrdinaryDiffEq.HSLDDRK64Type

Deprecated SHLDDRK64 scheme from 'D. Stanescu, W. G. Habashi. 2N-Storage Low Dissipation and Dispersion Runge-Kutta Schemes for Computational Acoustics'

HSLDDRK64: Low-Storage Method 6-stage, fourth order low-stage, low-dissipation, low-dispersion scheme. Fixed timestep only.

source
OrdinaryDiffEq.DGLDDRK73_CType
DGLDDRK73_C(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
              step_limiter! = OrdinaryDiffEq.trivial_limiter!,
              thread = OrdinaryDiffEq.False(),
              williamson_condition = true)

7-stage, third order low-storage low-dissipation, low-dispersion scheme for discontinuous Galerkin space discretizations applied to wave propagation problems. Optimized for PDE discretizations when maximum spatial step is small due to geometric features of computational domain. Fixed timestep only.

Like SSPRK methods, this method also takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

  • T. Toulorge, W. Desmet. Optimal Runge–Kutta Schemes for Discontinuous Galerkin Space Discretizations Applied to Wave Propagation Problems. Journal of Computational Physics, 231(4), pp 2067-2091, 2012. doi: https://doi.org/10.1016/j.jcp.2011.11.024
source
OrdinaryDiffEq.DGLDDRK84_CType
DGLDDRK84_C(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
              step_limiter! = OrdinaryDiffEq.trivial_limiter!,
              thread = OrdinaryDiffEq.False(),
              williamson_condition = true)

8-stage, fourth order low-storage low-dissipation, low-dispersion scheme for discontinuous Galerkin space discretizations applied to wave propagation problems. Optimized for PDE discretizations when maximum spatial step is small due to geometric features of computational domain. Fixed timestep only.

Like SSPRK methods, this method also takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

  • T. Toulorge, W. Desmet. Optimal Runge–Kutta Schemes for Discontinuous Galerkin Space Discretizations Applied to Wave Propagation Problems. Journal of Computational Physics, 231(4), pp 2067-2091, 2012. doi: https://doi.org/10.1016/j.jcp.2011.11.024
source
OrdinaryDiffEq.DGLDDRK84_FType
DGLDDRK84_F(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
              step_limiter! = OrdinaryDiffEq.trivial_limiter!,
              thread = OrdinaryDiffEq.False(),
              williamson_condition = true)

8-stage, fourth order low-storage low-dissipation, low-dispersion scheme for discontinuous Galerkin space discretizations applied to wave propagation problems. Optimized for PDE discretizations when the maximum spatial step size is not constrained. Fixed timestep only.

Like SSPRK methods, this method also takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

  • T. Toulorge, W. Desmet. Optimal Runge–Kutta Schemes for Discontinuous Galerkin Space Discretizations Applied to Wave Propagation Problems. Journal of Computational Physics, 231(4), pp 2067-2091, 2012. doi: https://doi.org/10.1016/j.jcp.2011.11.024
source
OrdinaryDiffEq.NDBLSRK124Type
NDBLSRK124(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
             step_limiter! = OrdinaryDiffEq.trivial_limiter!,
             thread = OrdinaryDiffEq.False(),
             williamson_condition = true)

12-stage, fourth order low-storage method with optimized stability regions for advection-dominated problems. Fixed timestep only.

Like SSPRK methods, this method also takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

  • Jens Niegemann, Richard Diehl, Kurt Busch. Efficient Low-Storage Runge–Kutta Schemes with Optimized Stability Regions. Journal of Computational Physics, 231, pp 364-372, 2012. doi: https://doi.org/10.1016/j.jcp.2011.09.003
source
OrdinaryDiffEq.NDBLSRK134Type
NDBLSRK134(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
             step_limiter! = OrdinaryDiffEq.trivial_limiter!,
             thread = OrdinaryDiffEq.False(),
             williamson_condition = true)

13-stage, fourth order low-storage method with optimized stability regions for advection-dominated problems. Fixed timestep only.

Like SSPRK methods, this method also takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

  • Jens Niegemann, Richard Diehl, Kurt Busch. Efficient Low-Storage Runge–Kutta Schemes with Optimized Stability Regions. Journal of Computational Physics, 231, pp 364-372, 2012. doi: https://doi.org/10.1016/j.jcp.2011.09.003
source
OrdinaryDiffEq.NDBLSRK144Type
NDBLSRK144(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
             step_limiter! = OrdinaryDiffEq.trivial_limiter!,
             thread = OrdinaryDiffEq.False(),
             williamson_condition = true)

14-stage, fourth order low-storage method with optimized stability regions for advection-dominated problems. Fixed timestep only.

Like SSPRK methods, this method also takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

  • Jens Niegemann, Richard Diehl, Kurt Busch. Efficient Low-Storage Runge–Kutta Schemes with Optimized Stability Regions. Journal of Computational Physics, 231, pp 364-372, 2012. doi: https://doi.org/10.1016/j.jcp.2011.09.003
source
OrdinaryDiffEq.CFRLDDRK64Type

M. Calvo, J. M. Franco, L. Randez. A New Minimum Storage Runge–Kutta Scheme for Computational Acoustics. Journal of Computational Physics, 201, pp 1-12, 2004. doi: https://doi.org/10.1016/j.jcp.2004.05.012

CFRLDDRK64: Low-Storage Method 6-stage, fourth order low-storage, low-dissipation, low-dispersion scheme. Fixed timestep only.

source
OrdinaryDiffEq.TSLDDRK74Type

Kostas Tselios, T. E. Simos. Optimized Runge–Kutta Methods with Minimal Dispersion and Dissipation for Problems arising from Computational Ccoustics. Physics Letters A, 393(1-2), pp 38-47, 2007. doi: https://doi.org/10.1016/j.physleta.2006.10.072

TSLDDRK74: Low-Storage Method 7-stage, fourth order low-storage low-dissipation, low-dispersion scheme with maximal accuracy and stability limit along the imaginary axes. Fixed timestep only.

source
OrdinaryDiffEq.CKLLSRK43_2Type

CKLLSRK43_2: Low-Storage Method 4-stage, third order low-storage scheme, optimised for compressible Navier–Stokes equations.

source
OrdinaryDiffEq.CKLLSRK54_3CType

CKLLSRK54_3C: Low-Storage Method 5-stage, fourth order low-storage scheme, optimised for compressible Navier–Stokes equations.

source
OrdinaryDiffEq.CKLLSRK95_4SType

CKLLSRK95_4S: Low-Storage Method 9-stage, fifth order low-storage scheme, optimised for compressible Navier–Stokes equations.

source
OrdinaryDiffEq.CKLLSRK95_4CType

CKLLSRK95_4C: Low-Storage Method 9-stage, fifth order low-storage scheme, optimised for compressible Navier–Stokes equations.

source
OrdinaryDiffEq.CKLLSRK95_4MType

CKLLSRK95_4M: Low-Storage Method 9-stage, fifth order low-storage scheme, optimised for compressible Navier–Stokes equations.

source
OrdinaryDiffEq.ParsaniKetchesonDeconinck3S32Type

Parsani, Matteo, David I. Ketcheson, and W. Deconinck. "Optimized explicit Runge–Kutta schemes for the spectral difference method applied to wave propagation problems." SIAM Journal on Scientific Computing 35.2 (2013): A957-A986. doi: https://doi.org/10.1137/120885899

ParsaniKetchesonDeconinck3S32: Low-Storage Method 3-stage, second order (3S) low-storage scheme, optimised for for the spectral difference method applied to wave propagation problems.

source
OrdinaryDiffEq.ParsaniKetchesonDeconinck3S82Type

Parsani, Matteo, David I. Ketcheson, and W. Deconinck. "Optimized explicit Runge–Kutta schemes for the spectral difference method applied to wave propagation problems." SIAM Journal on Scientific Computing 35.2 (2013): A957-A986. doi: https://doi.org/10.1137/120885899

ParsaniKetchesonDeconinck3S82: Low-Storage Method 8-stage, second order (3S) low-storage scheme, optimised for for the spectral difference method applied to wave propagation problems.

source
OrdinaryDiffEq.ParsaniKetchesonDeconinck3S53Type

Parsani, Matteo, David I. Ketcheson, and W. Deconinck. "Optimized explicit Runge–Kutta schemes for the spectral difference method applied to wave propagation problems." SIAM Journal on Scientific Computing 35.2 (2013): A957-A986. doi: https://doi.org/10.1137/120885899

ParsaniKetchesonDeconinck3S53: Low-Storage Method 5-stage, third order (3S) low-storage scheme, optimised for for the spectral difference method applied to wave propagation problems.

source
OrdinaryDiffEq.ParsaniKetchesonDeconinck3S173Type

Parsani, Matteo, David I. Ketcheson, and W. Deconinck. "Optimized explicit Runge–Kutta schemes for the spectral difference method applied to wave propagation problems." SIAM Journal on Scientific Computing 35.2 (2013): A957-A986. doi: https://doi.org/10.1137/120885899

ParsaniKetchesonDeconinck3S173: Low-Storage Method 17-stage, third order (3S) low-storage scheme, optimised for for the spectral difference method applied to wave propagation problems.

source
OrdinaryDiffEq.ParsaniKetchesonDeconinck3S94Type

Parsani, Matteo, David I. Ketcheson, and W. Deconinck. "Optimized explicit Runge–Kutta schemes for the spectral difference method applied to wave propagation problems." SIAM Journal on Scientific Computing 35.2 (2013): A957-A986. doi: https://doi.org/10.1137/120885899

ParsaniKetchesonDeconinck3S94: Low-Storage Method 9-stage, fourth order (3S) low-storage scheme, optimised for for the spectral difference method applied to wave propagation problems.

source
OrdinaryDiffEq.ParsaniKetchesonDeconinck3S184Type

Parsani, Matteo, David I. Ketcheson, and W. Deconinck. "Optimized explicit Runge–Kutta schemes for the spectral difference method applied to wave propagation problems." SIAM Journal on Scientific Computing 35.2 (2013): A957-A986. doi: https://doi.org/10.1137/120885899

ParsaniKetchesonDeconinck3S184: Low-Storage Method 18-stage, fourth order (3S) low-storage scheme, optimised for for the spectral difference method applied to wave propagation problems.

source
OrdinaryDiffEq.ParsaniKetchesonDeconinck3S105Type

Parsani, Matteo, David I. Ketcheson, and W. Deconinck. "Optimized explicit Runge–Kutta schemes for the spectral difference method applied to wave propagation problems." SIAM Journal on Scientific Computing 35.2 (2013): A957-A986. doi: https://doi.org/10.1137/120885899

ParsaniKetchesonDeconinck3S105: Low-Storage Method 10-stage, fifth order (3S) low-storage scheme, optimised for for the spectral difference method applied to wave propagation problems.

source
OrdinaryDiffEq.ParsaniKetchesonDeconinck3S205Type

Parsani, Matteo, David I. Ketcheson, and W. Deconinck. "Optimized explicit Runge–Kutta schemes for the spectral difference method applied to wave propagation problems." SIAM Journal on Scientific Computing 35.2 (2013): A957-A986. doi: https://doi.org/10.1137/120885899

ParsaniKetchesonDeconinck3S205: Low-Storage Method 20-stage, fifth order (3S) low-storage scheme, optimised for for the spectral difference method applied to wave propagation problems.

source
OrdinaryDiffEq.RDPK3Sp35Type
RDPK3Sp35(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
            step_limiter! = OrdinaryDiffEq.trivial_limiter!,
            thread = OrdinaryDiffEq.False())

A third-order, five-stage explicit Runge-Kutta method with embedded error estimator designed for spectral element discretizations of compressible fluid mechanics.

Like SSPRK methods, this method also takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

  • Ranocha, Dalcin, Parsani, Ketcheson (2021) Optimized Runge-Kutta Methods with Automatic Step Size Control for Compressible Computational Fluid Dynamics arXiv:2104.06836
source
OrdinaryDiffEq.RDPK3SpFSAL35Type
RDPK3SpFSAL35(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
                step_limiter! = OrdinaryDiffEq.trivial_limiter!,
                thread = OrdinaryDiffEq.False())

A third-order, five-stage explicit Runge-Kutta method with embedded error estimator using the FSAL property designed for spectral element discretizations of compressible fluid mechanics.

Like SSPRK methods, this method also takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

  • Ranocha, Dalcin, Parsani, Ketcheson (2021) Optimized Runge-Kutta Methods with Automatic Step Size Control for Compressible Computational Fluid Dynamics arXiv:2104.06836
source
OrdinaryDiffEq.RDPK3Sp49Type
RDPK3Sp49(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
            step_limiter! = OrdinaryDiffEq.trivial_limiter!,
            thread = OrdinaryDiffEq.False())

A fourth-order, nine-stage explicit Runge-Kutta method with embedded error estimator designed for spectral element discretizations of compressible fluid mechanics.

Like SSPRK methods, this method also takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

  • Ranocha, Dalcin, Parsani, Ketcheson (2021) Optimized Runge-Kutta Methods with Automatic Step Size Control for Compressible Computational Fluid Dynamics arXiv:2104.06836
source
OrdinaryDiffEq.RDPK3SpFSAL49Type
RDPK3SpFSAL49(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
                step_limiter! = OrdinaryDiffEq.trivial_limiter!,
                thread = OrdinaryDiffEq.False())

A fourth-order, nine-stage explicit Runge-Kutta method with embedded error estimator using the FSAL property designed for spectral element discretizations of compressible fluid mechanics.

Like SSPRK methods, this method also takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

  • Ranocha, Dalcin, Parsani, Ketcheson (2021) Optimized Runge-Kutta Methods with Automatic Step Size Control for Compressible Computational Fluid Dynamics arXiv:2104.06836
source
OrdinaryDiffEq.RDPK3Sp510Type
RDPK3Sp510(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
             step_limiter! = OrdinaryDiffEq.trivial_limiter!,
             thread = OrdinaryDiffEq.False())

A fifth-order, ten-stage explicit Runge-Kutta method with embedded error estimator designed for spectral element discretizations of compressible fluid mechanics.

Like SSPRK methods, this method also takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

  • Ranocha, Dalcin, Parsani, Ketcheson (2021) Optimized Runge-Kutta Methods with Automatic Step Size Control for Compressible Computational Fluid Dynamics arXiv:2104.06836
source
OrdinaryDiffEq.RDPK3SpFSAL510Type
RDPK3SpFSAL510(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
                 step_limiter! = OrdinaryDiffEq.trivial_limiter!,
                 thread = OrdinaryDiffEq.False())

A fifth-order, ten-stage explicit Runge-Kutta method with embedded error estimator using the FSAL property designed for spectral element discretizations of compressible fluid mechanics.

Like SSPRK methods, this method also takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

  • Ranocha, Dalcin, Parsani, Ketcheson (2021) Optimized Runge-Kutta Methods with Automatic Step Size Control for Compressible Computational Fluid Dynamics arXiv:2104.06836
source

SSP Optimized Runge-Kutta Methods

Missing docstring.

Missing docstring for KYK2014DGSSPRK_3S2. Check Documenter's build log for details.

OrdinaryDiffEq.SSPRK22Type
SSPRK22(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
          step_limiter! = OrdinaryDiffEq.trivial_limiter!,
          thread = OrdinaryDiffEq.False())

A second-order, two-stage explicit strong stability preserving (SSP) method. Fixed timestep only.

Like all SSPRK methods, this method takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

  • Shu, Chi-Wang, and Stanley Osher. "Efficient implementation of essentially non-oscillatory shock-capturing schemes." Journal of Computational Physics 77.2 (1988): 439-471. https://doi.org/10.1016/0021-9991(88)90177-5
source
OrdinaryDiffEq.SSPRK33Type
SSPRK33(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
          step_limiter! = OrdinaryDiffEq.trivial_limiter!,
          thread = OrdinaryDiffEq.False())

A third-order, three-stage explicit strong stability preserving (SSP) method. Fixed timestep only.

Like all SSPRK methods, this method takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

  • Shu, Chi-Wang, and Stanley Osher. "Efficient implementation of essentially non-oscillatory shock-capturing schemes." Journal of Computational Physics 77.2 (1988): 439-471. https://doi.org/10.1016/0021-9991(88)90177-5
source
OrdinaryDiffEq.SSPRK53Type
SSPRK53(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
          step_limiter! = OrdinaryDiffEq.trivial_limiter!,
          thread = OrdinaryDiffEq.False())

A third-order, five-stage explicit strong stability preserving (SSP) method. Fixed timestep only.

Like all SSPRK methods, this method takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

  • Ruuth, Steven. "Global optimization of explicit strong-stability-preserving Runge-Kutta methods." Mathematics of Computation 75.253 (2006): 183-207.
source
Missing docstring.

Missing docstring for KYKSSPRK42. Check Documenter's build log for details.

OrdinaryDiffEq.SSPRK53_2N1Type
SSPRK53_2N1(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
              step_limiter! = OrdinaryDiffEq.trivial_limiter!,
              thread = OrdinaryDiffEq.False())

A third-order, five-stage explicit strong stability preserving (SSP) low-storage method. Fixed timestep only.

Like all SSPRK methods, this method takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

  • Higueras and T. Roldán. "New third order low-storage SSP explicit Runge–Kutta methods". arXiv:1809.04807v1.
source
OrdinaryDiffEq.SSPRK53_2N2Type
SSPRK53_2N2(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
              step_limiter! = OrdinaryDiffEq.trivial_limiter!,
              thread = OrdinaryDiffEq.False())

A third-order, five-stage explicit strong stability preserving (SSP) low-storage method. Fixed timestep only.

Like all SSPRK methods, this method takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

  • Higueras and T. Roldán. "New third order low-storage SSP explicit Runge–Kutta methods". arXiv:1809.04807v1.
source
OrdinaryDiffEq.SSPRK53_HType
SSPRK53_H(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
            step_limiter! = OrdinaryDiffEq.trivial_limiter!,
            thread = OrdinaryDiffEq.False())

A third-order, five-stage explicit strong stability preserving (SSP) low-storage method. Fixed timestep only.

Like all SSPRK methods, this method takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

  • Higueras and T. Roldán. "New third order low-storage SSP explicit Runge–Kutta methods". arXiv:1809.04807v1.
source
OrdinaryDiffEq.SSPRK63Type
SSPRK63(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
          step_limiter! = OrdinaryDiffEq.trivial_limiter!,
          thread = OrdinaryDiffEq.False())

A third-order, six-stage explicit strong stability preserving (SSP) method. Fixed timestep only.

Like all SSPRK methods, this method takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

  • Ruuth, Steven. "Global optimization of explicit strong-stability-preserving Runge-Kutta methods." Mathematics of Computation 75.253 (2006): 183-207.
source
OrdinaryDiffEq.SSPRK73Type
SSPRK73(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
          step_limiter! = OrdinaryDiffEq.trivial_limiter!,
          thread = OrdinaryDiffEq.False())

A third-order, seven-stage explicit strong stability preserving (SSP) method. Fixed timestep only.

Like all SSPRK methods, this method takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

  • Ruuth, Steven. "Global optimization of explicit strong-stability-preserving Runge-Kutta methods." Mathematics of Computation 75.253 (2006): 183-207.
source
OrdinaryDiffEq.SSPRK83Type
SSPRK83(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
          step_limiter! = OrdinaryDiffEq.trivial_limiter!,
          thread = OrdinaryDiffEq.False())

A third-order, eight-stage explicit strong stability preserving (SSP) method. Fixed timestep only.

Like all SSPRK methods, this method takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

  • Ruuth, Steven. "Global optimization of explicit strong-stability-preserving Runge-Kutta methods." Mathematics of Computation 75.253 (2006): 183-207.
source
OrdinaryDiffEq.SSPRK43Type
SSPRK43(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
          step_limiter! = OrdinaryDiffEq.trivial_limiter!,
          thread = OrdinaryDiffEq.False())

A third-order, four-stage explicit strong stability preserving (SSP) method.

Like all SSPRK methods, this method takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

Optimal third-order explicit SSP method with four stages discovered by

  • J. F. B. M. Kraaijevanger. "Contractivity of Runge-Kutta methods." In: BIT Numerical Mathematics 31.3 (1991), pp. 482–528. DOI: 10.1007/BF01933264.

Embedded method constructed by

  • Sidafa Conde, Imre Fekete, John N. Shadid. "Embedded error estimation and adaptive step-size control for optimal explicit strong stability preserving Runge–Kutta methods." arXiv: 1806.08693

Efficient implementation (and optimized controller) developed by

  • Hendrik Ranocha, Lisandro Dalcin, Matteo Parsani, David I. Ketcheson (2021) Optimized Runge-Kutta Methods with Automatic Step Size Control for Compressible Computational Fluid Dynamics arXiv:2104.06836
source
OrdinaryDiffEq.SSPRK432Type
SSPRK432(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
           step_limiter! = OrdinaryDiffEq.trivial_limiter!,
           thread = OrdinaryDiffEq.False())

A third-order, four-stage explicit strong stability preserving (SSP) method.

Consider using SSPRK43 instead, which uses the same main method and an improved embedded method.

Like all SSPRK methods, this method takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

  • Gottlieb, Sigal, David I. Ketcheson, and Chi-Wang Shu. Strong stability preserving Runge-Kutta and multistep time discretizations. World Scientific, 2011. Example 6.1.
source
OrdinaryDiffEq.SSPRKMSVS43Type
SSPRKMSVS43(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
              step_limiter! = OrdinaryDiffEq.trivial_limiter!)

A third-order, four-step explicit strong stability preserving (SSP) linear multistep method. This method does not come with an error estimator and requires a fixed time step size.

Like all SSP methods, this method takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

References

  • Shu, Chi-Wang. "Total-variation-diminishing time discretizations." SIAM Journal on Scientific and Statistical Computing 9, no. 6 (1988): 1073-1084. DOI: 10.1137/0909073
source
OrdinaryDiffEq.SSPRKMSVS32Type
SSPRKMSVS32(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
              step_limiter! = OrdinaryDiffEq.trivial_limiter!)

A second-order, three-step explicit strong stability preserving (SSP) linear multistep method. This method does not come with an error estimator and requires a fixed time step size.

Like all SSP methods, this method takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

References

  • Shu, Chi-Wang. "Total-variation-diminishing time discretizations." SIAM Journal on Scientific and Statistical Computing 9, no. 6 (1988): 1073-1084. DOI: 10.1137/0909073
source
OrdinaryDiffEq.SSPRK932Type
SSPRK932(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
           step_limiter! = OrdinaryDiffEq.trivial_limiter!,
           thread = OrdinaryDiffEq.False())

A third-order, nine-stage explicit strong stability preserving (SSP) method.

Consider using SSPRK43 instead, which uses the same main method and an improved embedded method.

Like all SSPRK methods, this method takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

  • Gottlieb, Sigal, David I. Ketcheson, and Chi-Wang Shu. Strong stability preserving Runge-Kutta and multistep time discretizations. World Scientific, 2011.
source
OrdinaryDiffEq.SSPRK54Type
SSPRK54(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
          step_limiter! = OrdinaryDiffEq.trivial_limiter!,
          thread = OrdinaryDiffEq.False())

A fourth-order, five-stage explicit strong stability preserving (SSP) method. Fixed timestep only.

Like all SSPRK methods, this method takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

  • Ruuth, Steven. "Global optimization of explicit strong-stability-preserving Runge-Kutta methods." Mathematics of Computation 75.253 (2006): 183-207.
source
OrdinaryDiffEq.SSPRK104Type
SSPRK104(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
           step_limiter! = OrdinaryDiffEq.trivial_limiter!,
           thread = OrdinaryDiffEq.False())

A fourth-order, ten-stage explicit strong stability preserving (SSP) method. Fixed timestep only.

Like all SSPRK methods, this method takes optional arguments stage_limiter! and step_limiter!, where stage_limiter! and step_limiter! are functions of the form limiter!(u, integrator, p, t).

The argument thread determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False(), default) or use multiple threads (thread = OrdinaryDiffEq.True()) when Julia is started with multiple threads.

References

  • Ketcheson, David I. "Highly efficient strong stability-preserving Runge–Kutta methods with low-storage implementations." SIAM Journal on Scientific Computing 30.4 (2008): 2113-2136.
source