SummationByPartsOperators.jl API

SummationByPartsOperators.SummationByPartsOperatorsModule
SummationByPartsOperators

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}
}
source
SummationByPartsOperators.BeljaddLeFlochMishraParés2017Type
BeljaddLeFlochMishraParé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.
source
SummationByPartsOperators.BurgersNonperiodicSemidiscretizationType
BurgersNonperiodicSemidiscretization(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.

source
SummationByPartsOperators.BurgersPeriodicSemidiscretizationType
BurgersPeriodicSemidiscretization(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.

source
SummationByPartsOperators.ConstantFilterMethod
ConstantFilter(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.

source
SummationByPartsOperators.ConstantFilterMethod
ConstantFilter(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.

source
SummationByPartsOperators.CubicNonperiodicSemidiscretizationType
CubicNonperiodicSemidiscretization(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.

source
SummationByPartsOperators.CubicPeriodicSemidiscretizationType
CubicPeriodicSemidiscretization(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.

source
SummationByPartsOperators.DienerDorbandSchnetterTiglio2007Type
DienerDorbandSchnetterTiglio2007()

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.
source
SummationByPartsOperators.Holoborodko2008Type
Holoborodko2008()

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/
source
SummationByPartsOperators.MadayTadmor1989Type
MadayTadmor1989()

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.
source
SummationByPartsOperators.MatrixDerivativeOperatorType
MatrixDerivativeOperator{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::Matrix 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.

source
SummationByPartsOperators.Mattsson2012Type
Mattsson2012()

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.
source
SummationByPartsOperators.Mattsson2014Type
Mattsson2014()

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.
source
SummationByPartsOperators.Mattsson2017Type
Mattsson2017(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.

source
SummationByPartsOperators.MattssonNordström2004Type
MattssonNordströ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.
source
SummationByPartsOperators.MattssonSvärdShoeybi2008Type
MattssonSvä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.
source
SummationByPartsOperators.PeriodicUpwindOperatorsType
PeriodicUpwindOperators
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

source
SummationByPartsOperators.QuarticNonconvexPeriodicSemidiscretizationType
QuarticNonconvexPeriodicSemidiscretization(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.

source
SummationByPartsOperators.SharanBradyLivescu2022Type
SharanBradyLivescu2022(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.

source
SummationByPartsOperators.Tadmor1989Type
Tadmor1989()

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.
source
SummationByPartsOperators.Tadmor1993Type
Tadmor1993()

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.
source
SummationByPartsOperators.TadmorWaagan2012ConvergentType
TadmorWaagan2012Convergent()

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.
source
SummationByPartsOperators.TadmorWaagan2012StandardType
TadmorWaagan2012Standard()

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.
source
SummationByPartsOperators.UniformMesh1DType
UniformMesh1D(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.

source
SummationByPartsOperators.UpwindOperatorsType
UpwindOperators
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

source
SummationByPartsOperators.VariableLinearAdvectionNonperiodicSemidiscretizationType
VariableLinearAdvectionNonperiodicSemidiscretization(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.

source
SummationByPartsOperators.VariableLinearAdvectionPeriodicSemidiscretizationType
VariableLinearAdvectionPeriodicSemidiscretization(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.

source
SummationByPartsOperators.WaveEquationNonperiodicSemidiscretizationType
WaveEquationNonperiodicSemidiscretization(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).

source
LinearAlgebra.mul!Function
mul!(du, D::DerivativeOperator, u, α=true, β=false)

Efficient in-place version of du = α * D * u + β * du. Note that du must not be aliased with u.

source
PolynomialBases.compute_coefficients!Method
compute_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.

source
PolynomialBases.evaluate_coefficients!Method
evaluate_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.

source
PolynomialBases.evaluate_coefficientsMethod
evaluate_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.

source
PolynomialBases.integrateMethod
integrate(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.

source
PolynomialBases.integrateMethod
integrate(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.

source
PolynomialBases.integrateMethod
integrate(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.

source
SummationByPartsOperators.couple_continuouslyFunction
couple_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 LegendreDerivativeOperators, 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

source
SummationByPartsOperators.couple_discontinuouslyFunction
couple_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 LegendreDerivativeOperators, 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 properties
  • Val(:minus) or Val(:plus), resulting in upwind SBP operators

References

source
SummationByPartsOperators.derivative_leftMethod
derivative_left(D::AbstractNonperiodicDerivativeOperator, der_order)

Get a representation of the linear functional evaluation the Nth derivative at the left boundary node as (dense) vector.

source
SummationByPartsOperators.derivative_operatorFunction
derivative_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().

source
SummationByPartsOperators.derivative_rightMethod
derivative_right(D::AbstractNonperiodicDerivativeOperator, der_order)

Get a representation of the linear functional evaluation the Nth derivative at the right boundary node as (dense) vector.

source
SummationByPartsOperators.dissipation_operatorFunction
dissipation_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().

source
SummationByPartsOperators.dissipation_operatorMethod
dissipation_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().

source
SummationByPartsOperators.dissipation_operatorMethod
dissipation_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().

source
SummationByPartsOperators.fornbergMethod
fornberg(x::Vector{T}, m::Int) where {T}

Calculate the weights of a finite difference approximation of the mth 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.

source
SummationByPartsOperators.fourier_derivative_matrixFunction
fourier_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.

source
SummationByPartsOperators.legendre_derivative_operatorMethod
legendre_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.

source
SummationByPartsOperators.legendre_second_derivative_operatorMethod
legendre_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.

source
SummationByPartsOperators.mul_transpose_derivative_left!Method
mul_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!.

source
SummationByPartsOperators.mul_transpose_derivative_right!Method
mul_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!.

source
SummationByPartsOperators.periodic_central_derivative_coefficientsFunction
periodic_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()).

source
SummationByPartsOperators.periodic_derivative_coefficientsFunction
periodic_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()`.

source
SummationByPartsOperators.periodic_derivative_coefficientsMethod
periodic_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()`.

source
SummationByPartsOperators.periodic_derivative_operatorFunction
periodic_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()).

source
SummationByPartsOperators.periodic_derivative_operatorFunction
periodic_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.
source
SummationByPartsOperators.periodic_derivative_operatorMethod
periodic_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/
source
SummationByPartsOperators.semidiscretizeMethod
semidiscretize(u0func, semi::AbstractSemidiscretization, tspan)

Apply the semidiscretization semi to the initial data given by u0func and return an ODEProblem with time span tspan.

source
SummationByPartsOperators.upwind_operatorsMethod
upwind_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
source
SummationByPartsOperators.upwind_operatorsMethod
upwind_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
source
SummationByPartsOperators.var_coef_derivative_operatorFunction
var_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().

source
SummationByPartsOperators.xmaxFunction
xmax(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.

source
SummationByPartsOperators.xminFunction
xmin(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.

source