Docstrings
ElementaryPDESolutions.Polynomial
— Typestruct Polynomial{N,T}
A polynomial in N
variables with coefficients of type T
.
The functor interface is implemented, so that p(x)
evaluates the polynomial. For performance reasons, x
is expected to be a Tuple
.
Examples
A polynomial with a single term can be created by passing a pair mapping the order to the coefficient:
julia> Polynomial((1,1)=>2)
2xy
When multiple terms are present, they must be passed as vector (or tuple) of pairs:
julia> Polynomial([(1,1)=>2, (0,1)=>-1])
-y + 2xy
The spatial dimension is automatically inferred from the length of the order tuple:
julia> Polynomial((1,1,1)=>2)
2xyz
ElementaryPDESolutions.anisotropic_laplacian
— Methodanisotropic_laplacian(A::AbstractMatrix, P::Polynomial)
Evaluate the anisotropic Laplacian ∇ ⋅ (A ∇P)
.
ElementaryPDESolutions.brinkman_component_solver
— Methodbrinkman_component_solver(Q::Polynomial{N,T}, α) where {N,T}
Compute a polynomial vector potential P
satisfying the auxiliary vector PDE.
(Δ³ - α⁴Δ)P = Q
for the Brinkman (linearized Navier-Stokes) system.
ElementaryPDESolutions.convert_coefs
— Methodconvert_coefs(p::Polynomial, T)
Return a version of p
where the coefficients have been converted to type T
(is such a conversion is possible).
ElementaryPDESolutions.degree
— Methoddegree(p::Polynomial)
The largest degree of any monomial in p
.
ElementaryPDESolutions.derivative
— Methodderivative(p::Polynomial, i::Int)
Differentiate p
with respect to the i
th variable.
ElementaryPDESolutions.drop_zeros!
— Functiondrop_zeros!(q::Polynomial,tol=0,p=2)
Drop all coefficients in q
for which the abs(p) ≤ tol
.
ElementaryPDESolutions.gradient
— Methodgradient(p::Polynomial)
Return an N
-tuple of the derivatives of p
with respect to each variable.
ElementaryPDESolutions.is_homogeneous
— Methodis_homogeneous(p::Polynomial)
Return true
if p
is homogeneous, i.e. if all the monomials in p
have the same degree.
ElementaryPDESolutions.multiply_by_anisotropic_r
— Methodmultiply_by_anisotropic_anisotropic_r(A::AbstractMatrix{T}, p::Polynomial, k::Int = 2)
Multiply a polynomial p
by the polynomial r_A^k
, where r_A = |r^T A^{-1} r|
, r = (x1, x2, ..., x_n), and k
is an even positive integer.
ElementaryPDESolutions.multiply_by_anisotropic_β_r
— Methodmultiply_by_anisotropic_β_r(β::AbstractVector{T}, p::Polynomial, k::Int)
Multiply a polynomial p
by the polynomial (β ⋅ 𝐫)ᵏ, 𝐫 = (x1, x2, ..., x_n), and k
is a non-negative integer.
ElementaryPDESolutions.multiply_by_r
— Methodmultiply_by_r(p::Polynomial, k::Int = 2)
Multiply a polynomial p
by the polynomial r^k
, where r = |𝐱|
and k
is an even positive integer.
ElementaryPDESolutions.solve_anisotropic_advect
— Methodsolve_anisotropic_advect(β::AbstractVector, Q::Polynomial)
Return a polynomial P
satisfying the anisotropic advection equation β⋅∇P = Q
.
Examples
using StaticArrays
β = SVector{2,Rational{Int64}}(2 // 1, 1 // 1)
Q = Polynomial([(0, 0) => 2 // 1])
P = solve_anisotropic_advect(β, Q)
# output
2//5y + 4//5x
ElementaryPDESolutions.solve_anisotropic_advect_diffuse
— Methodsolve_anisotropic_advect_diffuse(A::SMatrix{N, N}, β::AbstractVector{T}, Q::Polynomial)
Return a polynomial P
satisfying the anisotropic advection-diffusion equation ∇ ⋅ (A ∇P) + β⋅∇P = Q
, A
a symmetric positive definite matrix.
Examples
using StaticArrays
A = SMatrix{2,2,Rational{Int64}}(2 // 1, 1 // 1, 1 // 1, 3 // 1)
β = SVector{2,Rational{Int64}}(2 // 1, 1 // 1)
Q = Polynomial([(0, 1) => 2 // 1])
P = solve_anisotropic_advect_diffuse(A, β, Q)
# output
-14//25y - 28//25x + 16//25xy + 9//25y² - 4//25x²
ElementaryPDESolutions.solve_anisotropic_laplace
— Methodsolve_anisotropic_laplace(A::AbstractMatrix{T}, Q::Polynomial)
Return a polynomial P
satisfying the anisotropic Laplace equation ∇ ⋅ (A ∇P) = Q
, A
an invertible matrix. Q
is required to be homogeneous. Inverse is anisotropic_laplacian
.
Examples
using StaticArrays
A = SMatrix{2,2,Rational{Int64}}(2 // 1, 1 // 1, 1 // 1, 3 // 1)
Q = Polynomial([(1, 1) => 2 // 1])
P = solve_anisotropic_laplace(A, Q)
# output
-3//400x⁴ + 11//100x³y + 11//150xy³ - 2//25x²y² - 1//300y⁴
ElementaryPDESolutions.solve_bilaplace
— Methodsolve_bilaplace(Q::Polynomial)
Compute a polynomial solution to Δ²P = Q
. Q
is required to be homogeneous.
Examples
julia> Q = Polynomial((1,0)=>1)
x
julia> P = solve_bilaplace(Q)
1//192x⁵ + 1//96x³y² + 1//192xy⁴
ElementaryPDESolutions.solve_brinkman
— Methodsolve_brinkman(Q::NTuple{N,Polynomial{N,T}};Re=1,α=1)
Compute a vector of polynomials U
and a polynomial P
satisfying the linearized unsteady Navier-Stokes equations, sometimes referred to as the Brinkman equations or the modified Stokes equations, (Δ - α²)U - Re ∇P = Q
with ∇⋅U = 0
. Each component of the polynomial Q
is required to be individually homogeneous.
The solutions are given by the expressions
u = (Δ + α²)(Δ - ∇∇⋅)g,
p = -1/Re (Δ² - α⁴)∇⋅g,
where the vector potential g satisfies
(Δ³ - α⁴Δ)g = Q.
Examples
julia> Q = (Polynomial([(2, 1) => 2 // 1]), Polynomial([(0, 2) => 4 // 1]))
(2//1x²y, 4//1y²)
julia> U, P = solve_brinkman(Q; Re=Rational(1), α=Rational(1))
((0//1y + xy + 5//24y³ - 5//8x²y, 0//1 + 4//1x + 1//2x² - 1//2y² + 5//8xy² + 11//24x³), -7//6y³ - 1//2x²y - 11//24x³y - 5//24xy³)
ElementaryPDESolutions.solve_elastodynamics
— Methodsolve_elastodynamics(Q::NTuple{N,Polynomial{N,T}};ρ=1,μ=1,ν=1/4,ω=1)
Compute a vector of polynomials U
satisfying -μ/(1-2ν) ∇(div U) - μ ΔU - μ k₂² U = Q
.
Examples
julia> Q = (Polynomial((2,1)=>1),Polynomial((1,0)=>1))
(x²y, x)
julia> P = solve_elastodynamics(Q;μ=1)
(-6//1y + x²y, -3//1x)
ElementaryPDESolutions.solve_elastostatic
— Methodsolve_elastostatic(Q::NTuple{N,Polynomial{N,T}};μ=1,ν=1)
Compute a vector of polynomials U
satisfying μ/(1-2ν) ∇(div U) + μΔU = Q
. Q
is required to be homogeneous.
Examples
julia> Q = (Polynomial((1,2)=>1), Polynomial((0,0)=>1))
(xy², 1)
julia> P = solve_elastostatic(Q;ν=1//2)
(-1//8xy + 1//480x⁵ + 1//32x³y² + 1//24xy⁴, 3//16x² + 1//16y² - 1//120y⁵ - 1//96x⁴y - 1//32x²y³)
ElementaryPDESolutions.solve_helmholtz
— Methodsolve_helmholtz(Q::Polynomial;k=1)
Return the unique polynomial P
satisfying ΔP + k²P = Q
.
Examples
julia> Q = Polynomial((1,2)=>1)
xy²
julia> P = solve_helmholtz(Q, k=1)
-2.0x + xy²
ElementaryPDESolutions.solve_laplace
— Methodsolve_laplace(Q::Polynomial)
Return a polynomial P
satisfying ΔP = Q
. Q
is required to be homogeneous.
Examples
julia> Q = Polynomial((1,0)=>1.0)
x
julia> P = solve_laplace(Q)
0.125xy² + 0.125x³
ElementaryPDESolutions.solve_maxwell
— Methodsolve_maxwell(J::NTuple{3,Polynomial{3,T}};ϵ=1,μ=1,ω=1)
Compute a pair of vectors of polynomials E
and H
satisfying the Maxwell system:
\[\begin{aligned} \mathrm{i}\omega\varepsilon\boldsymbol{E} + \operatorname{rot} \boldsymbol{H} &= \boldsymbol{J}, \qquad & -\mathrm{i}\omega\mu\boldsymbol{H} + \operatorname{rot}\boldsymbol{E} &= \boldsymbol{0}, \\ \varepsilon\operatorname{div}\boldsymbol{E} &= \rho, & \mu\operatorname{div}\boldsymbol{H} &= 0, \end{aligned}\]
with the sources being constrained by the charge conservation equation:
\[\begin{aligned} \operatorname{div}\boldsymbol{J} - \mathrm{i}\omega\rho &= 0. \end{aligned}\]
Returns the pair (E, H)
.
Examples
julia> J = (Polynomial((0,2,1) => 1), Polynomial((0,1,1) => 1), Polynomial((1,0,1) => 1))
(y²z, yz, xz)
julia> E, H = solve_maxwell(J);
julia> E
((-0.0 - 1.0im) + (0.0 + 2.0im)z + (-0.0 - 1.0im)y²z, (-0.0 - 1.0im)yz, (-0.0 - 1.0im) + (-0.0 - 1.0im)xz)
julia> H
(y, 2.0 + z - y², 2.0yz)
ElementaryPDESolutions.solve_stokes
— Methodsolve_stokes(Q::NTuple{N,Polynomial{N,T}};μ=1)
Compute a vector of polynomials U
and a polynomial P
satisfying μΔU - ∇P = Q
with ∇ ⋅ U = 0
. Q
is required to be homogeneous.
Examples
julia> Q = (Polynomial((1,0)=>1),Polynomial((0,0)=>1))
(x, 1)
julia> P = solve_stokes(Q;μ=Rational(1))
((-1//8xy + 1//16xy² + 1//48x³, 3//16x² + 1//16y² - 1//48y³ - 1//16x²y), -1//2y - 3//8x² - 1//8y²)