ImmersedPotentialFlowModule

This module implements an immersed Galerkin discretization of the Potential flow with (non)homogenous Dirichlet and Neumann boundary conditions.

Immersed potential flow around rotating cylinder

./examples/immersed_potentialflow/rotating_cylinder.jl

Immersed potential flow around torus

./examples/immersed_potentialflow/torus.jl

Index

Taiga.ImmersedPotentialFlowModule.ImmersedPotentialFlowType
ImmersedPotentialFlow{Dim, T, M <: GeometricMapping{Dim}, V <: ScalarSplineSpace{Dim, T}, W <: Field{Dim}}

Immersed Galerkin formulation of potential flow. The solution field can be either velocity potential or stream function.

Fields:

  • F::M: geometric mapping
  • S::V: constrained solution space
  • Δ::Partition{Dim, T}: partition
  • uʰ::W: solution field
  • ūʰ::W: function satisfying Dirichlet boundary condition
  • t::Function: traction vector
  • ϕ::AlgoimCallLevelSetFunction: level set function
source
Taiga.ImmersedPotentialFlowModule.LinearOperatorType
LinearOperator{Dim, T, S <: SparseMatrixCSC{T}, M <: AbstractMatrix{T}, V <: AbstractVector{T}}

Linear operator for ImmersedPotentialFlow model.

Fields:

  • C::S: sparse extraction operator
  • E::S: sparse extension operator
  • K::M: sparse stiffness matrix
  • L::S: sparse linear operator Eᵀ Cᵀ K C E
  • b::V: right hand side vector
source
Taiga.ImmersedPotentialFlowModule.PressureType
Pressure{Dim, T1, T2 <: Velocity{Dim}}

Pressure mapping for uniform far field flow.

Fields:

  • ρ::T1: density
  • p::T1: far field pressure
  • U::T1: far field velocity magnitude
  • velocity::T2: velocity field
source
Taiga.ImmersedPotentialFlowModule.VelocityType
Velocity{Dim, F <: PrimaryField, T1 <: GeometricMapping{Dim}, T2 <: Field{Dim}}

Velocity mapping. Supports potential as velocity potential or stream function.

Fields:

  • mapping::T1: geometric mapping
  • potential::T2: potential
source
Taiga.ImmersedPotentialFlowModule.assemble_bodyMethod
assemble_body(acc::ElementAccessor{Dim}, acc_immersed::ElementAccessor{Dim}, partition::Partition{Dim}, pullback::PullbackBody{Dim}, ϕ::AlgoimCallLevelSetFunction; show_progress::Bool = true)

Assemble unconstrained stiffness matrix.

Arguments:

  • acc: element accessor
  • acc_immersed: immersed element accessor
  • partition: partition
  • pullback: bilinear form pullback
  • ϕ: level set function
  • show_progress: boolean flag for formation progress bar
source
Taiga.ImmersedPotentialFlowModule.assemble_boundaryMethod
assemble_boundary(acc::ElementAccessor{Dim}, acc_immersed::ElementAccessor{Dim}, partition::Partition{Dim, T}, pullbacks_boundary::NTuple{N, PullbackBoundary{M}}, ϕ::AlgoimCallLevelSetFunction)

Assemble right hand side vector.

Arguments:

  • acc: element accessor
  • acc_immersed: immersed element accessor
  • partition: partition
  • pullbacks_boundary: tuple of 2Dim linear form pullbacks restricted to each boundary
  • ϕ: level set function
source
Taiga.apply_particular_solutionMethod
Taiga.apply_particular_solution(L::LinearOperator, model::ImmersedPotentialFlow, x₀::V)

Return particular solution.

Arguments:

  • L: linear operator
  • model: model
  • x₀: homogeneous solution
source
Taiga.forcing!Method
Taiga.forcing!(b::AbstractVector, L::LinearOperator, model::ImmersedPotentialFlow)

Update forcing vector.

Arguments:

  • b: cache vector
  • L: linear operator
  • model: model
source
Taiga.forcingMethod
Taiga.forcing(L::LinearOperator, model::ImmersedPotentialFlow)

Return foring vector.

source
Taiga.linear_operatorMethod
Taiga.linear_operator(model::ImmersedPotentialFlow; show_progress::Bool = true)

Construct linear operator.

source