SummationByPartsOperators.jl API
SummationByPartsOperators.SummationByPartsOperators
— ModuleSummationByPartsOperators
SummationByPartsOperators.jl is a Julia library of summation-by-parts (SBP) operators, which are discrete derivative operators developed to get provably stable semidiscretizations, paying special attention to boundary conditions. Discretizations included in this framework are finite difference, Fourier pseudospectral, continuous Galerkin, and discontinuous Galerkin methods. The main aim of SummationByPartsOperators.jl is to be useful for researchers and students to learn the basic concepts by providing a unified framework of all of these seemingly different discretizations. At the same time, the implementation is optimized to achieve good performance without sacrificing flexibility.
Check out the documentation for further information. Some noticeable functions to start with are derivative_operator
, legendre_derivative_operator
, periodic_derivative_operator
, fourier_derivative_operator
, dissipation_operator
, and grid
.
If you use this package for your research, please cite it using
@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}
}
SummationByPartsOperators.BeljaddLeFlochMishraParés2017
— TypeBeljaddLeFlochMishraParés2017()
Coefficients of the periodic operators given in
- Beljadid, LeFloch, Mishra, Parés (2017) Schemes with Well-Controlled Dissipation. Hyperbolic Systems in Nonconservative Form. Communications in Computational Physics 21.4, pp. 913-946.
SummationByPartsOperators.BurgersNonperiodicSemidiscretization
— TypeBurgersNonperiodicSemidiscretization(D, Di, split_form, left_bc, right_bc)
A semidiscretization of Burgers' equation $\partial_t u(t,x) + \partial_x \frac{u(t,x)^2}{2} = 0$ with boundary conditions left_bc(t)
, right_bc(t)
.
D
is a first-derivative SBP operator, Di
an associated dissipation operator or nothing
, and split_form::Union{Val(true), Val(false)}
determines whether the canonical split form or the conservative form is used.
SummationByPartsOperators.BurgersPeriodicSemidiscretization
— TypeBurgersPeriodicSemidiscretization(D, Di, split_form)
A semidiscretization of Burgers' equation $\partial_t u(t,x) + \partial_x \frac{u(t,x)^2}{2} = 0$ with periodic boundary conditions.
D
is a first-derivative SBP operator, Di
an associated dissipation operator or nothing
, and split_form::Union{Val(true), Val(false)}
determines whether the canonical split form or the conservative form is used.
SummationByPartsOperators.ConstantFilter
— TypeConstantFilter
Represents the action of a modal filter on values in a nodal basis with fixed strength.
SummationByPartsOperators.ConstantFilter
— MethodConstantFilter(D::FourierDerivativeOperator, filter)
Create a modal filter with constant parameters adapted to the Fourier derivative operator D
with parameters given by the filter function filter
.
SummationByPartsOperators.ConstantFilter
— MethodConstantFilter(D::LegendreDerivativeOperator, filter, TmpEltype=T)
Create a modal filter with constant parameters adapted to the Legendre derivative operator D
with parameters given by the filter function filter
.
SummationByPartsOperators.CubicNonperiodicSemidiscretization
— TypeCubicNonperiodicSemidiscretization(D, Di, split_form, left_bc, right_bc)
A semidiscretization of the cubic conservation law $\partial_t u(t,x) + \partial_x u(t,x)^3 = 0$ with nonperiodic boundary conditions left_bc(t)
, right_bc(t)
.
D
is a first-derivative SBP operator, Di
an associated dissipation operator or nothing
, and split_form::Union{Val(true), Val(false)}
determines whether the canonical split form or the conservative form is used.
SummationByPartsOperators.CubicPeriodicSemidiscretization
— TypeCubicPeriodicSemidiscretization(D, Di, split_form)
A semidiscretization of the cubic conservation law $\partial_t u(t,x) + \partial_x u(t,x)^3 = 0$ with periodic boundary conditions.
D
is a first-derivative SBP operator, Di
an associated dissipation operator or nothing
, and split_form::Union{Val(true), Val(false)}
determines whether the canonical split form or the conservative form is used.
SummationByPartsOperators.DerivativeCoefficientRow
— TypeDerivativeCoefficientRow{T,Start,Length}
A struct representing a row in the boundary block of an SBP derivative operator with scalar type T
.
SummationByPartsOperators.DerivativeCoefficients
— TypeDerivativeCoefficients
The coefficients of a derivative operator on a nonperiodic grid.
SummationByPartsOperators.DerivativeOperator
— TypeDerivativeOperator
A derivative operator on a nonperiodic finite difference grid. See derivative_operator
.
SummationByPartsOperators.DienerDorbandSchnetterTiglio2007
— TypeDienerDorbandSchnetterTiglio2007()
Coefficients of the SBP operators given in
- Diener, Dorband, Schnetter, Tiglio (2007) Optimized high-order derivative and dissipation operators satisfying summation by parts, and applications in three-dimensional multi-block evolutions. Journal of Scientific Computing 32.1, pp. 109-145.
See also (second- and fourth-order operators)
- Mattsson, Nordström (2004) Summation by parts operators for finite difference approximations of second derivatives. Journal of Computational Physics 199, pp. 503-540.
The dissipation operators proposed by Diener, Dorband, Schnetter, Tiglio (2007) for the diagonal-norm operators are the same as the ones of
- Mattsson, Svärd, Nordström (2004) Stable and Accurate Artificial Dissipation. Journal of Scientific Computing 21.1, pp. 57-79.
SummationByPartsOperators.DissipationOperator
— TypeDissipationOperator
A dissipation operator on a nonperiodic finite difference grid. See dissipation_operator
.
SummationByPartsOperators.ExponentialFilter
— TypeExponentialFilter
Represents the exponential filter function σ(η) = exp(-α*η^p)
.
SummationByPartsOperators.FactorisationWrapper
— TypeFactorisationWrapper
A small wrapper around a a factorisation fact
, allowing to represent multiplication by the inverse of fact
.
SummationByPartsOperators.FastMode
— TypeFastMode()
A (probably) faster execution mode that might depend on packages such as LoopVectorization.jl.
SummationByPartsOperators.Fornberg1998
— TypeFornberg1998()
Coefficients of the periodic operators given in
- Fornberg (1998) Calculation of Weights in Finite Difference Formulas. SIAM Rev. 40.3, pp. 685-691.
SummationByPartsOperators.FourierConstantViscosity
— TypeFourierConstantViscosity
Fourier viscosity operator with constant coefficients for the periodic 1st derivative Fourier operator.
SummationByPartsOperators.FourierDerivativeOperator
— TypeFourierDerivativeOperator{T}
A derivative operator on a periodic grid with real scalar type T
computing the first derivative using a spectral Fourier expansion via real discrete Fourier transforms.
See also fourier_derivative_operator
.
SummationByPartsOperators.FourierDerivativeOperator
— MethodFourierDerivativeOperator(xmin::T, xmax::T, N::Integer) where {T<:Real}
Construct the FourierDerivativeOperator
on a uniform grid between xmin
and xmax
using N
nodes and N÷2+1
complex Fourier modes.
See also fourier_derivative_operator
.
SummationByPartsOperators.FourierDerivativeOperator2D
— TypeFourierDerivativeOperator2D{T<:Real}
A derivative operator on a two-dimensional periodic grid with scalar type T
computing the first derivatives using a spectral Fourier expansion via real discrete Fourier transforms.
SummationByPartsOperators.FourierDerivativeOperator2D
— MethodFourierDerivativeOperator2D(xmin, xmax, Nx, ymin, ymax, Ny)
Construct the FourierDerivativeOperator
on a uniform grid between xmin
and xmax
using Nx
nodes and ymin
and ymax
using Ny
nodes.
SummationByPartsOperators.GlaubitzNordströmÖffner2023
— TypeGlaubitzNordströmÖffner2023()
Function space SBP (FSBP) operators given in
- Glaubitz, Nordström, Öffner (2023) Summation-by-parts operators for general function spaces. SIAM Journal on Numerical Analysis 61, 2, pp. 733-754.
See also
- Glaubitz, Nordström, Öffner (2024) An optimization-based construction procedure for function space based summation-by-parts operators on arbitrary grids. arXiv, arXiv:2405.08770v1.
SummationByPartsOperators.Holoborodko2008
— TypeHoloborodko2008()
Coefficients of the periodic operators given in
- Holoborodko (2008) Smooth Noise Robust Differentiators. http://www.holoborodko.com/pavel/numerical-methods/numerical-derivative/smooth-low-noise-differentiators/
SummationByPartsOperators.LanczosLowNoise
— TypeLanczosLowNoise()
Coefficients of the periodic operators given in
- Holoborodko (2008) Smooth Noise Robust Differentiators. http://www.holoborodko.com/pavel/numerical-methods/numerical-derivative/lanczos-low-noise-differentiators/
SummationByPartsOperators.LegendreDerivativeOperator
— TypeLegendreDerivativeOperator{T<:Real}
A derivative operator on a nonperiodic Lobatto-Legendre grid with scalar type T
computing the first derivative using a Legendre expansion.
SummationByPartsOperators.LegendreDerivativeOperator
— MethodLegendreDerivativeOperator(xmin::T, xmax::T, N::Int) where {T<:Real}
Construct the LegendreDerivativeOperator
on a uniform grid between xmin
and xmax
using N
nodes and N-1
Legendre modes.
SummationByPartsOperators.LegendreSecondDerivativeOperator
— TypeLegendreSecondDerivativeOperator{T<:Real}
A derivative operator on a nonperiodic Lobatto-Legendre grid with scalar type T
computing the second derivative using a Legendre expansion.
SummationByPartsOperators.LinearlyCombinedDerivativeOperators
— TypeLinearlyCombinedDerivativeOperators
Form linear combinations of several derivative operators lazily.
SummationByPartsOperators.MadayTadmor1989
— TypeMadayTadmor1989()
Coefficients of the Fourier spectral viscosity given in
- Maday, Tadmor (1989) Analysis of the Spectral Vanishing Viscosity Method for Periodic Conservation Laws. SIAM Journal on Numerical Analysis 26.4, pp. 854-870.
SummationByPartsOperators.MatrixDerivativeOperator
— TypeMatrixDerivativeOperator{T <: Real}
MatrixDerivativeOperator(xmin, xmax, nodes, weights, D, accuracy_order, source)
A derivative operator on a nonperiodic grid with scalar type T
computing a derivative as matrix vector product. This type is designed to make it easy to experiment with new operators given in matrix form.
An instance of this type can be constructed by passing the endpoints xmin
, xmax
of the desired grid as well as the nodes
, weights
, and the derivative operator D::AbstractMatrix
on a reference interval, assuming that the nodes
contain the boundary points of the reference interval. source
is the source of coefficients and can be nothing
for experimentation.
SummationByPartsOperators.Mattsson2012
— TypeMattsson2012()
Coefficients of the SBP operators given in
- Mattsson (2012) Summation by Parts Operators for Finite Difference Approximations of Second-Derivatives with Variable Coefficients. Journal of Scientific Computing 51, pp. 650-682.
SummationByPartsOperators.Mattsson2014
— TypeMattsson2014()
Coefficients of the SBP operators given in
- Mattsson (2014) Diagonal-norm summation by parts operators for finite difference approximations of third and fourth derivatives. Journal of Computational Physics 274, pp. 432-454.
SummationByPartsOperators.Mattsson2017
— TypeMattsson2017(version::Symbol)
Coefficients of the upwind SBP operators given in
- Mattsson (2017) Diagonal-norm upwind SBP operators. Journal of Computational Physics 335, pp. 283-310.
You can choose between the different versions :central
, :plus
, and :minus
.
SummationByPartsOperators.MattssonAlmquistCarpenter2014Extended
— TypeMattssonAlmquistCarpenter2014Extended()
Coefficients of the extended SBP operators given in
- Mattsson, Almquist, Carpenter (2014) Optimal diagonal-norm SBP operators. Journal of Computational Physics 264, pp. 91-111.
SummationByPartsOperators.MattssonAlmquistCarpenter2014Optimal
— TypeMattssonAlmquistCarpenter2014Optimal()
Coefficients of the optimal SBP operators with nonuniform grid given in
- Mattsson, Almquist, Carpenter (2014) Optimal diagonal-norm SBP operators. Journal of Computational Physics 264, pp. 91-111.
SummationByPartsOperators.MattssonAlmquistVanDerWeide2018Accurate
— TypeMattssonAlmquistVanDerWeide2018Accurate()
Coefficients of the optimized SBP operators with nonuniform grid given in
- Mattsson, Almquist, van der Weide (2018) Boundary optimized diagonal-norm SBP operators. Journal of Computational Physics 374, pp. 1261-1266.
SummationByPartsOperators.MattssonAlmquistVanDerWeide2018Minimal
— TypeMattssonAlmquistVanDerWeide2018Minimal()
Coefficients of the optimized SBP operators with nonuniform grid given in
- Mattsson, Almquist, van der Weide (2018) Boundary optimized diagonal-norm SBP operators. Journal of Computational Physics 374, pp. 1261-1266.
SummationByPartsOperators.MattssonNordström2004
— TypeMattssonNordström2004()
Coefficients of the SBP operators given in
- Mattsson, Nordström (2004) Summation by parts operators for finite difference approximations of second derivatives. Journal of Computational Physics 199, pp. 503-540.
SummationByPartsOperators.MattssonSvärdNordström2004
— TypeMattssonSvärdNordström2004()
Coefficients of the SBP operators given in
- Mattsson, Svärd, Nordström (2004) Stable and Accurate Artificial Dissipation. Journal of Scientific Computing 21.1, pp. 57-79.
SummationByPartsOperators.MattssonSvärdShoeybi2008
— TypeMattssonSvärdShoeybi2008()
Coefficients of the SBP operators given in
- Mattsson, Svärd, Shoeybi (2008) Stable and accurate schemes for the compressible Navier-Stokes equations. Journal of Computational Physics 227, pp. 2293-2316.
SummationByPartsOperators.PeriodicDerivativeCoefficients
— TypePeriodicDerivativeCoefficients
The coefficients of a derivative operator on a periodic grid.
SummationByPartsOperators.PeriodicDerivativeOperator
— TypePeriodicDerivativeOperator
A derivative operator on a uniform periodic grid. See periodic_derivative_operator
and periodic_central_derivative_operator
.
SummationByPartsOperators.PeriodicDissipationOperator
— TypePeriodicDissipationOperator
A dissipation operator on a periodic finite difference grid. See dissipation_operator
.
SummationByPartsOperators.PeriodicUpwindOperators
— TypePeriodicUpwindOperators
PeriodicUpwindOperators(D_minus, D_central, D_plus)
A struct bundling the individual operators available for periodic upwind SBP operators. The individual operators are available as D.minus
, D.plus
(and optionally D.central
, if provided), where D::PeriodicUpwindOperators
.
The combined struct behaves as much as possible as an operator itself as long as no ambiguities arise. For example, upwind operators need to use the same grid and mass matrix, so mass_matrix
, grid
, xmin
, xmax
etc. are available but mul!
is not.
It is recommended to construct an instance of PeriodicUpwindOperators
using upwind_operators
. An instance can also be constructed manually by passing the operators in the order D_minus, D_central, D_plus
.
See also upwind_operators
, UpwindOperators
SummationByPartsOperators.QuarticNonconvexPeriodicSemidiscretization
— TypeQuarticNonconvexPeriodicSemidiscretization(D, Di, split_form)
A semidiscretization of the quartic nonconvex conservation law $\partial_t u(t,x) + \partial_x ( u(t,x)^4 - 10 u(t,x)^2 + 3 u(t,x) ) = 0$ with periodic boundary conditions.
D
is a first-derivative SBP operator, Di
an associated dissipation operator or nothing
, and split_form::Union{Val(true), Val(false)}
determines whether the canonical split form or the conservative form is used.
SummationByPartsOperators.SafeMode
— TypeSafeMode()
A safe execution mode relying only on basic functionality of Julia.
SummationByPartsOperators.SharanBradyLivescu2022
— TypeSharanBradyLivescu2022(alpha_left, alpha_right)
Coefficients of the cut-cell SBP operators given in
- Sharan, Brady, Livescu (2022) High-order dimensionally-split Cartesian embedded boundary method for non-dissipative schemes. Journal of Computational Physics 464, 111341.
Here, alpha_left
* Δx is the spacing between the left endpoint and the second node, while alpha_right
* Δx is the spacing between the right endpoint and the last node.
SummationByPartsOperators.SourceOfCoefficients
— TypeSourceOfCoefficients
All sources of coefficients (articles) are subtypes of this abstract type.
SummationByPartsOperators.Tadmor1989
— TypeTadmor1989()
Coefficients of the Fourier spectral viscosity given in
- Tadmor (1989) Convergence of Spectral Methods for Nonlinear Conservation Laws. SIAM Journal on Numerical Analysis 26.1, pp. 30-44.
SummationByPartsOperators.Tadmor1993
— TypeTadmor1993()
Coefficients of the Fourier super spectral viscosity given in
- Tadmor (1993) Super Viscosity and Spectral Approximations of Nonlinear Conservation Laws. Numerical Methods for Fluid Dynamics IV, pp. 69-82.
SummationByPartsOperators.TadmorWaagan2012Convergent
— TypeTadmorWaagan2012Convergent()
Coefficients of the Fourier spectral viscosity given in
- Tadmor, Waagan (2012) Adaptive Spectral Viscosity for Hyperbolic Conservation Laws. SIAM Journal on Scientific Computing 34.2, pp. A993-A1009.
See also
- Schochet (1990) The Rate of Convergence of Spectral-Viscosity Methods for Periodic Scalar Conservation Laws. SIAM Journal on Numerical Analysis 27.5, pp. 1142-1159.
SummationByPartsOperators.TadmorWaagan2012Standard
— TypeTadmorWaagan2012Standard()
Coefficients of the Fourier spectral viscosity given in
- Tadmor, Waagan (2012) Adaptive Spectral Viscosity for Hyperbolic Conservation Laws. SIAM Journal on Scientific Computing 34.2, pp. A993-A1009.
SummationByPartsOperators.ThreadedMode
— TypeThreadedMode()
An execution mode using multiple threads and possibly further optimizations, cf. FastMode
.
SummationByPartsOperators.UniformMesh1D
— TypeUniformMesh1D(xmin::Real, xmax::Real, Nx::Integer)
UniformMesh1D(; xmin::Real, xmax::Real, Nx::Integer)
A uniform mesh in one space dimension of Nx
cells between xmin
and xmax
.
SummationByPartsOperators.UniformPeriodicMesh1D
— TypeUniformPeriodicMesh1D(xmin::Real, xmax::Real, Nx::Integer)
UniformPeriodicMesh1D(; xmin::Real, xmax::Real, Nx::Integer)
A uniform periodic mesh in one space dimension of Nx
cells between xmin
and xmax
.
SummationByPartsOperators.UpwindOperators
— TypeUpwindOperators
UpwindOperators(D_minus, D_central, D_plus)
A struct bundling the individual operators available for non-periodic upwind SBP operators. The individual operators are available as D.minus
, D.plus
(and optionally D.central
, if provided), where D::UpwindOperators
.
The combined struct behaves as much as possible as an operator itself as long as no ambiguities arise. For example, upwind operators need to use the same grid and mass matrix, so mass_matrix
, grid
, xmin
, xmax
etc. are available but mul!
is not.
It is recommended to construct an instance of UpwindOperators
using upwind_operators
. An instance can also be constructed manually by passing the operators in the order D_minus, D_central, D_plus
.
See also upwind_operators
, PeriodicUpwindOperators
SummationByPartsOperators.VarCoefDerivativeCoefficients
— TypeVarCoefDerivativeCoefficients
The coefficients of a variable coefficient derivative operator on a nonperiodic grid.
SummationByPartsOperators.VarCoefDerivativeOperator
— TypeVarCoefDerivativeOperator
A dissipation operator on a nonperiodic finite difference grid.
SummationByPartsOperators.VariableLinearAdvectionNonperiodicSemidiscretization
— TypeVariableLinearAdvectionNonperiodicSemidiscretization(D, Di, a, split_form,
left_bc, right_bc)
A semidiscretization of the linear advection equation $\partial_t u(t,x) + \partial_x ( a(x) u(t,x) ) = 0$ with boundary conditions left_bc(t)
, right_bc(t)
.
D
is an SBP derivative operator, Di
an associated dissipation operator or nothing
, a(x)
the variable coefficient, and split_form::Union{Val(false), Val(true)}
determines whether the canonical split form or the conservative form should be used.
SummationByPartsOperators.VariableLinearAdvectionPeriodicSemidiscretization
— TypeVariableLinearAdvectionPeriodicSemidiscretization(D, Di, a, split_form)
A semidiscretization of the linear advection equation $\partial_t u(t,x) + \partial_x ( a(x) u(t,x) ) = 0$ with periodic boundary conditions.
D
is a periodic SBP derivative operator, Di
an associated dissipation operator or nothing
, a(x)
the variable coefficient, and split_form::Union{Val(false), Val(true)}
determines whether the canonical split form or the conservative form should be used.
SummationByPartsOperators.WaveEquationNonperiodicSemidiscretization
— TypeWaveEquationNonperiodicSemidiscretization(D, left_bc, right_bc)
A semidiscretization of the linear wave equation $\partial_t^2 u(t,x) = \partial_x^2 u(t,x)$.
D
is assumed to be a second-derivative SBP operator and the boundary conditions can be Val(:HomogeneousNeumann)
, Val(:HomogeneousDirichlet)
, or Val(:NonReflecting)
.
SummationByPartsOperators.WilliamsDuru2024
— TypeWilliamsDuru2024(version::Symbol)
Coefficients of the upwind SBP operators given in
- Williams, Duru (2024) Full-Spectrum Dispersion Relation Preserving Summation-By-Parts Operators. SIAM Journal on Numerical Analysis 62.4, pp. 1565-1588.
You can choose between the different versions :central
, :plus
, and :minus
.
LinearAlgebra.mul!
— Functionmul!(du, D::DerivativeOperator, u, α=true, β=false)
Efficient in-place version of du = α * D * u + β * du
. Note that du
must not be aliased with u
.
PolynomialBases.compute_coefficients!
— Methodcompute_coefficients!(uval, u, D::AbstractDerivativeOperator)
Compute the nodal values of the function u
at the grid associated to the derivative operator D
and stores the result in uval
.
PolynomialBases.compute_coefficients
— Methodcompute_coefficients(u, D::AbstractDerivativeOperator)
Compute the nodal values of the function u
at the grid associated to the derivative operator D
.
PolynomialBases.evaluate_coefficients!
— Methodevaluate_coefficients!(xplot, uplot, u, D::AbstractDerivativeOperator)
Evaluates the nodal coefficients u
at a grid associated to the derivative operator D
and stores the result in xplot, uplot
. Returns xplot, uplot
, where xplot
contains the nodes and uplot
the corresponding values of u
.
PolynomialBases.evaluate_coefficients
— Methodevaluate_coefficients(u, D::AbstractDerivativeOperator)
Evaluates the nodal coefficients u
at a grid associated to the derivative operator D
. Returns xplot, uplot
, where xplot
contains the nodes and uplot
the corresponding values of u
.
PolynomialBases.integrate
— Methodintegrate(func, u, D::DerivativeOperator)
Map the function func
to the coefficients u
and integrate with respect to the quadrature rule associated with the SBP derivative operator D
.
PolynomialBases.integrate
— Methodintegrate(func, u, D::PeriodicDerivativeOperator)
Map the function func
to the coefficients u
and integrate with respect to the quadrature rule associated with the periodic derivative operator D
.
PolynomialBases.integrate
— Methodintegrate(func, u, D::AbstractPeriodicDerivativeOperator)
Map the function func
to the coefficients u
and integrate with respect to the quadrature rule associated with the derivative operator D
.
SummationByPartsOperators.accuracy_order
— Functionaccuracy_order(D)
Return the order of accuracy of a derivative operator D
. For SBP finite difference operators, this refers to the interior order of accuracy.
SummationByPartsOperators.couple_continuously
— Functioncouple_continuously(D, mesh)
Return a derivative operator corresponding to a continuous coupling of D
on the cells of the given mesh
as in (nodal) continuous Galerkin (CG) methods. If the underlying SBP operators are LegendreDerivativeOperator
s, these are CG spectral element methods (CGSEM). However, a continuous coupling of arbitrary SBP operators is supported.
The mesh
can be a UniformMesh1D
or a UniformPeriodicMesh1D
.
References
- Ranocha, Mitsotakis, Ketcheson (2021). A Broad Class of Conservative Numerical Methods for Dispersive Wave Equations. DOI: 10.4208/cicp.OA-2020-0119
SummationByPartsOperators.couple_discontinuously
— Functioncouple_discontinuously(D, mesh, [coupling=Val(:central)])
Return a derivative operator corresponding to a discontinuous coupling of D
on the cells of the given mesh
as in (nodal) discontinuous Galerkin (CG) methods. If the underlying SBP operators are LegendreDerivativeOperator
s, these are DG spectral element methods (DGSEM). However, a discontinuous coupling of arbitrary SBP operators is supported.
The mesh
can be a UniformMesh1D
or a UniformPeriodicMesh1D
. The coupling
can be
Val(:central)
(default), resulting in classical SBP propertiesVal(:minus)
orVal(:plus)
, resulting in upwind SBP operators
References
- Ranocha, Mitsotakis, Ketcheson (2021). A Broad Class of Conservative Numerical Methods for Dispersive Wave Equations. DOI: 10.4208/cicp.OA-2020-0119
SummationByPartsOperators.derivative_left
— Methodderivative_left(D::AbstractNonperiodicDerivativeOperator, der_order)
Get a representation of the linear functional evaluation the N
th derivative at the left boundary node as (dense) vector.
SummationByPartsOperators.derivative_left
— Methodderivative_left(D::DerivativeOperator, u, der_order::Val{N})
Compute the N
-th derivative of the function given by the coefficients u
at the left boundary of the grid.
SummationByPartsOperators.derivative_operator
— Functionderivative_operator(source_of_coefficients,
derivative_order, accuracy_order,
xmin, xmax, N, mode=FastMode())
derivative_operator(source_of_coefficients;
derivative_order, accuracy_order,
xmin, xmax, N, mode=FastMode())
Create a DerivativeOperator
approximating the derivative_order
-th derivative on a grid between xmin
and xmax
with N
grid points up to order of accuracy accuracy_order
. with coefficients given by source_of_coefficients
. The evaluation of the derivative can be parallelized using threads by choosing mode=ThreadedMode()
.
SummationByPartsOperators.derivative_order
— Functionderivative_order(D)
Return the order of the derivative associated to the derivative operator D
. For example, it will return 1
for a first-derivative SBP operator.
SummationByPartsOperators.derivative_right
— Methodderivative_right(D::AbstractNonperiodicDerivativeOperator, der_order)
Get a representation of the linear functional evaluation the N
th derivative at the right boundary node as (dense) vector.
SummationByPartsOperators.derivative_right
— Methodderivative_right(D::DerivativeOperator, u, der_order::Val{N})
Compute the N
-th derivative of the function given by the coefficients u
at the right boundary of the grid.
SummationByPartsOperators.dissipation_operator
— Functiondissipation_operator(source_of_coefficients, order, xmin, xmax, N,
left_weights, right_weights, mode=FastMode())
Create a negative semidefinite DissipationOperator
using undivided differences approximating a weighted order
-th derivative on a grid between xmin
and xmax
with N
grid points up to order of accuracy 2 with coefficients given by source_of_coefficients
. The norm matrix is given by left_weights
and right_weights
. The evaluation of the derivative can be parallelized using threads by choosing mode=ThreadedMode()
.
SummationByPartsOperators.dissipation_operator
— Methoddissipation_operator(D::PeriodicDerivativeOperator;
strength=one(eltype(D)),
order=accuracy_order(D),
mode=D.coefficients.mode)
Create a negative semidefinite DissipationOperator
using undivided differences approximating a order
-th derivative with strength strength
adapted to the derivative operator D
. The evaluation of the derivative can be parallelized using threads by choosing mode=ThreadedMode()
.
SummationByPartsOperators.dissipation_operator
— Methoddissipation_operator([source_of_coefficients=MattssonSvärdNordström2004()],
D::DerivativeOperator{T};
strength=one(T),
order::Int=accuracy_order(D),
mode=D.coefficients.mode)
Create a negative semidefinite DissipationOperator
using undivided differences approximating a weighted order
-th derivative adapted to the derivative operator D
with coefficients given in source_of_coefficients
. The evaluation of the derivative can be parallelized using threads by choosing mode=ThreadedMode()
.
SummationByPartsOperators.fornberg
— Methodfornberg(x::Vector{T}, m::Int) where {T}
Calculate the weights of a finite difference approximation of the m
th derivative with maximal order of accuracy at 0
using the nodes x
, see Fornberg (1998) Calculation of Weights in Finite Difference Formulas SIAM Rev. 40.3, pp. 685-691.
SummationByPartsOperators.fourier_derivative_matrix
— Functionfourier_derivative_matrix(N, xmin::Real=0.0, xmax::Real=2π)
Compute the Fourier derivative matrix with respect to the corresponding nodal basis using N
nodes, see Kopriva (2009) Implementing Spectral Methods for PDEs, Algorithm 18.
SummationByPartsOperators.fourier_derivative_operator
— Methodfourier_derivative_operator(xmin::Real, xmax::Real, N::Integer)
fourier_derivative_operator(; xmin::Real, xmax::Real, N::Integer)
Construct the FourierDerivativeOperator
on a uniform grid between xmin
and xmax
using N
nodes and N÷2+1
complex Fourier modes.
SummationByPartsOperators.function_space_operator
— Functionfunction_space_operator(basis_functions, nodes, source;
derivative_order = 1, accuracy_order = 0,
opt_alg = Optim.LBFGS(), options = Optim.Options(g_tol = 1e-14, iterations = 10000),
verbose = false)
Construct an operator that represents a first-derivative operator in a function space spanned by the basis_functions
, which is an iterable of functions. The operator is constructed on the interval [x_min, x_max]
with the nodes nodes
, where x_min
is taken as the minimal value in nodes
and x_max
the maximal value. Note that the nodes
will be sorted internally. The accuracy_order
is the order of the accuracy of the operator, which can optionally be passed, but does not have any effect on the operator. The operator is constructed solving an optimization problem with Optim.jl. You can specify the optimization algorithm and options for the optimization problem with the keyword arguments opt_alg
and options
respectively, see also the documentation of Optim.jl
The operator that is returned follows the general interface. Currently, it is wrapped in a MatrixDerivativeOperator
, but this might change in the future. In order to use this function, the package Optim
must be loaded.
See also GlaubitzNordströmÖffner2023
.
This function requires at least Julia 1.9.
This is an experimental feature and may change in future releases.
SummationByPartsOperators.grid
— Functiongrid(D)
Return the grid associated to a derivative operator D
.
SummationByPartsOperators.left_boundary_weight
— Functionleft_boundary_weight(D)
Return the left-boundary weight of the (diagonal) mass matrix M
associated to the derivative operator D
.
SummationByPartsOperators.legendre_derivative_operator
— Methodlegendre_derivative_operator(xmin::Real, xmax::Real, N::Integer)
legendre_derivative_operator(; xmin::Real, xmax::Real, N::Integer)
Construct the LegendreDerivativeOperator
on a uniform grid between xmin
and xmax
using N
nodes and N-1
Legendre modes.
SummationByPartsOperators.legendre_second_derivative_operator
— Methodlegendre_second_derivative_operator(xmin::Real, xmax::Real, N::Integer)
legendre_second_derivative_operator(; xmin::Real, xmax::Real, N::Integer)
Construct the LegendreDerivativeOperator
on a uniform grid between xmin
and xmax
using N
nodes and N-1
Legendre modes.
SummationByPartsOperators.mass_matrix
— Methodmass_matrix(D::Union{DerivativeOperator,VarCoefDerivativeOperator})
Create the diagonal mass matrix for the SBP derivative operator D
.
SummationByPartsOperators.mul_transpose_derivative_left!
— Methodmul_transpose_derivative_left!(u, D::DerivativeOperator, der_order::Val{N}, α=true, β=false)
Set the grid function u
to α
times the transposed N
-th derivative functional applied to u
plus β
times u
in the domain of the N
-th derivative functional at the left boundary of the grid. Thus, the coefficients α, β
have the same meaning as in mul!
.
SummationByPartsOperators.mul_transpose_derivative_right!
— Methodmul_transpose_derivative_right!(u, D::DerivativeOperator, der_order::Val{N}, α=true, β=false)
Set the grid function u
to α
times the transposed N
-th derivative functional applied to u
plus β
times u
in the domain of the N
-th derivative functional at the right boundary of the grid. Thus, the coefficients α, β
have the same meaning as in mul!
.
SummationByPartsOperators.periodic_central_derivative_coefficients
— Functionperiodic_central_derivative_coefficients(derivative_order, accuracy_order, T=Float64, mode=FastMode())
Create the PeriodicDerivativeCoefficients
approximating the derivative_order
-th derivative with an order of accuracy accuracy_order
and scalar type T
. The evaluation of the derivative can be parallelized using threads by choosing mode=ThreadedMode())
.
SummationByPartsOperators.periodic_central_derivative_operator
— Functionperiodic_central_derivative_operator(derivative_order, accuracy_order,
xmin, xmax, N, mode=FastMode())
Create a PeriodicDerivativeOperator
approximating the derivative_order
-th derivative on a uniform grid between xmin
and xmax
with N
grid points up to order of accuracy accuracy_order
. The evaluation of the derivative can be parallelized using threads by choosing mode=ThreadedMode()
.
SummationByPartsOperators.periodic_derivative_coefficients
— Functionperiodic_derivative_coefficients(derivative_order, accuracy_order,
left_offset=-(accuracy_order+1)÷2,
T=Float64, mode=FastMode())
Create the PeriodicDerivativeCoefficients
approximating the derivative_order
-th derivative with an order of accuracy accuracy_order
and scalar type T
where the leftmost grid point used is determined by left_offset
. The evaluation of the derivative can be parallelized using threads by choosing mode=ThreadedMode()`.
SummationByPartsOperators.periodic_derivative_coefficients
— Methodperiodic_derivative_coefficients(source::Holoborodko2008, derivative_order, accuracy_order;
T=Float64, mode=FastMode(),
stencil_width=accuracy_order+3)
Create the PeriodicDerivativeCoefficients
approximating the derivative_order
-th derivative with an order of accuracy accuracy_order
and scalar type T
given by Holoborodko2008
. The evaluation of the derivative can be parallelized using threads by choosing mode=ThreadedMode()`.
SummationByPartsOperators.periodic_derivative_coefficients
— Methodperiodic_derivative_coefficients(source::LanczosLowNoise, derivative_order, accuracy_order;
T=Float64, mode=FastMode(),
stencil_width=accuracy_order+3)
Create the PeriodicDerivativeCoefficients
approximating the derivative_order
-th derivative with an order of accuracy accuracy_order
and scalar type T
for (super) Lanczos low Noise filters, see LanczosLowNoise
. The evaluation of the derivative can be parallelized using threads by choosing mode = ThreadedMode()
.
SummationByPartsOperators.periodic_derivative_operator
— Functionperiodic_derivative_operator(derivative_order, accuracy_order, grid,
left_offset=-(accuracy_order+1)÷2, mode=FastMode())
Create a PeriodicDerivativeOperator
approximating the derivative_order
-th derivative on the uniform grid
up to order of accuracy accuracy_order
where the leftmost grid point used is determined by left_offset
. The evaluation of the derivative can be parallelized using threads by choosing mode=ThreadedMode())
.
SummationByPartsOperators.periodic_derivative_operator
— Functionperiodic_derivative_operator(derivative_order, accuracy_order,
xmin, xmax, N,
left_offset=-(accuracy_order+1)÷2,
mode=FastMode())
periodic_derivative_operator(; derivative_order, accuracy_order,
xmin, xmax, N,
left_offset=-(accuracy_order+1)÷2,
mode=FastMode())
Create a PeriodicDerivativeOperator
approximating the derivative_order
-th derivative on a uniform grid between xmin
and xmax
with N
grid points up to order of accuracy accuracy_order
where the leftmost grid point used is determined by left_offset
. The evaluation of the derivative can be parallelized using threads by choosing mode=ThreadedMode())
.
Examples
julia> periodic_derivative_operator(derivative_order=1, accuracy_order=2,
xmin=0.0, xmax=1.0, N=11)
Periodic first-derivative operator of order 2 on a grid in [0.0, 1.0] using 11 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.
SummationByPartsOperators.periodic_derivative_operator
— Methodperiodic_derivative_operator(source::Holoborodko2008,
derivative_order, accuracy_order,
xmin, xmax, N; mode=FastMode(), kwargs...)
periodic_derivative_operator(source::Holoborodko2008;
derivative_order, accuracy_order,
xmin, xmax, N, mode=FastMode(), kwargs...)
Create a PeriodicDerivativeOperator
approximating the derivative_order
-th derivative on a uniform grid between xmin
and xmax
with N
grid points up to order of accuracy accuracy_order
where the leftmost grid point used is determined by left_offset
. The evaluation of the derivative can be parallelized using threads by choosing mode=ThreadedMode()
.
Examples
julia> periodic_derivative_operator(Holoborodko2008(), derivative_order=1, accuracy_order=2,
xmin=0.0, xmax=1.0, N=11)
Periodic first-derivative operator of order 2 on a grid in [0.0, 1.0] using 11 nodes,
stencils with 2 nodes to the left, 2 nodes to the right, and coefficients of Holoborodko (2008)
Smooth Noise Robust Differentiators.
http://www.holoborodko.com/pavel/numerical-methods/numerical-derivative/smooth-low-noise-differentiators/
SummationByPartsOperators.periodic_derivative_operator
— Methodperiodic_derivative_operator(source::LanczosLowNoise;
derivative_order = 1,
accuracy_order,
stencil_width = accuracy_order + 3,
xmin, xmax, N,
mode = FastMode())
Create a PeriodicDerivativeOperator
approximating the derivative_order
-th derivative on a uniform grid between xmin
and xmax
with N
grid points up to order of accuracy accuracy_order
where stencil_width
determines the number of nodes used for the finite difference stencil. The evaluation of the derivative can be parallelized using threads by choosing mode = ThreadedMode()
.
Examples
julia> D = periodic_derivative_operator(LanczosLowNoise();
accuracy_order = 2,
xmin = 0.0, xmax = 1.0, N = 10)
Periodic first-derivative operator of order 2 on a grid in [0.0, 1.0] using 10 nodes,
stencils with 2 nodes to the left, 2 nodes to the right, and coefficients of Holoborodko (2008)
Smooth Noise Robust Differentiators.
http://www.holoborodko.com/pavel/numerical-methods/numerical-derivative/lanczos-low-Noise-differentiators/
julia> Matrix(D)
10×10 Matrix{Float64}:
0.0 1.0 2.0 0.0 0.0 0.0 0.0 0.0 -2.0 -1.0
-1.0 0.0 1.0 2.0 0.0 0.0 0.0 0.0 0.0 -2.0
-2.0 -1.0 0.0 1.0 2.0 0.0 0.0 0.0 0.0 0.0
0.0 -2.0 -1.0 0.0 1.0 2.0 0.0 0.0 0.0 0.0
0.0 0.0 -2.0 -1.0 0.0 1.0 2.0 0.0 0.0 0.0
0.0 0.0 0.0 -2.0 -1.0 0.0 1.0 2.0 0.0 0.0
0.0 0.0 0.0 0.0 -2.0 -1.0 0.0 1.0 2.0 0.0
0.0 0.0 0.0 0.0 0.0 -2.0 -1.0 0.0 1.0 2.0
2.0 0.0 0.0 0.0 0.0 0.0 -2.0 -1.0 0.0 1.0
1.0 2.0 0.0 0.0 0.0 0.0 0.0 -2.0 -1.0 0.0
julia> D = periodic_derivative_operator(LanczosLowNoise();
accuracy_order = 2,
xmin = 0 // 1, xmax = 1 // 1, N = 10,
stencil_width = 7)
Periodic first-derivative operator of order 2 on a grid in [0//1, 1//1] using 10 nodes,
stencils with 3 nodes to the left, 3 nodes to the right, and coefficients of Holoborodko (2008)
Smooth Noise Robust Differentiators.
http://www.holoborodko.com/pavel/numerical-methods/numerical-derivative/lanczos-low-Noise-differentiators/
SummationByPartsOperators.right_boundary_weight
— Functionright_boundary_weight(D)
Return the left-boundary weight of the (diagonal) mass matrix M
associated to the derivative operator D
.
SummationByPartsOperators.scale_by_inverse_mass_matrix!
— Functionscale_by_inverse_mass_matrix!(u, D)
Scale the vector u
by the inverse of the mass matrix associated to the derivative operator D
.
SummationByPartsOperators.scale_by_mass_matrix!
— Functionscale_by_mass_matrix!(u, D)
Scale the vector u
by the mass matrix associated to the derivative operator D
.
SummationByPartsOperators.semidiscretize
— Methodsemidiscretize(u0func, semi::AbstractSemidiscretization, tspan)
Apply the semidiscretization semi
to the initial data given by u0func
and return an ODEProblem
with time span tspan
.
SummationByPartsOperators.source_of_coefficients
— Methodsource_of_coefficients(D)
Return the source of coefficients of the derivative operator D
. If you use the operator D
for your research, please cite this source in addition to SummationByPartsOperators
.
SummationByPartsOperators.upwind_operators
— Methodupwind_operators(source_type, args...; derivative_order = 1, kwargs...)
Create UpwindOperators
from the given source type. The positional arguments args...
and keyword arguments kwargs...
are passed directly to derivative_operator
.
Examples
julia> D = upwind_operators(Mattsson2017, accuracy_order = 2,
xmin = 0//1, xmax = 9//1, N = 10)
Upwind SBP first-derivative operators of order 2 on a grid in [0//1, 9//1] using 10 nodes
and coefficients of Mattsson2017
julia> D.minus
SBP first-derivative operator of order 2 on a grid in [0//1, 9//1] using 10 nodes
and coefficients of Mattsson (2017)
Diagonal-norm upwind SBP operators.
Journal of Computational Physics 335, pp. 283-310.
(upwind coefficients minus)
julia> D.plus
SBP first-derivative operator of order 2 on a grid in [0//1, 9//1] using 10 nodes
and coefficients of Mattsson (2017)
Diagonal-norm upwind SBP operators.
Journal of Computational Physics 335, pp. 283-310.
(upwind coefficients plus)
julia> Matrix(D.central)
10×10 Matrix{Rational{Int64}}:
-2//1 3//1 -1//1 0//1 0//1 0//1 0//1 0//1 0//1 0//1
-3//5 0//1 4//5 -1//5 0//1 0//1 0//1 0//1 0//1 0//1
1//4 -1//1 0//1 1//1 -1//4 0//1 0//1 0//1 0//1 0//1
0//1 1//4 -1//1 0//1 1//1 -1//4 0//1 0//1 0//1 0//1
0//1 0//1 1//4 -1//1 0//1 1//1 -1//4 0//1 0//1 0//1
0//1 0//1 0//1 1//4 -1//1 0//1 1//1 -1//4 0//1 0//1
0//1 0//1 0//1 0//1 1//4 -1//1 0//1 1//1 -1//4 0//1
0//1 0//1 0//1 0//1 0//1 1//4 -1//1 0//1 1//1 -1//4
0//1 0//1 0//1 0//1 0//1 0//1 1//5 -4//5 0//1 3//5
0//1 0//1 0//1 0//1 0//1 0//1 0//1 1//1 -3//1 2//1
SummationByPartsOperators.upwind_operators
— Methodupwind_operators(periodic_derivative_operator;
derivative_order = 1, accuracy_order,
xmin, xmax, N,
mode = FastMode()))
Create PeriodicUpwindOperators
from operators constructed by periodic_derivative_operator
. The keyword arguments are passed directly to periodic_derivative_operator
.
Examples
julia> D = upwind_operators(periodic_derivative_operator, accuracy_order = 2,
xmin = 0//1, xmax = 8//1, N = 8)
Upwind SBP first-derivative operators of order 2 on a grid in [0//1, 7//1] using 8 nodes
and coefficients of Fornberg1998
julia> D.minus
Periodic first-derivative operator of order 2 on a grid in [0//1, 8//1] using 8 nodes,
stencils with 2 nodes to the left, 0 nodes to the right, and coefficients of Fornberg (1998)
Calculation of Weights in Finite Difference Formulas.
SIAM Rev. 40.3, pp. 685-691.
julia> D.plus
Periodic first-derivative operator of order 2 on a grid in [0//1, 8//1] using 8 nodes,
stencils with 0 nodes to the left, 2 nodes to the right, and coefficients of Fornberg (1998)
Calculation of Weights in Finite Difference Formulas.
SIAM Rev. 40.3, pp. 685-691.
julia> Matrix(D.central)
8×8 Matrix{Rational{Int64}}:
0//1 1//1 -1//4 0//1 0//1 0//1 1//4 -1//1
-1//1 0//1 1//1 -1//4 0//1 0//1 0//1 1//4
1//4 -1//1 0//1 1//1 -1//4 0//1 0//1 0//1
0//1 1//4 -1//1 0//1 1//1 -1//4 0//1 0//1
0//1 0//1 1//4 -1//1 0//1 1//1 -1//4 0//1
0//1 0//1 0//1 1//4 -1//1 0//1 1//1 -1//4
-1//4 0//1 0//1 0//1 1//4 -1//1 0//1 1//1
1//1 -1//4 0//1 0//1 0//1 1//4 -1//1 0//1
SummationByPartsOperators.var_coef_derivative_operator
— Functionvar_coef_derivative_operator(source_of_coefficients, derivative_order, accuracy_order,
xmin, xmax, N, left_weights, right_weights, bfunc,
mode=FastMode())
Create a VarCoefDerivativeOperator
approximating a derivative_order
-th derivative with variable coefficients bfunc
on a grid between xmin
and xmax
with N
grid points up to order of accuracy accuracy_order
with coefficients given by source_of_coefficients
. The evaluation of the derivative can be parallelized using threads by choosing mode=ThreadedMode()
.
SummationByPartsOperators.xmax
— Functionxmax(D)
Return the right boundary xmax
of the domain specified when constructing the derivative operator D
. Note that this might be different from the rightmost node of the grid
of D
when not all boundary nodes are included, e.g., for periodic derivative operators.
SummationByPartsOperators.xmin
— Functionxmin(D)
Return the left boundary xmin
of the domain specified when constructing the derivative operator D
. Note that this might be different from the leftmost node of the grid
of D
when not all boundary nodes are included, e.g., for periodic derivative operators.