SummationByPartsOperators.jl: A Julia library of provably stable discretization techniques with mimetic properties
The Julia library SummationByPartsOperators.jl provides a unified interface of different discretization approaches including finite difference, Fourier pseudospectral, continuous Galerkin, and discontinuous Galerkin methods. This unified interface is based on the notion of summation-by-parts (SBP) operators. Originally developed for finite difference methods, SBP operators are discrete derivative operators designed specifically to get provably stable (semi-) discretizations, mimicking energy/entropy estimates from the continuous level discretely and paying special attention to boundary conditions.
SummationByPartsOperators.jl is mainly written to be useful for both students learning the basic concepts and researchers developing new numerical algorithms based on SBP operators. Thus, this package uses Julia's multiple dispatch and strong type system to provide a unified framework of all of these seemingly different discretizations while being reasonably optimized at the same time, achieving good performance without sacrificing flexibility.
SummationByPartsOperators.jl is a registered Julia package. Thus, you can install it from the Julia REPL via
julia> using Pkg; Pkg.add("SummationByPartsOperators")
If you want to update SummationByPartsOperators.jl, you can use
julia> using Pkg; Pkg.update("SummationByPartsOperators")
As usual, if you want to update SummationByPartsOperators.jl and all other packages in your current project, you can execute
julia> using Pkg; Pkg.update()
A brief list of notable changes is available in NEWS.md
.
Compute the derivative on a periodic domain using a central finite difference operator.
julia> using SummationByPartsOperators
julia> using Plots: plot, plot!
julia> D = periodic_derivative_operator(derivative_order=1, accuracy_order=2,
xmin=0.0, xmax=2.0, N=20)
Periodic first-derivative operator of order 2 on a grid in [0.0, 2.0] using 20 nodes,
stencils with 1 nodes to the left, 1 nodes to the right, and coefficients of Fornberg (1998)
Calculation of Weights in Finite Difference Formulas.
SIAM Rev. 40.3, pp. 685-691.
julia> x = grid(D); u = sinpi.(x);
julia> plot(x, D * u, label="numerical")
julia> plot!(x, π .* cospi.(x), label="analytical")
You should see a plot like the following.
Compute the derivative on a bounded domain using an SBP finite difference operator.
julia> using SummationByPartsOperators
julia> using Plots: plot, plot!
julia> D = derivative_operator(MattssonNordström2004(), derivative_order=1, accuracy_order=2,
xmin=0.0, xmax=1.0, N=21)
SBP first-derivative operator of order 2 on a grid in [0.0, 1.0] using 21 nodes
and coefficients of Mattsson, Nordström (2004)
Summation by parts operators for finite difference approximations of second
derivatives.
Journal of Computational Physics 199, pp. 503-540.
julia> x = grid(D); u = exp.(x);
julia> plot(x, D * u, label="numerical")
julia> plot!(x, exp.(x), label="analytical")
You should see a plot like the following.
The following derivative operators are implemented as "lazy"/matrix-free
operators, i.e. no large (size of the computational grid) matrix is formed
explicitly. They are linear operators and implement the same interface as
matrices in Julia (at least partially). In particular, *
and mul!
are
supported.
-
periodic_derivative_operator(; derivative_order, accuracy_order, xmin, xmax, N)
These are classical central finite difference operators using
N
nodes on the interval[xmin, xmax]
. -
periodic_derivative_operator(Holoborodko2008(); derivative_order, accuracy_order, xmin, xmax, N)
These are central finite difference operators using
N
nodes on the interval[xmin, xmax]
and the coefficients of Pavel Holoborodko. -
fourier_derivative_operator(; xmin, xmax, N)
Fourier derivative operators are implemented using the fast Fourier transform of FFTW.jl.
All of these periodic derivative operators support multiplication and addition
such that polynomials and rational functions of them can be represented efficiently,
e.g. to solve elliptic problems of the form u = (D^2 + I) \ f
.
-
derivative_operator(source_of_coefficients; derivative_order, accuracy_order, xmin, xmax, N)
Finite difference SBP operators for first and second derivatives can be obtained by using
MattssonNordström2004()
assource_of_coefficients
. Other sources of coefficients are implemented as well. To obtain a full list of all operators, usesubtypes(SourceOfCoefficients)
. -
legendre_derivative_operator(; xmin, xmax, N)
Use Lobatto Legendre polynomial collocation schemes on
N
, i.e. polynomials of degreeN-1
, implemented via PolynomialBases.jl.
Additionally, some artificial dissipation/viscosity operators are implemented.
The most basic usage is Di = dissipation_operator(D)
,
where D
can be a (periodic, Fourier, Legendre, SBP FD) derivative
operator. Use ?dissipation_operator
for more details.
SBP operators on bounded domains can be coupled continuously or discontinuously
to obtain CG//DG-type methods. You need to create an appropriate mesh
and
a basic operator D
that should be used on each element.
Then, global CG/DG operators are constructed lazily/matrix-free by calling
couple_continuously(D, mesh)
or
couple_discontinuously(D, mesh, coupling::Union{Val{:plus}, Val{:central}, Val{:minus}}=Val(:central))
.
Choosing coupling=Val(:central)
yields a classical SBP operator; the other two
coupling
types result in upwind SBP operators. Currently, only uniform meshes
UniformMesh1D(xmin::Real, xmax::Real, Nx::Integer)
UniformPeriodicMesh1D(xmin::Real, xmax::Real, Nx::Integer)
are implemented.
Sometimes, it can be convenient to obtain an explicit (sparse, banded) matrix form of the operators. Therefore, some conversion functions are supplied, e.g.
julia> using SummationByPartsOperators
julia> D = derivative_operator(MattssonNordström2004(),
derivative_order=1, accuracy_order=2,
xmin=0.0, xmax=1.0, N=5)
SBP first-derivative operator of order 2 on a grid in [0.0, 1.0] using 5 nodes
and coefficients of Mattsson, Nordström (2004)
Summation by parts operators for finite difference approximations of second
derivatives.
Journal of Computational Physics 199, pp. 503-540.
julia> Matrix(D)
5×5 Array{Float64,2}:
-4.0 4.0 0.0 0.0 0.0
-2.0 0.0 2.0 0.0 0.0
0.0 -2.0 0.0 2.0 0.0
0.0 0.0 -2.0 0.0 2.0
0.0 0.0 0.0 -4.0 4.0
julia> using SparseArrays
julia> sparse(D)
5×5 SparseMatrixCSC{Float64, Int64} with 10 stored entries:
-4.0 4.0 ⋅ ⋅ ⋅
-2.0 ⋅ 2.0 ⋅ ⋅
⋅ -2.0 ⋅ 2.0 ⋅
⋅ ⋅ -2.0 ⋅ 2.0
⋅ ⋅ ⋅ -4.0 4.0
julia> using BandedMatrices
julia> BandedMatrix(D)
5×5 BandedMatrix{Float64,Array{Float64,2},Base.OneTo{Int64}}:
-4.0 4.0 ⋅ ⋅ ⋅
-2.0 0.0 2.0 ⋅ ⋅
⋅ -2.0 0.0 2.0 ⋅
⋅ ⋅ -2.0 0.0 2.0
⋅ ⋅ ⋅ -4.0 4.0
The latest documentation is available
online
and under docs/src
.
Some additional examples can be found in the directory
notebooks
.
In particular, examples of complete discretizations of
the linear advection equation,
the heat equation,
and the wave equation are available.
Further examples are supplied as
tests.
If you use SummationByPartsOperators.jl for your research, please cite it using the bibtex entry
@article{ranocha2021sbp,
title={{SummationByPartsOperators.jl}: {A} {J}ulia library of provably stable
semidiscretization techniques with mimetic properties},
author={Ranocha, Hendrik},
journal={Journal of Open Source Software},
year={2021},
month={08},
doi={10.21105/joss.03454},
volume={6},
number={64},
pages={3454},
publisher={The Open Journal},
url={https://github.com/ranocha/SummationByPartsOperators.jl}
}
This project is licensed under the MIT license (see LICENSE.md). Since it is an open-source project, we are very happy to accept contributions from the community. Please refer to CONTRIBUTING.md for more details.