PDE-Specialized Explicit Runge-Kutta Methods
Low Storage Explicit Runge-Kutta Methods
OrdinaryDiffEq.CarpenterKennedy2N54
— TypeCarpenterKennedy2N54(; 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} }
OrdinaryDiffEq.SHLDDRK64
— TypeSHLDDRK64(; 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
Missing docstring for SHLDDRK52
. Check Documenter's build log for details.
Missing docstring for SHLDDRK_2N
. Check Documenter's build log for details.
OrdinaryDiffEq.HSLDDRK64
— TypeDeprecated 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.
OrdinaryDiffEq.DGLDDRK73_C
— TypeDGLDDRK73_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
OrdinaryDiffEq.DGLDDRK84_C
— TypeDGLDDRK84_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
OrdinaryDiffEq.DGLDDRK84_F
— TypeDGLDDRK84_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
OrdinaryDiffEq.NDBLSRK124
— TypeNDBLSRK124(; 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
OrdinaryDiffEq.NDBLSRK134
— TypeNDBLSRK134(; 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
OrdinaryDiffEq.NDBLSRK144
— TypeNDBLSRK144(; 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
OrdinaryDiffEq.CFRLDDRK64
— TypeM. 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.
OrdinaryDiffEq.TSLDDRK74
— TypeKostas 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.
OrdinaryDiffEq.CKLLSRK43_2
— TypeCKLLSRK43_2: Low-Storage Method 4-stage, third order low-storage scheme, optimised for compressible Navier–Stokes equations.
OrdinaryDiffEq.CKLLSRK54_3C
— TypeCKLLSRK54_3C: Low-Storage Method 5-stage, fourth order low-storage scheme, optimised for compressible Navier–Stokes equations.
OrdinaryDiffEq.CKLLSRK95_4S
— TypeCKLLSRK95_4S: Low-Storage Method 9-stage, fifth order low-storage scheme, optimised for compressible Navier–Stokes equations.
OrdinaryDiffEq.CKLLSRK95_4C
— TypeCKLLSRK95_4C: Low-Storage Method 9-stage, fifth order low-storage scheme, optimised for compressible Navier–Stokes equations.
OrdinaryDiffEq.CKLLSRK95_4M
— TypeCKLLSRK95_4M: Low-Storage Method 9-stage, fifth order low-storage scheme, optimised for compressible Navier–Stokes equations.
OrdinaryDiffEq.CKLLSRK54_3C_3R
— TypeCKLLSRK543C3R: Low-Storage Method 5-stage, fourth order low-storage scheme, optimised for compressible Navier–Stokes equations.
OrdinaryDiffEq.CKLLSRK54_3M_3R
— TypeCKLLSRK543M3R: Low-Storage Method 5-stage, fourth order low-storage scheme, optimised for compressible Navier–Stokes equations.
OrdinaryDiffEq.CKLLSRK54_3N_3R
— TypeCKLLSRK543N3R: Low-Storage Method 5-stage, fourth order low-storage scheme, optimised for compressible Navier–Stokes equations.
OrdinaryDiffEq.CKLLSRK85_4C_3R
— TypeCKLLSRK854C3R: Low-Storage Method 8-stage, fifth order low-storage scheme, optimised for compressible Navier–Stokes equations.
OrdinaryDiffEq.CKLLSRK85_4M_3R
— TypeCKLLSRK854M3R: Low-Storage Method 8-stage, fifth order low-storage scheme, optimised for compressible Navier–Stokes equations.
OrdinaryDiffEq.CKLLSRK85_4P_3R
— TypeCKLLSRK854P3R: Low-Storage Method 8-stage, fifth order low-storage scheme, optimised for compressible Navier–Stokes equations.
OrdinaryDiffEq.CKLLSRK54_3N_4R
— TypeCKLLSRK543N4R: Low-Storage Method 5-stage, fourth order low-storage scheme, optimised for compressible Navier–Stokes equations.
OrdinaryDiffEq.CKLLSRK54_3M_4R
— TypeCKLLSRK543M4R: Low-Storage Method 5-stage, fourth order low-storage scheme, optimised for compressible Navier–Stokes equations.
OrdinaryDiffEq.CKLLSRK65_4M_4R
— TypeCKLLSRK654M4R: Low-Storage Method 6-stage, fifth order low-storage scheme, optimised for compressible Navier–Stokes equations.
OrdinaryDiffEq.CKLLSRK85_4FM_4R
— TypeCKLLSRK854FM4R: Low-Storage Method 8-stage, fifth order low-storage scheme, optimised for compressible Navier–Stokes equations.
OrdinaryDiffEq.CKLLSRK75_4M_5R
— TypeCKLLSRK754M5R: Low-Storage Method 7-stage, fifth order low-storage scheme, optimised for compressible Navier–Stokes equations.
OrdinaryDiffEq.ParsaniKetchesonDeconinck3S32
— TypeParsani, 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.
OrdinaryDiffEq.ParsaniKetchesonDeconinck3S82
— TypeParsani, 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.
OrdinaryDiffEq.ParsaniKetchesonDeconinck3S53
— TypeParsani, 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.
OrdinaryDiffEq.ParsaniKetchesonDeconinck3S173
— TypeParsani, 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.
OrdinaryDiffEq.ParsaniKetchesonDeconinck3S94
— TypeParsani, 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.
OrdinaryDiffEq.ParsaniKetchesonDeconinck3S184
— TypeParsani, 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.
OrdinaryDiffEq.ParsaniKetchesonDeconinck3S105
— TypeParsani, 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.
OrdinaryDiffEq.ParsaniKetchesonDeconinck3S205
— TypeParsani, 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.
OrdinaryDiffEq.RDPK3Sp35
— TypeRDPK3Sp35(; 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
OrdinaryDiffEq.RDPK3SpFSAL35
— TypeRDPK3SpFSAL35(; 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
OrdinaryDiffEq.RDPK3Sp49
— TypeRDPK3Sp49(; 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
OrdinaryDiffEq.RDPK3SpFSAL49
— TypeRDPK3SpFSAL49(; 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
OrdinaryDiffEq.RDPK3Sp510
— TypeRDPK3Sp510(; 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
OrdinaryDiffEq.RDPK3SpFSAL510
— TypeRDPK3SpFSAL510(; 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
SSP Optimized Runge-Kutta Methods
Missing docstring for KYK2014DGSSPRK_3S2
. Check Documenter's build log for details.
OrdinaryDiffEq.SSPRK22
— TypeSSPRK22(; 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
OrdinaryDiffEq.SSPRK33
— TypeSSPRK33(; 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
OrdinaryDiffEq.SSPRK53
— TypeSSPRK53(; 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.
Missing docstring for KYKSSPRK42
. Check Documenter's build log for details.
OrdinaryDiffEq.SSPRK53_2N1
— TypeSSPRK53_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.
OrdinaryDiffEq.SSPRK53_2N2
— TypeSSPRK53_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.
OrdinaryDiffEq.SSPRK53_H
— TypeSSPRK53_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.
OrdinaryDiffEq.SSPRK63
— TypeSSPRK63(; 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.
OrdinaryDiffEq.SSPRK73
— TypeSSPRK73(; 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.
OrdinaryDiffEq.SSPRK83
— TypeSSPRK83(; 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.
OrdinaryDiffEq.SSPRK43
— TypeSSPRK43(; 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
OrdinaryDiffEq.SSPRK432
— TypeSSPRK432(; 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.
OrdinaryDiffEq.SSPRKMSVS43
— TypeSSPRKMSVS43(; 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
OrdinaryDiffEq.SSPRKMSVS32
— TypeSSPRKMSVS32(; 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
OrdinaryDiffEq.SSPRK932
— TypeSSPRK932(; 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.
OrdinaryDiffEq.SSPRK54
— TypeSSPRK54(; 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.
OrdinaryDiffEq.SSPRK104
— TypeSSPRK104(; 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.