PoissonModule

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

Poisson on horseshoe domain

./examples/poisson/horseshoe.jl

Index

Taiga.FastDiagonalizationMethod
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.

source
Taiga.FastDiagonalizationMethod
Taiga.FastDiagonalization(::Type{<:Wachspress}, model::Poisson; niter::Int)

Constructs a fast diagonalization preconditioner using niter of Wachspress algorithm to approximate patch data.

source
Taiga.PoissonModule.LinearOperatorType
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 matrix
  • b::V: rhs vector
  • C::B: Kronecker product extraction operator
source
Taiga.PoissonModule.MatrixfreeLinearOperatorType
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: model
  • acc::A: patch accessor
  • C::B: Kronecker product extraction operator
  • pullback_body::Function: pullback for the body terms
  • sumfact_cache_u::V: sum factorization vector cache
source
Taiga.PoissonModule.PoissonType
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 mapping
  • S::V: solution space
  • Δ::Partition{Dim,T}: partition of the space
  • uʰ::Field: solution field
  • ūʰ::W: scalar solution field satisfying Dirichlet boundary conditions
  • κ::Function: conductivity matrix of size Dim × Dim as function of spatial coordinates
  • f::Function: body force as function of spatial coordinates
  • t::Function: traction as function of spatial coordinates
source
Taiga.PoissonModule.PoissonMethod
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 mapping
  • S::V: scalar function space
  • C::ScalarSplineSpaceConstraints{Dim}: function space constraints
  • uʰ::W: scalar solution field
  • ūʰ::W: scalar solution field satisfying Dirichlet boundary conditions
  • κ::Function: conductivity matrix as a function of spatial coordinates
  • f::Function: body force as function of spatial coordinates
  • t::Function: traction vector as function of spatial coordinates
  • matrixfree::Bool: boolean flag for matrix-free forward problem, default is false
source
Taiga.apply_particular_solutionMethod
apply_particular_solution(L::LinearOperator, x₀::V) where {V<:AbstractVector}

Applies particular solution to homogeneous solution vector.

source
Taiga.forcing!Method
forcing!(b::AbstractVector, L::LinearOperator)

Updates forcing vector b.

source