PotentialFlowModule
This module implements a Galerkin discretization of the potential flow with (non)homogenous Dirichlet and Neumann boundary conditions.

./examples/immersed_potentialflow/cylinder.jl
Index
Taiga.FastDiagonalization — Method
Taiga.FastDiagonalization(model::PotentialFlow{Dim, T}; method::Type{<:DataApproximationMethod} = Wachspress, kwargs)Potential flow fast diagonalization preconditioner.
Taiga.PotentialFlowModule.LinearOperator — Type
LinearOperator{Dim, T, B <: KroneckerProduct{T}, M <: AbstractMatrix{T}, V <: AbstractVector{T}}Linear operator for PotentialFlow model.
Fields:
C::B: Kronecker product extraction operatorK::M: sparse stiffness matrixb::V: forcing vectorcache₁::Vector{T}:mul!cachecache₂::Vector{T}:mul!cache
Taiga.PotentialFlowModule.PotentialFlow — Type
PotentialFlow{Dim, T, M <: GeometricMapping{Dim}, V <: ScalarSplineSpace{Dim, T}, W <: Field{Dim}}Galerkin formulation of potential flow. The solution field can be either velocity potential or stream function.
Fields:
F::M: geometric mappingS::V: constrained solution spaceΔ::Partition{Dim, T}: partitionuʰ::W: solution fieldūʰ::W: function satisfying Dirichlet boundary conditionst::Function: traction vector
Taiga.PotentialFlowModule.Pressure — Type
Pressure{Dim, T1, T2 <: Velocity{Dim}}Pressure mapping for uniform far field flow.
Fields:
ρ::T1: densityp::T1: far field pressureU::T1: far field velocity magnitudevelocity::T2: velocity field
Taiga.PotentialFlowModule.PullbackBody — Type
PullbackBody{Dim, T <: GeometricMapping{Dim}}Pullback of bilinear form.
Fields:
mapping::T: geometric mapping
Taiga.PotentialFlowModule.PullbackBoundary — Type
PullbackBoundary{Dim, T <: GeometricMapping{Dim}}Pullback of linear form.
Fields:
mapping::T: geometric mappingtraction::Function: traction vectorside::Int: side restriction (side <= 2Dim)
Taiga.PotentialFlowModule.Velocity — Type
Velocity{Dim, T1 <: GeometricMapping{Dim}, T2 <: Field{Dim}}Velocity mapping.
Fields:
mapping::T1: geometric mappingpotential::T2: potential field
Taiga.PotentialFlowModule.assemble_body — Method
assemble_body(acc::ElementAccessor{Dim}, partition::Partition{Dim}, pullback::PullbackBody{Dim})Assemble unconstrained stiffness matrix.
Arguments:
acc: element accessorpartition: partitionpullback: bilinear form pullback
Taiga.PotentialFlowModule.assemble_boundary — Method
assemble_boundary(acc::ElementAccessor{Dim}, partition::Partition{Dim, T}, pullbacks_boundary::NTuple{N, PullbackBoundary{M}})Assemble right hand side vector.
Arguments:
acc: element accessorpartition: partitionpullbacks_boundary: tuple of2Dimlinear form pullbacks restricted to each boundary
Taiga.apply_particular_solution — Method
Taiga.apply_particular_solution(L::LinearOperator, model::PotentialFlow, x₀::V)Return particular solution.
Arguments:
L: linear operatormodel: modelx₀: homogeneous solution
Taiga.forcing! — Method
Taiga.forcing!(b::AbstractVector, L::LinearOperator, model::PotentialFlow)Update forcing vector.
Arguments:
b: cache vectorL: linear operatormodel: model
Taiga.forcing — Method
Taiga.forcing(L::LinearOperator, model::PotentialFlow)Return forcing vector.
Taiga.linear_operator — Method
Taiga.linear_operator(model::PotentialFlow)Construct linear operator.