SummationByPartsOperators.jl API
SummationByPartsOperators.SummationByPartsOperators — ModuleSummationByPartsOperators.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 noticable 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
@misc{ranocha2021sbp,
title={{SummationByPartsOperators.jl}: {A} {J}ulia library of provably stable
semidiscretization techniques with mimetic properties},
author={Ranocha, Hendrik},
year={2021},
howpublished={\url{https://github.com/ranocha/SummationByPartsOperators.jl},
doi={10.5281/zenodo.4773575}
}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.BurgersNonperiodicSemidiscretisation — TypeBurgersNonperiodicSemidiscretisation(D, Di, split_form, left_bc, right_bc)A semidiscretisation 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.BurgersPeriodicSemidiscretisation — TypeBurgersPeriodicSemidiscretisation(D, Di, split_form)A semidiscretisation 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 — TypeConstantFilterRepresents 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.CubicNonperiodicSemidiscretisation — TypeCubicNonperiodicSemidiscretisation(D, Di, split_form, left_bc, right_bc)A semidiscretisation 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.CubicPeriodicSemidiscretisation — TypeCubicPeriodicSemidiscretisation(D, Di, split_form)A semidiscretisation 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 — TypeDerivativeCoefficientsThe coefficients of a derivative operator on a nonperiodic grid.
SummationByPartsOperators.DerivativeOperator — TypeDerivativeOperatorA derivative operator on a nonperiodic finite difference grid.
SummationByPartsOperators.DissipationOperator — TypeDissipationOperatorA dissipation operator on a nonperiodic finite difference grid.
SummationByPartsOperators.ExponentialFilter — TypeExponentialFilterRepresents the exponential filter function σ(η) = exp(-α*η^p).
SummationByPartsOperators.FactorisationWrapper — TypeFactorisationWrapperA small wrapper around a a factorisation fact, allowing to represent multiplication by the inverse of fact.
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 — TypeFourierConstantViscosityFourier 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.
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.
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.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.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.MadayTadmor1989 — TypeMadayTadmor1989Coefficients 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.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 fiite 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 — TypePeriodicDerivativeCoefficientsThe coefficients of a derivative operator on a periodic grid.
SummationByPartsOperators.PeriodicDerivativeOperator — TypePeriodicDerivativeOperatorA derivative operator on a uniform periodic grid.
SummationByPartsOperators.PeriodicDissipationOperator — TypePeriodicDissipationOperatorA dissipation operator on a periodic finite difference grid.
SummationByPartsOperators.QuarticNonconvexPeriodicSemidiscretisation — TypeQuarticNonconvexPeriodicSemidiscretisation(D, Di, split_form)A semidiscretisation 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.SourceOfCoefficients — TypeSourceOfCoefficientsAll sources of coefficients (articles) are subtypes of this abstract type.
SummationByPartsOperators.SumOfDerivativeOperators — TypeSumOfDerivativeOperatorsSum several derivative operators lazily.
SummationByPartsOperators.Tadmor1989 — TypeTadmor1989Coefficients 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 — TypeTadmor1993Coefficients 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 — TypeTadmorWaagan2012ConvergentCoefficients 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 — TypeTadmorWaagan2012StandardCoefficients 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.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.VarCoefDerivativeCoefficients — TypeVarCoefDerivativeCoefficientsThe coefficients of a variable coefficient derivative operator on a nonperiodic grid.
SummationByPartsOperators.VarCoefDerivativeOperator — TypeVarCoefDerivativeOperatorA dissipation operator on a nonperiodic finite difference grid.
SummationByPartsOperators.VariableLinearAdvectionNonperiodicSemidiscretisation — TypeVariableLinearAdvectionNonperiodicSemidiscretisation(D, Di, a, split_form,
left_bc, right_bc)A semidiscretisation 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.WaveEquationNonperiodicSemidiscretisation — TypeWaveEquationNonperiodicSemidiscretisation(D, left_bc, right_bc)A semidiscretisation 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).
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.add_transpose_derivative_left! — Methodadd_transpose_derivative_left!(u, D::DerivativeOperator, der_order::Val{N}, α)Add α times the transposed N-th derivative functional to the grid function u at the left boundary of the grid.
SummationByPartsOperators.add_transpose_derivative_right! — Methodadd_transpose_derivative_right!(u, D::DerivativeOperator, der_order::Val{N}, α)Add α times the transposed N-th derivative functional to the grid function u at the right boundary of the grid.
SummationByPartsOperators.derivative_left — Methodderivative_left(D::AbstractNonperiodicDerivativeOperator, der_order)Get a representation of the linear functional evaluation the Nth 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, parallel=Val{:serial}())
derivative_operator(source_of_coefficients;
derivative_order, accuracy_order,
xmin, xmax, N, parallel=Val{:serial}())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 parallised using threads by chosing parallel=Val{:threads}()).
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 Nth 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, parallel=Val{:serial}())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 parallised using threads by chosing parallel=Val{:threads}()).
SummationByPartsOperators.dissipation_operator — Methoddissipation_operator(D::DerivativeOperator; kwargs...)Create a negative semidefinite DissipationOperator using undivided differences approximating a weighted order-th derivative adapted to the derivative operator D. The evaluation of the derivative can be parallised using threads by chosing parallel=Val{:threads}()).
SummationByPartsOperators.dissipation_operator — Methoddissipation_operator(D::PeriodicDerivativeOperator;
strength=one(eltype(D)),
order=accuracy_order(D),
parallel=D.coefficients.parallel)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 parallised using threads by chosing parallel=Val{:threads}()).
SummationByPartsOperators.dissipation_operator — Methoddissipation_operator(source_of_coefficients, D::DerivativeOperator{T};
strength=one(T),
order::Int=accuracy_order(D),
parallel=D.coefficients.parallel)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 parallised using threads by chosing parallel=Val{:threads}()).
SummationByPartsOperators.fornberg — Methodfornberg(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.
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.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.periodic_central_derivative_coefficients — Functionperiodic_central_derivative_coefficients(derivative_order, accuracy_order, T=Float64, parallel=Val{:serial}())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 parallised using threads by chosing parallel=Val{:threads}()).
SummationByPartsOperators.periodic_central_derivative_operator — Functionperiodic_central_derivative_operator(derivative_order, accuracy_order,
grid, parallel=Val{:serial}())Create a PeriodicDerivativeOperator approximating the derivative_order-th derivative on the uniform grid up to order of accuracy accuracy_order. The evaluation of the derivative can be parallised using threads by chosing parallel=Val{:threads}()).
SummationByPartsOperators.periodic_central_derivative_operator — Functionperiodic_central_derivative_operator(derivative_order, accuracy_order,
xmin, xmax, N, parallel=Val{:serial}())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 parallised using threads by chosing parallel=Val{:threads}()).
SummationByPartsOperators.periodic_derivative_coefficients — Functionperiodic_derivative_coefficients(derivative_order, accuracy_order, left_offset=-(accuracy_order+1)÷2, T=Float64, parallel=Val{:serial}())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 parallised using threads by chosing parallel=Val{:threads}())`.
SummationByPartsOperators.periodic_derivative_coefficients — Methodperiodic_derivative_coefficients(source::Holoborodko2008, derivative_order, accuracy_order;
T=Float64, parallel=Val{:serial}(),
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 parallised using threads by chosing parallel=Val{:threads}())`.
SummationByPartsOperators.periodic_derivative_operator — Functionperiodic_derivative_operator(derivative_order, accuracy_order,
xmin, xmax, N,
left_offset=-(accuracy_order+1)÷2,
parallel=Val{:serial}())
periodic_derivative_operator(; derivative_order, accuracy_order,
xmin, xmax, N,
left_offset=-(accuracy_order+1)÷2,
parallel=Val{:serial}())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 parallised using threads by chosing parallel=Val{:threads}()).
Examples
julia> periodic_derivative_operator(derivative_order=1, accuracy_order=2,
xmin=0.0, xmax=1.0, N=11)
Periodic 1st derivative operator of order 2 {T=Float64, Parallel=Val{:serial}}
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 from
Fornberg (1998)
Calculation of Weights in Finite Difference Formulas.
SIAM Rev. 40.3, pp. 685-691.
SummationByPartsOperators.periodic_derivative_operator — Functionperiodic_derivative_operator(derivative_order, accuracy_order, grid,
left_offset=-(accuracy_order+1)÷2, parallel=Val{:serial}())Create a PeriodicDerivativeOperator approximating the derivative_order-th derivative on thr 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 parallised using threads by chosing parallel=Val{:threads}()).
SummationByPartsOperators.periodic_derivative_operator — Methodperiodic_derivative_operator(source::Holoborodko2008,
derivative_order, accuracy_order,
xmin, xmax, N; parallel=Val{:serial}(), kwargs...)
periodic_derivative_operator(source::Holoborodko2008;
derivative_order, accuracy_order,
xmin, xmax, N, parallel=Val{:serial}(), 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 parallised using threads by chosing parallel=Val{:threads}()).
Examples
julia> periodic_derivative_operator(Holoborodko2008(), derivative_order=1, accuracy_order=2,
xmin=0.0, xmax=1.0, N=11)
Periodic 1st derivative operator of order 2 {T=Float64, Parallel=Val{:serial}}
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 from
Holoborodko (2008)
Smooth Noise Robust Differentiators.
http://www.holoborodko.com/pavel/numerical-methods/numerical-derivative/smooth-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.semidiscretize — Methodsemidiscretize(u0func, semidisc::AbstractSemidiscretisation, tspan)Apply the semidiscretisation semidisc 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.var_coef_derivative_operator — Functionvar_coef_derivative_operator(source_of_coefficients, derivative_order, accuracy_order, xmin, xmax, N, left_weights, right_weights, bfunc, parallel=Val{:serial}())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 parallised using threads by chosing parallel=Val{:threads}()).