Array Pools

Immersa.StructuralStateType
StructuralState{T,A<:AbstractVector{T}}

Container for the kinematic state of a deformable structure, storing displacements, velocities, and accelerations in generalised (structural) coordinates.

Fields

  • χ::A : Structural displacements.
  • ζ::A : Structural velocities.
  • ζdot::A : Structural accelerations.

Constructor

StructuralState{T}(backend, n)

Allocates a zero-initialised state with n structural degrees of freedom on the given computation backend (e.g. CPU()).

Arguments

  • backend : Computation backend (CPU() or GPU device).
  • n::Int : Number of structural degrees of freedom.

Returns

A StructuralState{T} with all fields initialised to zero.

source
Immersa.StructureBCType
StructureBC{T}

Represents a single structural boundary condition that prescribes the value of a specific degree of freedom as a function of (index, step, time).

Fields

  • index::Int : Global DOF index to constrain.
  • value::FunctionWrapper{T,Tuple{Int,Int,T}} : Callable (index, i, t) -> T returning the prescribed value at step i and time t.

Constructor

StructureBC{T}(index, value)

Wraps value in a FunctionWrapper for type-stable dispatch.

source
Immersa.GeometricNonlinearBodyType
GeometricNonlinearBody{N,T,S,V,B} <: AbstractBody

A deformable body described by a geometrically nonlinear (co-rotational) beam model.

The body is discretised as a chain of 2-node Euler–Bernoulli beam elements in N dimensions, each carrying axial, transverse, and rotational degrees of freedom. An optional prescribed sub-body handles any rigidly attached or prescribed-motion points (e.g. a clamped root or a moving base).

Type parameters

  • N : Spatial dimension (typically 2).
  • T : Scalar type (e.g. Float64).
  • S : Vector type for per-element scalar properties.
  • V : Vector type for per-node position vectors (SVector{N,T}).
  • B : Type of the prescribed sub-body.

Fields

  • xref::V : Reference (undeformed) positions of the deforming nodes.
  • ds0::S : Reference element arc-lengths.
  • m::S : Per-element mass per unit length.
  • kb::S : Per-element bending stiffness.
  • ke::S : Per-element extensional stiffness.
  • bcs::Vector{StructureBC{T}}: Boundary conditions on structural DOFs.
  • prescribed::B : Sub-body for any prescribed-motion points (default NothingBody()).
source
Immersa.deforming_point_rangeFunction
deforming_point_range(body::GeometricNonlinearBody) -> UnitRange{Int}

Return the index range 1:n covering the deforming nodes in the combined point array.

source
Immersa.structure_var_countFunction
structure_var_count(body) -> Int

Return the number of structural degrees of freedom for body.

  • For any AbstractPrescribedBody, returns 0 (no structural DOFs).
  • For a 2-D GeometricNonlinearBody, returns 3n where n is the number of deforming nodes (two translations + one rotation per node).
source
Immersa.prescribed_point_rangeFunction
prescribed_point_range(body::GeometricNonlinearBody) -> UnitRange{Int}

Return the index range covering the prescribed (non-deforming) points that follow the deforming nodes in the combined point array.

source
Immersa.GeometricNonlinearBodyOperatorsType
GeometricNonlinearBodyOperators{TM,TK,TQ,S,TKh1,TKh2}

Mutable container holding the finite-element operators for a GeometricNonlinearBody. These are assembled and updated during the simulation as the body deforms.

Fields

  • M::TM : Global consistent mass matrix.
  • K::TK : Global tangent stiffness matrix (material + geometric).
  • Q::TQ : Force-projection operator (fluid → structure), implemented as a LinearMap.
  • β0::S : Reference element orientations (angles of undeformed elements).
  • Khat_temp::TKh1 : Workspace for the effective stiffness matrix K + (4/dt²) M.
  • Khat::TKh2 : LU factorisation of Khat_temp, used in the Newmark solve.
  • Fint::S : Internal (elastic) force vector in global coordinates.
  • f_temp::S : Temporary vector for force transformations.
source
Immersa.structural_operatorsFunction
structural_operators(backend::CPU, body::GeometricNonlinearBody{N,T}) -> GeometricNonlinearBodyOperators

Allocate and initialise all finite-element operators for body on the given backend.

Computes the reference element orientations β0, allocates the global mass, stiffness, and internal-force arrays, builds the force-projection operator Q, and returns a fully constructed GeometricNonlinearBodyOperators.

Arguments

  • backend::CPU : Computation backend.
  • body::GeometricNonlinearBody{N,T} : The nonlinear structural body.

Returns

A GeometricNonlinearBodyOperators ready for use in FSI time-stepping.

source
Immersa.Q_matrixFunction
Q_matrix(::CPU, body::GeometricNonlinearBody{N,T}) -> LinearMap{T}

Build the force-projection operator Q that maps fluid-frame nodal forces to generalised structural forces using consistent finite-element shape functions.

The operator is returned as a LinearMap (matrix-free, in-place). Rows corresponding to constrained DOFs (boundary conditions) are zeroed.

Arguments

  • ::CPU : Computation backend.
  • body::GeometricNonlinearBody{N,T} : The nonlinear structural body.

Returns

A LinearMap{T} of size n × n where n = structure_var_count(body).

source
Immersa.init_structure_operators!Function
init_structure_operators!(ops, body, points, state, dt)

Perform the initial assembly of all structural operators at the start of a simulation.

Assembles the mass matrix M, tangent stiffness K, effective stiffness , and internal force vector Fint from the current geometry.

Arguments

  • ops::GeometricNonlinearBodyOperators : Operator container to fill.
  • body::GeometricNonlinearBody : The structural body.
  • points::BodyPoints : Current body-point positions and spacings.
  • state::StructuralState : Current structural state.
  • dt : Time step size.
source
Immersa.update_structure_operators!Function
update_structure_operators!(ops, body, points, state, dt)

Re-assemble geometry-dependent structural operators after the body has deformed.

Unlike init_structure_operators!, this skips the mass matrix (which depends only on the reference configuration) and updates only the tangent stiffness K, effective stiffness , and internal forces Fint.

Arguments

  • ops::GeometricNonlinearBodyOperators : Operator container to update.
  • body::GeometricNonlinearBody : The structural body.
  • points::BodyPoints : Current body-point positions and spacings.
  • state::StructuralState : Current structural state.
  • dt : Time step size.
source
Immersa.update_Khat!Function
update_Khat!(ops, dt)

Form and factorise the effective stiffness matrix used in the Newmark time integrator:

K̂ = K + (4 / dt²) M

The result is stored as an in-place LU factorisation in ops.Khat.

Arguments

  • ops::GeometricNonlinearBodyOperators : Operator container (modified in place).
  • dt : Time step size.
source
Immersa.update_M!Function
update_M!(ops, body, points, state)

Assemble the global consistent mass matrix from element contributions.

Each element uses the standard Euler–Bernoulli beam mass matrix (axial + transverse inertia) with mass per unit length m and current element spacing Δs. Rows and columns corresponding to boundary conditions are zeroed.

Arguments

  • ops::GeometricNonlinearBodyOperators : Operator container (modified in place).
  • body::GeometricNonlinearBody : The structural body.
  • points::BodyPoints : Current body-point positions and spacings.
  • state::StructuralState : Current structural state (unused; kept for interface consistency).
source
Immersa.update_K!Function
update_K!(ops, body, points, state)

Assemble the global tangent stiffness matrix (material + geometric contributions).

For each element, the material stiffness is computed from the co-rotational strain–displacement matrix B and the constitutive matrix CL (bending + extensional). The geometric stiffness accounts for the normal force Nf and bending moments Mf1, Mf2 in the deformed configuration. Rows and columns corresponding to boundary conditions are zeroed with a unit diagonal.

Arguments

  • ops::GeometricNonlinearBodyOperators : Operator container (modified in place).
  • body::GeometricNonlinearBody : The structural body.
  • points::BodyPoints : Current body-point positions and spacings.
  • state::StructuralState : Current structural state (provides rotations θ).
source
Immersa.update_Fint!Function
update_Fint!(ops, body, points, state)

Assemble the global internal (elastic) force vector from element contributions.

For each element, the local generalised forces — normal force Nf and bending moments Mf1, Mf2 — are computed in the co-rotational frame and mapped to global DOFs via the strain–displacement matrix B. Constrained DOFs are zeroed.

Arguments

  • ops::GeometricNonlinearBodyOperators : Operator container (modified in place).
  • body::GeometricNonlinearBody : The structural body.
  • points::BodyPoints : Current body-point positions and spacings.
  • state::StructuralState : Current structural state (provides rotations θ).
source
Immersa.structure_to_fluid_displacement!Function
structure_to_fluid_displacement!(x_fluid, x_structure, body, ops)

Map structural DOFs to fluid-frame (immersed-boundary) displacements.

Extracts the translational components from the interleaved structural vector x_structure (which stores [u₁, v₁, θ₁, u₂, v₂, θ₂, …]) and writes them as SVector{N,T} entries into the deforming-point range of x_fluid.

Arguments

  • x_fluid : Fluid-frame displacement or velocity array (modified in place).
  • x_structure : Structural DOF vector.
  • body::GeometricNonlinearBody{N,T} : The structural body.
  • ops::GeometricNonlinearBodyOperators : Structural operators (unused; kept for dispatch).

Returns

The modified x_fluid.

source
Immersa.fluid_to_structure_force!Function
fluid_to_structure_force!(f_structure, f_fluid, body, ops)

Map fluid-frame body forces to generalised structural forces via the projection operator Q.

The SVector{N,T} entries in f_fluid are unpacked into an interleaved temporary vector, then multiplied by ops.Q to produce the structural force vector f_structure.

Arguments

  • f_structure : Output structural force vector (modified in place).
  • f_fluid : Fluid-frame body force array (Vector{SVector{N,T}}).
  • body::GeometricNonlinearBody{N,T} : The structural body (unused; kept for dispatch).
  • ops::GeometricNonlinearBodyOperators : Structural operators (provides Q and workspace).

Returns

The modified f_structure.

source
Immersa.update_structure!Function
update_structure!(points, def, body, ops, i, t)

Synchronise the immersed-boundary point data with the current structural state.

Converts structural displacements and velocities to fluid-frame quantities, updates the deforming node positions (reference + displacement), velocities, and recomputes the element arc-lengths ds from the deformed geometry.

Arguments

  • points::BodyPoints : Body-point data (modified in place).
  • def::StructuralState : Current structural state.
  • body::GeometricNonlinearBody : The structural body.
  • ops::GeometricNonlinearBodyOperators : Structural operators.
  • i : Current time step index.
  • t : Current simulation time.
source
Immersa.update_structure_bc!Function
update_structure_bc!(χ, body, i, t)

Apply structural boundary conditions by overwriting constrained entries in the displacement vector χ with their prescribed values at step i and time t.

Arguments

  • χ : Structural displacement vector (modified in place).
  • body::GeometricNonlinearBody : The structural body (provides bcs).
  • i : Current time step index.
  • t : Current simulation time.
source