Restrictions

ExtendableFEM.jl provides functionality to apply various restrictions to your finite element problems through the AbstractRestriction type system. Restrictions enforce additional constraints on your solution that cannot be easily expressed through standard boundary conditions or the weak formulation.

Built-in Restrictions

Linear Functional Restriction

LinearFunctionalRestriction allows you to constrain a linear functional of a finite element unknown to a specific value. This is useful, for example, for enforcing that the solution has mean zero, or that an integral constraint is satisfied.

Documentation:

ExtendableFEM.LinearFunctionalRestrictionMethod
LinearFunctionalRestriction(
    O::LinearOperator{Tv};
    value
) -> LinearFunctionalRestriction{Int64}

Construct a linear functional restriction for a finite element unknown from a given LinearOperator.

Arguments

  • O::LinearOperator{Tv}: The linear operator representing the functional to be restricted.
  • value::T: The target value for the linear functional (default: 0).
source
ExtendableFEM.MassRestrictionMethod
MassRestriction(
    u::Unknown;
    kernel,
    value,
    operator,
    Tv
) -> LinearFunctionalRestriction{Int64, Float64}

Construct a mass/integral restriction for a finite element unknown.

This restriction enforces that the integral of the given unknown (or of an operator applied to it) attains a given value. Internally, it creates a LinearOperator using the provided kernel and operator, and wraps it in a LinearFunctionalRestriction.

Arguments

  • u::Unknown: The unknown whose mean value is to be restricted.
  • kernel: Kernel function for the linear operator (default: ExtendableFEMBase.constant_one_kernel).
  • value::T: The target integral/mass value (default: 0).
  • operator: Operator to apply to the unknown before restriction (default: Identity).
  • Tv: Value type for the restriction (default: Float64).
source
ExtendableFEM.ZeroMeanValueRestrictionMethod
ZeroMeanValueRestriction(
    u::Unknown;
    kernel,
    operator,
    Tv
) -> LinearFunctionalRestriction{Int64, Float64}

Construct a zero mean value restriction for a finite element unknown.

This restriction enforces that the mean value of the given unknown (or of an operator applied to it) vanishes. Internally, it creates a LinearOperator using the provided kernel and operator, and wraps it in a LinearFunctionalRestriction.

Arguments

  • u::Unknown: The unknown whose mean value is to be restricted.
  • kernel: Kernel function for the linear operator (default: ExtendableFEMBase.constant_one_kernel).
  • operator: Operator to apply to the unknown before restriction (default: Identity).
  • Tv: Value type for the restriction (default: Float64).
source

Example Usage

# Restrict the mean value of unknown u to 0
restriction = ZeroMeanValueRestriction(u)

# Restrict the mean value to 1.0 with a specific operator
restriction = MassRestriction(u, value = 1.0, operator = MyCustomOperator)

# Assign to the problem
assign_restriction!(PD, restriction)

Coupled DOFs Restriction

CoupledDofsRestriction enables coupling between different degrees of freedom (DOFs) in your system. This is particularly useful for implementing periodic boundary conditions or other constraints where DOFs need to be related to each other. Compared to manipulating the system matrix directly via operators (e.g. CombineDofs), CoupledDofsRestriction is much faster for large systems.

Documentation:

ExtendableFEM.CoupledDofsRestrictionMethod
CoupledDofsRestriction(matrix::AbstractMatrix)

Creates an restriction that couples dofs together.

The coupling is stored in a CSC Matrix `matrix`, s.t.,

dofᵢ = Σⱼ Aⱼᵢ dofⱼ (col-wise)

The matrix can be obtained from, e.g., `get_periodic_coupling_matrix`.
source

Example Usage

# Create a coupling matrix for periodic boundary conditions
coupling_matrix = get_periodic_coupling_matrix(...)
restriction = CoupledDofsRestriction(coupling_matrix)
assign_restriction!(PD, restriction)

Boundary Data Restriction

BoundaryDataRestriction enforces prescribed boundary values for a finite element unknown. It can be used for both homogeneous (zero) and non-homogeneous boundary conditions, interpolating the boundary data as needed.

Documentation:

ExtendableFEM.BoundaryDataRestrictionType
BoundaryDataRestriction(
    u::Unknown;
    ...
) -> BoundaryDataRestriction{BOT} where BOT<:(HomogeneousData{Unknown{IT}} where IT)
BoundaryDataRestriction(
    u::Unknown,
    data;
    kwargs...
) -> Union{BoundaryDataRestriction{BOT} where BOT<:(HomogeneousData{Unknown{IT}} where IT), BoundaryDataRestriction{BOT} where BOT<:(InterpolateBoundaryData{Unknown{IT}} where IT)}

Construct a boundary data restriction for an unknown.

This restriction enforces prescribed boundary values for the given unknown, either as homogeneous data (similar to a HomogeneousBoundaryData operator) or interpolated from a provided data function (similar to an InterpolateBoundaryOperator).

Arguments

  • u::Unknown: The unknown whose boundary values are to be restricted.
  • data: A function with signature data!(result, qpinfo) that computes the desired boundary values at each quadrature point. The qpinfo argument provides information such as global coordinates (qpinfo.x). If no

data is given, homogeneous boundary data is assembled

  • kwargs...: Additional keyword arguments passed to the boundary operator constructors.
source

Example Usage

# Homogeneous boundary data (default)
restriction = BoundaryDataRestriction(u)

# Inhomogeneous boundary data using a function
restriction = BoundaryDataRestriction(u, data = (result, qpinfo) -> result .= sin(qpinfo.x[1]))

assign_restriction!(PD, restriction)

AbstractRestriction API

All restrictions are subtypes of AbstractRestriction. The following functions are available for interacting with restrictions:

  • assemble!(restriction, solution, SC; kwargs...): Assemble the internal matrices and vectors for the restriction.
  • restriction_matrix(restriction): Get the restriction matrix.
  • restriction_rhs(restriction): Get the right-hand side vector for the restriction.
  • fixed_dofs(restriction): Get the list of fixed degrees of freedom affected by the restriction.
  • is_timedependent(restriction): Returns whether the restriction is time-dependent.