Multifield API
Weak-form DSL and multifield assembly interface.
Weak-Form DSL Operators
LowLevelFEM.Grad — FunctionGrad(P)Create a weak-form DSL gradient operator applied to P.
Arguments
P: Field descriptor (Problem) used in the weak form.
Returns
OpApplied: Operator application object for∇.
Example
K = ∫(Grad(Pu) ⋅ Grad(Pu); Ω="solid")LowLevelFEM.Div — FunctionDiv(P)Create a weak-form DSL divergence operator applied to P.
Arguments
P: Field descriptor (Problem) used in the weak form.
Returns
OpApplied: Operator application object for∇⋅.
Example
A = ∫(Div(Pu) ⋅ Div(Pu); Ω="domain")LowLevelFEM.Curl — FunctionCurl(P)Create a weak-form DSL curl operator applied to P.
Arguments
P: Field descriptor (Problem) used in the weak form.
Returns
OpApplied: Operator application object for curl.
Example
A = ∫(Curl(Pu) ⋅ Curl(Pu); Ω="domain")LowLevelFEM.SymGrad — FunctionSymGrad(P)Create a weak-form DSL symmetric-gradient operator applied to P.
Arguments
P: Field descriptor (Problem) used in the weak form.
Returns
OpApplied: Operator application object forε(u).
Example
K = ∫(SymGrad(Pu) ⋅ C ⋅ SymGrad(Pu); Ω="solid")LowLevelFEM.Id — FunctionId(P)Create a weak-form DSL identity operator applied to P.
Arguments
P: Field descriptor (Problem) used in the weak form.
Returns
OpApplied: Operator application object for identity mapping.
Example
M = ∫(Id(Pu) ⋅ Id(Pu); Ω="solid")LowLevelFEM.TensorDiv — FunctionTensorDiv(P)Create a weak-form DSL tensor-divergence operator applied to P.
Arguments
P: Field descriptor (Problem) used in the weak form.
Returns
OpApplied: Operator application object for tensor divergence.
Example
A = ∫(TensorDiv(Pσ) ⋅ TensorDiv(Pσ); Ω="solid")LowLevelFEM.Adv — FunctionAdv(P)Create a weak-form DSL advection operator applied to P.
Arguments
P: Field descriptor (Problem) used in the weak form.
Returns
OpApplied: Operator application object for advection terms.
Example
A = ∫(Adv(Pu) ⋅ Id(Pu); Ω="domain")LowLevelFEM.AxialGrad — FunctionAxialGrad(P)Create a weak-form DSL axial gradient operator applied to P.
Arguments
P: Field descriptor (Problem) used in the weak form.
Returns
OpApplied: Operator application object representingt ⋅ ∇u.
Example
```julia K = ∫(AxialGrad(Pu) ⋅ (E*A) ⋅ AxialGrad(Pu); Ω="truss")
LowLevelFEM.ε — FunctionSymGrad(P)Create a weak-form DSL symmetric-gradient operator applied to P.
Arguments
P: Field descriptor (Problem) used in the weak form.
Returns
OpApplied: Operator application object forε(u).
Example
K = ∫(SymGrad(Pu) ⋅ C ⋅ SymGrad(Pu); Ω="solid")Weak-Form Integration
LowLevelFEM.∫ — Function∫(expr::WeakExpr; Ω=nothing, Γ=nothing, weight=nothing)Assemble a finite element operator from a weak-form expression.
Examples
Diffusion
K = ∫( Grad(Pu) ⋅ Grad(Pu); Ω="domain")Elasticity
K = ∫( SymGrad(Pu) ⋅ C ⋅ SymGrad(Pu); Ω="solid")Tensor chain
K = ∫( Grad(Pu) ⋅ F' ⋅ S ⋅ F ⋅ Grad(Pu); Ω="solid")Elastic support
K = ∫( Id(Pu) ⋅ [kx 0; 0 ky] ⋅ Id(Pu); Γ="boundary")or
K = ∫( Pu ⋅ [kx 0; 0 ky] ⋅ Pu; Γ="boundary")Mixed formulation
A = ∫( Div(Pu) ⋅ Pp )
B = ∫( Pp ⋅ Div(Pu) )Linear form
f = ∫(Pu ⋅ g)With operator chain
f = ∫(Pu ⋅ A ⋅ g)With coefficient
f = ∫(PT ⋅ PT * h, Γ="right")Arguments
expr
Weak-form expression composed of operators and coefficients.
Keyword arguments
Ω
Volume physical group name.
Γ
Boundary physical group name.
Returns
SystemMatrix or ScalarField, VectorField, TensorField
LowLevelFEM.∫Ω — Function∫Ω(name, expr)Convenience wrapper for volume integration on physical group name.
Arguments
name: Gmsh physical group name used as domainΩ.expr::WeakExpr: Weak-form expression to assemble.
Returns
SystemMatrix: Assembled matrix over the selected volume.
Example
K = ∫Ω("solid", Grad(Pu) ⋅ Grad(Pu))LowLevelFEM.∫Γ — Function∫Γ(name, expr)Convenience wrapper for boundary integration on physical group name.
Arguments
name: Gmsh physical group name used as boundaryΓ.expr::WeakExpr: Weak-form expression to assemble.
Returns
SystemMatrix: Assembled matrix over the selected boundary.
Example
KΓ = ∫Γ("loaded_boundary", Id(Pu) ⋅ Id(Pu))Multifield Solver
LowLevelFEM.solveField — FunctionsolveField(K::SystemMatrix,
f::Union{ScalarField,VectorField,TensorField};
support::Vector{BoundaryCondition}=BoundaryCondition[],
iterative=false,
reltol::Real = sqrt(eps()),
maxiter::Int = K.model.non * K.model.dim,
preconditioner = Identity(),
ordering=true)multifield version:
solveField(K::SystemMatrix,
F::SystemVector;
support::Vector{BoundaryCondition}=BoundaryCondition[])Solves a linear static field problem with Dirichlet boundary conditions.
The algebraic system
K * u = f
is solved for the unknown field u, where prescribed values are enforced solver-side via boundary conditions.
The solution is obtained by partitioning the degrees of freedom into free and constrained sets. On constrained DOFs the values prescribed by support override the solution, while on free DOFs the reduced system is solved:
K_ff * u_f = f_f − K_fc * u_c
The function supports both direct and iterative linear solvers and works uniformly for ScalarField and VectorField unknowns.
Arguments
K::SystemMatrixGlobal system matrix.f::Union{ScalarField,VectorField}Right-hand-side vector.support::Vector{BoundaryCondition}(keyword, default =BoundaryCondition[]) Dirichlet-type boundary conditions.iterative::Bool(keyword, default =false) Iftrue, use an iterative solver (conjugate gradient).reltol::Real(keyword, default =sqrt(eps())) Relative tolerance for the iterative solver.maxiter::Int(keyword, default =K.model.non * K.model.dim) Maximum number of iterations for the iterative solver.preconditioner(keyword, default =Identity()) Preconditioner for the iterative solver.ordering::Bool(keyword, default =true) Iffalse, use an explicit LU factorization without reordering.
Returns
u::Union{ScalarField,VectorField}Solution field with prescribed values enforced on constrained DOFs.
Notes
- Boundary conditions are applied inside the solver; the input field
fis not modified. - The algorithm itself is agnostic to the physical meaning of the field (scalar, vector, tensor), as long as
Kandfare dimensionally consistent. - When
iterative = true, the system is solved using conjugate gradient on the reduced matrixK_ff.
solveField(K::SystemMatrix,
F::SystemVector;
support::Vector{BoundaryCondition}=BoundaryCondition[])single field version:
solveField(K::SystemMatrix,
f::Union{ScalarField,VectorField,TensorField};
support::Vector{BoundaryCondition}=BoundaryCondition[],
iterative=false,
reltol::Real = sqrt(eps()),
maxiter::Int = K.model.non * K.model.dim,
preconditioner = Identity(),
ordering=true)Solve the linear system
K * x = Ffor a multifield finite element problem.
Dirichlet boundary conditions are imposed through BoundaryCondition objects.
Arguments
K
Global system matrix (SystemMatrix).
F
Global RHS vector (SystemVector or ScalarField or VectorField)
Keyword arguments
support
Vector of boundary conditions.
Returns
Tuple of fields corresponding to the Problems stored in K.
Example
(u, p) = solveField(K, F)LowLevelFEM.solveEigenFields — FunctionsolveEigenFields(K::SystemMatrix, M::SystemMatrix; n=6, fmin=0.0)Solve eigenproblem for multifield system and return field-wise Eigen objects.
Usage: umodes, pmodes = solveEigenFields(K, M)