PoissonModule
This module implements a Galerkin discretization of the Poisson equation with (non)homogenous Dirichlet and Neumann boundary conditions.

./examples/poisson/horseshoe.jl
Index
Taiga.FastDiagonalization — Method
Taiga.FastDiagonalization(model::Poisson{Dim,T}; method::Type{<:DataApproximationMethod}=Wachspress, kwargs...) where {Dim,T}Constructs a FastDiagonalization preconditioner using a particular data approximation method. Poisson supports Wachspress and ModalSplines data approximation methods.
Taiga.FastDiagonalization — Method
Taiga.FastDiagonalization(::Type{<:ModalSplines}, model::Poisson; spaces::NTuple{Dim, SplineSpace{T}})Constructs a fast diagonalization preconditioner using rank 1 modal splines approximation of model patch data on univariate spline spaces in spaces.
Taiga.FastDiagonalization — Method
Taiga.FastDiagonalization(::Type{<:Wachspress}, model::Poisson; niter::Int)Constructs a fast diagonalization preconditioner using niter of Wachspress algorithm to approximate patch data.
Taiga.PoissonModule.LinearOperator — Type
LinearOperator{Dim,T,M<:AbstractMatrix{T},V<:AbstractVector{T},B<:KroneckerProduct{T}} <: Taiga.LinearOperator{Dim,T}Linear operator for the homogeneous system in Poisson model.
Fields:
K::M: stiffness matrixb::V: rhs vectorC::B: Kronecker product extraction operator
Taiga.PoissonModule.LinearOperatorApproximation — Type
LinearOperatorApproximation{Dim,T} <: Taiga.LinearOperatorApproximation{Dim,T}Linear operator approximation for the homogeneous system in Poisson model.
Fields:
L::KroneckerProductAggregate{T}: Kronecker product aggregate
Taiga.PoissonModule.MatrixfreeLinearOperator — Type
MatrixfreeLinearOperator{Dim, T, M <: Poisson{Dim, T}, A <: Accessor{Dim}, B <: KroneckerProduct{T}, V <: VectorSumfactoryCache{Dim}} <: Taiga.LinearOperator{Dim,T}Matrix-free linear operator for the homogeneous system in Poisson model.
Fields:
model::M: modelacc::A: patch accessorC::B: Kronecker product extraction operatorpullback_body::Function: pullback for the body termssumfact_cache_u::V: sum factorization vector cache
Taiga.PoissonModule.Poisson — Type
Poisson{Dim,T,M<:GeometricMapping{Dim},V<:ScalarSplineSpace{Dim,T},W<:Field{Dim}} <: Model{Dim,T}Strong form
-∇ ⋅ (κ∇u) = f in Ω
u = ū on Γ₁
∇u ⋅ n = t ⋅ n on Γ₂Weak form
Find the homogeneous solution uₒ ∈ Sʰ, s.t. ∀v ∈ Sʰ
∫ ∇(uₒ - ū) ⋅ ∇v dΩ - ∫ (t ⋅ n) v dΓ₁ = - ∫ f v dΓ
for the homogeneous problem uₒ = 0 on Γ. The solution
u ∈ ū ⊕ Sʰ satisfying nonhomogeneous Dirichlet boundary
conditions is u = uₒ + ū.Note
The linear operators Taiga.PoissonModule.LinearOperator and Taiga.PoissonModule.MatrixfreeLinearOperator correspond to the homogeneous problem on the constrained space. The solution of the nonhomogeneous problem is obtained by Taiga.PoissonModule.apply_particular_solution.
Fields:
F::M: geometric mappingS::V: solution spaceΔ::Partition{Dim,T}: partition of the spaceuʰ::Field: solution fieldūʰ::W: scalar solution field satisfying Dirichlet boundary conditionsκ::Function: conductivity matrix of sizeDim × Dimas function of spatial coordinatesf::Function: body force as function of spatial coordinatest::Function: traction as function of spatial coordinates
Taiga.PoissonModule.Poisson — Method
Poisson(F::M, S::V, C::ScalarSplineSpaceConstraints{Dim}, uʰ::W, ūʰ::W, κ::Function, f::Function, t::Function; matrixfree=false) where {Dim,T,M<:GeometricMapping{Dim},V<:ScalarFunctionSpace{Dim,T},W<:Field{Dim}}Construct Poisson model.
Arguments:
F::M: geometric mappingS::V: scalar function spaceC::ScalarSplineSpaceConstraints{Dim}: function space constraintsuʰ::W: scalar solution fieldūʰ::W: scalar solution field satisfying Dirichlet boundary conditionsκ::Function: conductivity matrix as a function of spatial coordinatesf::Function: body force as function of spatial coordinatest::Function: traction vector as function of spatial coordinatesmatrixfree::Bool: boolean flag for matrix-free forward problem, default isfalse
Taiga.apply_particular_solution — Method
apply_particular_solution(model::Poisson, x₀::V) where {V<:AbstractVector}Applies particular solution to homogeneous solution vector.
Taiga.apply_particular_solution — Method
apply_particular_solution(L::LinearOperator, x₀::V) where {V<:AbstractVector}Applies particular solution to homogeneous solution vector.
Taiga.forcing! — Method
forcing!(b::AbstractVector, L::LinearOperator)Updates forcing vector b.
Taiga.forcing — Method
forcing(L::LinearOperator)Returns forcing vector.
Taiga.linear_operator — Method
linear_operator(model::Poisson; matrixfree::Bool=false, constrained::Bool=true)Returns a linear operator corresponding for model. Depending on the matrixfree keyword either PoissonModule.LinearOperator or PoissonModule.MatrixfreeLinearOperator is returned. For PoissonModule.MatrixfreeLinearOperator the constrained flag controls if the extraction operator is built into the matrix-free linear operator.
Taiga.linear_operator_approximation — Method
linear_operator_approximation(model::Poisson; method::Type{<:DataApproximationMethod})Returns PoissonModule.LinearOperatorApproximation for Poisson model using one of Taiga.DataApproximationMethod.