OrdinaryDiffEqLowStorageRK
Low-storage Runge-Kutta methods are specialized explicit schemes designed to minimize memory requirements while maintaining high-order accuracy. These methods are essential for large-scale computational fluid dynamics and wave propagation problems where memory constraints are critical.
Key Properties
Low-storage RK methods provide:
- Drastically reduced memory usage (typically 2-4 registers vs 7-10 for standard RK)
- High-order accuracy comparable to standard RK methods
- Preservation of important stability properties (low dissipation/dispersion)
- Scalability to very large PDE discretizations
When to Use Low-Storage RK Methods
These methods are recommended for:
- Large-scale PDE discretizations where memory is the limiting factor
- Computational fluid dynamics and wave propagation simulations
- High-performance computing applications with memory constraints
- GPU computations where memory bandwidth is critical
- Compressible flow simulations and aerodynamics
- Seismic wave propagation and acoustic simulations
- Problems with millions or billions of unknowns
Memory Efficiency Comparison
Registers refer to the number of copies of the u0
vector that must be stored in memory during integration:
- Standard Tsit5: ~9 registers (copies of the state vector)
- Low-storage methods: 2-4 registers (copies of the state vector)
- Practical example: If
u0
is from a PDE semi-discretization requiring 2 GB, then Tsit5 needs 18 GB of working memory, while a 2-register method only needs 4 GB and can achieve the same order - Trade-off: These methods achieve memory reduction by being less computationally efficient, trading compute performance for lower memory requirements
Solver Selection Guide
General-purpose low-storage
CarpenterKennedy2N54
: Fourth-order, 5-stage, excellent general choiceRDPK3Sp510
: Fifth-order with only 3 registers, very memory efficient
Wave propagation optimized
ORK256
: Second-order, 5-stage, optimized for wave equationsCFRLDDRK64
: Low-dissipation and low-dispersion variantTSLDDRK74
: Seventh-order for high accuracy wave propagation
Discontinuous Galerkin optimized
DGLDDRK73_C
: Optimized for DG discretizations (constrained)DGLDDRK84_C
,DGLDDRK84_F
: Fourth-order DG variants
Specialized high-order
NDBLSRK124
,NDBLSRK144
: Multi-stage fourth-order methodsSHLDDRK64
: Low dissipation and dispersion propertiesRK46NL
: Six-stage fourth-order method
Computational fluid dynamics
- Carpenter-Kennedy-Lewis series (
CKLLSRK*
): Optimized for Navier-Stokes equations - Parsani-Ketcheson-Deconinck series (
ParsaniKetcheson*
): CFD-optimized variants - Ranocha-Dalcin-Parsani-Ketcheson series (
RDPK*
): Modern CFD methods
Performance Considerations
- Use only when memory-bound: Standard RK methods are often more efficient when memory is not limiting
- Best for large systems: Most beneficial for problems with >10⁶ unknowns
- GPU acceleration: Particularly effective on memory-bandwidth limited hardware
Installation
To be able to access the solvers in OrdinaryDiffEqLowStorageRK
, you must first install them use the Julia package manager:
using Pkg
Pkg.add("OrdinaryDiffEqLowStorageRK")
This will only install the solvers listed at the bottom of this page. If you want to explore other solvers for your problem, you will need to install some of the other libraries listed in the navigation bar on the left.
Example usage
using OrdinaryDiffEqLowStorageRK
function lorenz!(du, u, p, t)
du[1] = 10.0 * (u[2] - u[1])
du[2] = u[1] * (28.0 - u[3]) - u[2]
du[3] = u[1] * u[2] - (8 / 3) * u[3]
end
u0 = [1.0; 0.0; 0.0]
tspan = (0.0, 100.0)
prob = ODEProblem(lorenz!, u0, tspan)
sol = solve(prob, CarpenterKennedy2N54())
Full list of solvers
OrdinaryDiffEqLowStorageRK.ORK256
— TypeORK256(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False(),
williamson_condition = true)
Explicit Runge-Kutta Method. A second-order, five-stage method for wave propagation equations. Fixed timestep only.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.williamson_condition
: allows for an optimization that allows fusing broadcast expressions with the function callf
. However, it only works forArray
types.
References
Matteo Bernardini, Sergio Pirozzoli. A General Strategy for the Optimization of Runge-Kutta Schemes for Wave Propagation Phenomena. Journal of Computational Physics, 228(11), pp 4182-4199, 2009. doi: https://doi.org/10.1016/j.jcp.2009.02.032
OrdinaryDiffEqLowStorageRK.DGLDDRK73_C
— TypeDGLDDRK73_C(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False(),
williamson_condition = true)
Explicit Runge-Kutta Method. 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.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.williamson_condition
: allows for an optimization that allows fusing broadcast expressions with the function callf
. However, it only works forArray
types.
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
OrdinaryDiffEqLowStorageRK.CarpenterKennedy2N54
— TypeCarpenterKennedy2N54(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False(),
williamson_condition = true)
Explicit Runge-Kutta Method. A fourth-order, five-stage low-storage method of Carpenter and Kennedy (free 3rd order Hermite interpolant). Fixed timestep only. Designed for hyperbolic PDEs (stability properties).
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.williamson_condition
: allows for an optimization that allows fusing broadcast expressions with the function callf
. However, it only works forArray
types.
References
@article{carpenter1994fourth, title={Fourth-order 2N-storage Runge-Kutta schemes}, author={Carpenter, Mark H and Kennedy, Christopher A}, year={1994} }
OrdinaryDiffEqLowStorageRK.NDBLSRK124
— TypeNDBLSRK124(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False(),
williamson_condition = true)
Explicit Runge-Kutta Method. 12-stage, fourth order low-storage method with optimized stability regions for advection-dominated problems. Fixed timestep only.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.williamson_condition
: allows for an optimization that allows fusing broadcast expressions with the function callf
. However, it only works forArray
types.
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
OrdinaryDiffEqLowStorageRK.NDBLSRK144
— TypeNDBLSRK144(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False(),
williamson_condition = true)
Explicit Runge-Kutta Method. 14-stage, fourth order low-storage method with optimized stability regions for advection-dominated problems. Fixed timestep only.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.williamson_condition
: allows for an optimization that allows fusing broadcast expressions with the function callf
. However, it only works forArray
types.
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
OrdinaryDiffEqLowStorageRK.CFRLDDRK64
— TypeCFRLDDRK64(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Low-Storage Method 6-stage, fourth order low-storage, low-dissipation, low-dispersion scheme. Fixed timestep only.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
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
OrdinaryDiffEqLowStorageRK.TSLDDRK74
— TypeTSLDDRK74(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. 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.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
Kostas Tselios, T. E. Simos. Optimized Runge–Kutta Methods with Minimal Dispersion and Dissipation for Problems arising from Computational Acoustics. Physics Letters A, 393(1-2), pp 38-47, 2007. doi: https://doi.org/10.1016/j.physleta.2006.10.072
OrdinaryDiffEqLowStorageRK.DGLDDRK84_C
— TypeDGLDDRK84_C(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False(),
williamson_condition = true)
Explicit Runge-Kutta Method. 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.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.williamson_condition
: allows for an optimization that allows fusing broadcast expressions with the function callf
. However, it only works forArray
types.
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
OrdinaryDiffEqLowStorageRK.DGLDDRK84_F
— TypeDGLDDRK84_F(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False(),
williamson_condition = true)
Explicit Runge-Kutta Method. 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.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.williamson_condition
: allows for an optimization that allows fusing broadcast expressions with the function callf
. However, it only works forArray
types.
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
OrdinaryDiffEqLowStorageRK.SHLDDRK64
— TypeSHLDDRK64(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False(),
williamson_condition = true)
Explicit Runge-Kutta Method. A fourth-order, six-stage low-storage method. Fixed timestep only.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.williamson_condition
: allows for an optimization that allows fusing broadcast expressions with the function callf
. However, it only works forArray
types.
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 }
OrdinaryDiffEqLowStorageRK.RK46NL
— TypeRK46NL(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. 6-stage, fourth order low-stage, low-dissipation, low-dispersion scheme. Fixed timestep only.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
Julien Berland, Christophe Bogey, Christophe Bailly. Low-Dissipation and Low-Dispersion Fourth-Order Runge-Kutta Algorithm. Computers & Fluids, 35(10), pp 1459-1463, 2006. doi: https://doi.org/10.1016/j.compfluid.2005.04.003
OrdinaryDiffEqLowStorageRK.ParsaniKetchesonDeconinck3S32
— TypeParsaniKetchesonDeconinck3S32(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Low-Storage Method 3-stage, second order (3S) low-storage scheme, optimized the spectral difference method applied to wave propagation problems.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
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
OrdinaryDiffEqLowStorageRK.ParsaniKetchesonDeconinck3S82
— TypeParsaniKetchesonDeconinck3S82(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Low-Storage Method 8-stage, second order (3S) low-storage scheme, optimized for the spectral difference method applied to wave propagation problems.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
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
OrdinaryDiffEqLowStorageRK.ParsaniKetchesonDeconinck3S53
— TypeParsaniKetchesonDeconinck3S53(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Low-Storage Method 5-stage, third order (3S) low-storage scheme, optimized for the spectral difference method applied to wave propagation problems.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
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
OrdinaryDiffEqLowStorageRK.ParsaniKetchesonDeconinck3S173
— TypeParsaniKetchesonDeconinck3S173(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Low-Storage Method 17-stage, third order (3S) low-storage scheme, optimized for the spectral difference method applied to wave propagation problems.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
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
OrdinaryDiffEqLowStorageRK.ParsaniKetchesonDeconinck3S94
— TypeParsaniKetchesonDeconinck3S94(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Low-Storage Method 9-stage, fourth order (3S) low-storage scheme, optimized for the spectral difference method applied to wave propagation problems.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
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
OrdinaryDiffEqLowStorageRK.ParsaniKetchesonDeconinck3S184
— TypeParsaniKetchesonDeconinck3S184(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Low-Storage Method 18-stage, fourth order (3S) low-storage scheme, optimized for the spectral difference method applied to wave propagation problems.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
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
OrdinaryDiffEqLowStorageRK.ParsaniKetchesonDeconinck3S105
— TypeParsaniKetchesonDeconinck3S105(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Low-Storage Method 10-stage, fifth order (3S) low-storage scheme, optimized for the spectral difference method applied to wave propagation problems.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
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
OrdinaryDiffEqLowStorageRK.ParsaniKetchesonDeconinck3S205
— TypeParsaniKetchesonDeconinck3S205(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Low-Storage Method 20-stage, fifth order (3S) low-storage scheme, optimized for the spectral difference method applied to wave propagation problems.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
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
OrdinaryDiffEqLowStorageRK.CKLLSRK43_2
— TypeCKLLSRK43_2(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Low-Storage Method 4-stage, third order low-storage scheme, optimized for compressible Navier–Stokes equations.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
@article{kennedy2000low, title={Low-storage, explicit Runge–Kutta schemes for the compressible Navier–Stokes equations}, author={Kennedy, Christopher A and Carpenter, Mark H and Lewis, R Michael}, journal={Applied numerical mathematics}, volume={35}, number={3}, pages={177–219}, year={2000}, publisher={Elsevier}}
OrdinaryDiffEqLowStorageRK.CKLLSRK54_3C
— TypeCKLLSRK54_3C(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Low-Storage Method 5-stage, fourth order low-storage scheme, optimized for compressible Navier–Stokes equations.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
@article{kennedy2000low, title={Low-storage, explicit Runge–Kutta schemes for the compressible Navier–Stokes equations}, author={Kennedy, Christopher A and Carpenter, Mark H and Lewis, R Michael}, journal={Applied numerical mathematics}, volume={35}, number={3}, pages={177–219}, year={2000}, publisher={Elsevier}}
OrdinaryDiffEqLowStorageRK.CKLLSRK95_4S
— TypeCKLLSRK95_4S(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Low-Storage Method 9-stage, fifth order low-storage scheme, optimized for compressible Navier–Stokes equations.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
@article{kennedy2000low, title={Low-storage, explicit Runge–Kutta schemes for the compressible Navier–Stokes equations}, author={Kennedy, Christopher A and Carpenter, Mark H and Lewis, R Michael}, journal={Applied numerical mathematics}, volume={35}, number={3}, pages={177–219}, year={2000}, publisher={Elsevier}}
OrdinaryDiffEqLowStorageRK.CKLLSRK95_4C
— TypeCKLLSRK95_4C(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Low-Storage Method 9-stage, fifth order low-storage scheme, optimized for compressible Navier–Stokes equations.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
@article{kennedy2000low, title={Low-storage, explicit Runge–Kutta schemes for the compressible Navier–Stokes equations}, author={Kennedy, Christopher A and Carpenter, Mark H and Lewis, R Michael}, journal={Applied numerical mathematics}, volume={35}, number={3}, pages={177–219}, year={2000}, publisher={Elsevier}}
OrdinaryDiffEqLowStorageRK.CKLLSRK95_4M
— TypeCKLLSRK95_4M(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Low-Storage Method 9-stage, fifth order low-storage scheme, optimized for compressible Navier–Stokes equations.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
@article{kennedy2000low, title={Low-storage, explicit Runge–Kutta schemes for the compressible Navier–Stokes equations}, author={Kennedy, Christopher A and Carpenter, Mark H and Lewis, R Michael}, journal={Applied numerical mathematics}, volume={35}, number={3}, pages={177–219}, year={2000}, publisher={Elsevier}}
OrdinaryDiffEqLowStorageRK.CKLLSRK54_3C_3R
— TypeCKLLSRK54_3C_3R(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Low-Storage Method 5-stage, fourth order low-storage scheme, optimized for compressible Navier–Stokes equations.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
@article{kennedy2000low, title={Low-storage, explicit Runge–Kutta schemes for the compressible Navier–Stokes equations}, author={Kennedy, Christopher A and Carpenter, Mark H and Lewis, R Michael}, journal={Applied numerical mathematics}, volume={35}, number={3}, pages={177–219}, year={2000}, publisher={Elsevier}}
OrdinaryDiffEqLowStorageRK.CKLLSRK54_3M_3R
— TypeCKLLSRK54_3M_3R(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Low-Storage Method 5-stage, fourth order low-storage scheme, optimized for compressible Navier–Stokes equations.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
@article{kennedy2000low, title={Low-storage, explicit Runge–Kutta schemes for the compressible Navier–Stokes equations}, author={Kennedy, Christopher A and Carpenter, Mark H and Lewis, R Michael}, journal={Applied numerical mathematics}, volume={35}, number={3}, pages={177–219}, year={2000}, publisher={Elsevier}}
OrdinaryDiffEqLowStorageRK.CKLLSRK54_3N_3R
— TypeCKLLSRK54_3N_3R(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Low-Storage Method 5-stage, fourth order low-storage scheme, optimized for compressible Navier–Stokes equations.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
@article{kennedy2000low, title={Low-storage, explicit Runge–Kutta schemes for the compressible Navier–Stokes equations}, author={Kennedy, Christopher A and Carpenter, Mark H and Lewis, R Michael}, journal={Applied numerical mathematics}, volume={35}, number={3}, pages={177–219}, year={2000}, publisher={Elsevier}}
OrdinaryDiffEqLowStorageRK.CKLLSRK85_4C_3R
— TypeCKLLSRK85_4C_3R(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Low-Storage Method 8-stage, fifth order low-storage scheme, optimized for compressible Navier–Stokes equations.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
@article{kennedy2000low, title={Low-storage, explicit Runge–Kutta schemes for the compressible Navier–Stokes equations}, author={Kennedy, Christopher A and Carpenter, Mark H and Lewis, R Michael}, journal={Applied numerical mathematics}, volume={35}, number={3}, pages={177–219}, year={2000}, publisher={Elsevier}}
OrdinaryDiffEqLowStorageRK.CKLLSRK85_4M_3R
— TypeCKLLSRK85_4M_3R(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Low-Storage Method 8-stage, fifth order low-storage scheme, optimized for compressible Navier–Stokes equations.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
@article{kennedy2000low, title={Low-storage, explicit Runge–Kutta schemes for the compressible Navier–Stokes equations}, author={Kennedy, Christopher A and Carpenter, Mark H and Lewis, R Michael}, journal={Applied numerical mathematics}, volume={35}, number={3}, pages={177–219}, year={2000}, publisher={Elsevier}}
OrdinaryDiffEqLowStorageRK.CKLLSRK85_4P_3R
— TypeCKLLSRK85_4P_3R(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Low-Storage Method 8-stage, fifth order low-storage scheme, optimized for compressible Navier–Stokes equations.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
@article{kennedy2000low, title={Low-storage, explicit Runge–Kutta schemes for the compressible Navier–Stokes equations}, author={Kennedy, Christopher A and Carpenter, Mark H and Lewis, R Michael}, journal={Applied numerical mathematics}, volume={35}, number={3}, pages={177–219}, year={2000}, publisher={Elsevier}}
OrdinaryDiffEqLowStorageRK.CKLLSRK54_3N_4R
— TypeCKLLSRK54_3N_4R(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Low-Storage Method 5-stage, fourth order low-storage scheme, optimized for compressible Navier–Stokes equations.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
@article{kennedy2000low, title={Low-storage, explicit Runge–Kutta schemes for the compressible Navier–Stokes equations}, author={Kennedy, Christopher A and Carpenter, Mark H and Lewis, R Michael}, journal={Applied numerical mathematics}, volume={35}, number={3}, pages={177–219}, year={2000}, publisher={Elsevier}}
OrdinaryDiffEqLowStorageRK.CKLLSRK54_3M_4R
— TypeCKLLSRK54_3M_4R(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Low-Storage Method 5-stage, fourth order low-storage scheme, optimized for compressible Navier–Stokes equations.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
@article{kennedy2000low, title={Low-storage, explicit Runge–Kutta schemes for the compressible Navier–Stokes equations}, author={Kennedy, Christopher A and Carpenter, Mark H and Lewis, R Michael}, journal={Applied numerical mathematics}, volume={35}, number={3}, pages={177–219}, year={2000}, publisher={Elsevier}}
OrdinaryDiffEqLowStorageRK.CKLLSRK65_4M_4R
— TypeCKLLSRK65_4M_4R(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. 6-stage, fifth order low-storage scheme, optimized for compressible Navier–Stokes equations.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
@article{kennedy2000low, title={Low-storage, explicit Runge–Kutta schemes for the compressible Navier–Stokes equations}, author={Kennedy, Christopher A and Carpenter, Mark H and Lewis, R Michael}, journal={Applied numerical mathematics}, volume={35}, number={3}, pages={177–219}, year={2000}, publisher={Elsevier}}
OrdinaryDiffEqLowStorageRK.CKLLSRK85_4FM_4R
— TypeCKLLSRK85_4FM_4R(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Low-Storage Method 8-stage, fifth order low-storage scheme, optimized for compressible Navier–Stokes equations.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
@article{kennedy2000low, title={Low-storage, explicit Runge–Kutta schemes for the compressible Navier–Stokes equations}, author={Kennedy, Christopher A and Carpenter, Mark H and Lewis, R Michael}, journal={Applied numerical mathematics}, volume={35}, number={3}, pages={177–219}, year={2000}, publisher={Elsevier}}
OrdinaryDiffEqLowStorageRK.CKLLSRK75_4M_5R
— TypeCKLLSRK75_4M_5R(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. CKLLSRK754M5R: Low-Storage Method 7-stage, fifth order low-storage scheme, optimized for compressible Navier–Stokes equations.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
@article{kennedy2000low, title={Low-storage, explicit Runge–Kutta schemes for the compressible Navier–Stokes equations}, author={Kennedy, Christopher A and Carpenter, Mark H and Lewis, R Michael}, journal={Applied numerical mathematics}, volume={35}, number={3}, pages={177–219}, year={2000}, publisher={Elsevier}}
OrdinaryDiffEqLowStorageRK.RDPK3Sp35
— TypeRDPK3Sp35(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. A third-order, five-stage method with embedded error estimator designed for spectral element discretizations of compressible fluid mechanics.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) 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
OrdinaryDiffEqLowStorageRK.RDPK3SpFSAL35
— TypeRDPK3SpFSAL35(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. A third-order, five-stage method with embedded error estimator using the FSAL property designed for spectral element discretizations of compressible fluid mechanics.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) 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
OrdinaryDiffEqLowStorageRK.RDPK3Sp49
— TypeRDPK3Sp49(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. A fourth-order, nine-stage method with embedded error estimator designed for spectral element discretizations of compressible fluid mechanics.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) 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
OrdinaryDiffEqLowStorageRK.RDPK3SpFSAL49
— TypeRDPK3SpFSAL49(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. A fourth-order, nine-stage method with embedded error estimator using the FSAL property designed for spectral element discretizations of compressible fluid mechanics.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) 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
OrdinaryDiffEqLowStorageRK.RDPK3Sp510
— TypeRDPK3Sp510(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. A fifth-order, ten-stage method with embedded error estimator designed for spectral element discretizations of compressible fluid mechanics.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) 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
OrdinaryDiffEqLowStorageRK.RDPK3SpFSAL510
— TypeRDPK3SpFSAL510(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. A fifth-order, ten-stage method with embedded error estimator using the FSAL property designed for spectral element discretizations of compressible fluid mechanics.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) 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
OrdinaryDiffEqLowStorageRK.HSLDDRK64
— TypeHSLDDRK64(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False(),
williamson_condition = true)
Explicit Runge-Kutta Method. Low-Storage Method 6-stage, fourth order low-stage, low-dissipation, low-dispersion scheme. Fixed timestep only.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.williamson_condition
: allows for an optimization that allows fusing broadcast expressions with the function callf
. However, it only works forArray
types.
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 }
OrdinaryDiffEqLowStorageRK.NDBLSRK134
— TypeNDBLSRK134(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False(),
williamson_condition = true)
Explicit Runge-Kutta Method. 13-stage, fourth order low-storage method with optimized stability regions for advection-dominated problems. Fixed timestep only.
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.williamson_condition
: allows for an optimization that allows fusing broadcast expressions with the function callf
. However, it only works forArray
types.
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
OrdinaryDiffEqLowStorageRK.SHLDDRK_2N
— TypeSHLDDRK_2N(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Low dissipation and dispersion Runge-Kutta schemes for computational acoustics
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
@article{stanescu19982n, title={2N-storage low dissipation and dispersion Runge-Kutta schemes for computational acoustics}, author={Stanescu, D and Habashi, WG}, journal={Journal of Computational Physics}, volume={143}, number={2}, pages={674–681}, year={1998}, publisher={Elsevier}}
OrdinaryDiffEqLowStorageRK.SHLDDRK52
— TypeSHLDDRK52(; stage_limiter! = OrdinaryDiffEq.trivial_limiter!,
step_limiter! = OrdinaryDiffEq.trivial_limiter!,
thread = OrdinaryDiffEq.False())
Explicit Runge-Kutta Method. Low dissipation and dispersion Runge-Kutta schemes for computational acoustics
Keyword Arguments
stage_limiter!
: function of the formlimiter!(u, integrator, p, t)
step_limiter!
: function of the formlimiter!(u, integrator, p, t)
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
@article{stanescu19982n, title={2N-storage low dissipation and dispersion Runge-Kutta schemes for computational acoustics}, author={Stanescu, D and Habashi, WG}, journal={Journal of Computational Physics}, volume={143}, number={2}, pages={674–681}, year={1998}, publisher={Elsevier}}