Assembly Patterns

The definition and assembly of an operator is essentially based on Assembly Patterns and a Assembly Type to choose the geometry information needed for providing quadrature and dof handling. The assembly pattern then basically evaluates the function operators and action for the ansatz and test functions and does the quadrature-weighted accumulation into matrices or vectors that represent the operators.

Each pattern comes with a number of arguments/quantities with associated Function Operators as well as one of the Assembly Type that states whether the form is evaluated over CELLS, FACES order BFACES (see above). Important note: this assembly type is relative to the grid of the first argument of the pattern. If this argument already lives ONFACES and the pattern is also ONFACES, it will ultimatively assemble on the faces of the faces (that are the edges of the grid with these faces). Moreover, patterns can have an Action that allow to make the evaluations parameter-, region- and/or function-dependent. Each pattern then has usually on to three implementation that writes into FEMatrix or FEVector (where e.g. a subset of arguments is fixed) or evaluates the pattern in the given FEVectorBlocks.

The patterns are used to assembly the PDE operators defined in a PDE Description.

GradientRobustMultiPhysics.AssemblyPatternType
mutable struct AssemblyPattern{APT<:AssemblyPatternType, T<:Real, AT<:AssemblyType, Tv<:Real, Ti<:Integer, ActionType<:Union{AbstractAction, GradientRobustMultiPhysics.AbstractNonlinearFormHandler}}

each assembly pattern has one of the assembly pattern types (APT) that trigger different assemblies for the involved finite element spaces, operators and an assigned action. The assembly type (AT) determines if the assembly takes place on cells, faces or edges etc. (relatively to the assembly type of the first argument of the pattern)

source

The following table lists all available assembly patterns, their constuctor names and how they can be used for assembly or evaluations.

AssemblyPatternTypeconstructorevaluateassembles into matrixassembles into vector
APT_ItemIntegratorItemIntegratoryesnono
APT_LinearFormLinearFormnonoyes
APT_BilinearFormBilinearFormnoyesno
APT_NonlinearFormNonlinearFormnoyesyes

Evaluations of the other AssemblyPatterns may be possible in a future update, but currently have to be performed by maintaining a duplicate of the pattern rewritten as an ItemIntegrator or by matrix vector mutliplications.

Constructor details

Below all assembly pattern types, constructor functions and evaluate/assembly functions are detailed. (For more info on the ItemIntegrator also see Item Integrators.)

GradientRobustMultiPhysics.ItemIntegratorFunction
ItemIntegrator(
    operators;
    ...
) -> AssemblyPattern{ItemIntegrator, Float64, ON_CELLS, Float64, Int32, NoAction}
ItemIntegrator(
    operators,
    action;
    T,
    AT,
    regions,
    name
) -> AssemblyPattern{ItemIntegrator, Float64, ON_CELLS, Float64, Int32}

Creates an ItemIntegrator assembly pattern based on:

  • operators : operators that should be evaluated for the coressponding FESpace (last one refers to test function)
  • action : an Action with kernel of interface (result, input, kwargs) that takes input (= all but last operator evaluations) and computes result to be dot-producted with test function evaluation (if no action is specified, the full input vector is dot-producted with the test function operator evaluation)

Optional arguments:

  • T : expected NumberType for evaluation output
  • AT : specifies on which entities of the grid the ItemINtegrator is evaluated
  • regions : specifies in which regions the operator should assemble, default [0] means all regions
  • name : name for this LinearForm that is used in print messages
source
GradientRobustMultiPhysics.L2DifferenceIntegratorMethod
L2DifferenceIntegrator(
    ncomponents::Int64,
    operator;
    AT,
    T,
    name,
    quadorder,
    regions
) -> AssemblyPattern{ItemIntegrator, T, AT, Float64, Int32, DefaultUserAction{T1, Ti, dx, dt, di, dl, ndim, KernelType}} where {T<:Real, AT<:AssemblyType, T1, Ti, dx, dt, di, dl, ndim, KernelType}

Creates an ItemIntegrator that computes the L2 norm difference between two arguments evalauted with the same operator (or with different operators if operator is an array) where ncomponents is the expected length of each operator evaluation. Note that all arguments in an evaluation call need to be defined on the same grid !

source
GradientRobustMultiPhysics.L2ErrorIntegratorFunction
L2ErrorIntegrator(
    compare_data::GradientRobustMultiPhysics.AbstractUserDataType;
    ...
) -> AssemblyPattern{ItemIntegrator, T, AT, Float64, Int32, DefaultUserAction{T1, Ti, dx, dt, di, dl, ndim, KernelType}} where {T<:Real, AT<:AssemblyType, T1, Ti, dx, dt, di, dl, ndim, KernelType}
L2ErrorIntegrator(
    compare_data::GradientRobustMultiPhysics.AbstractUserDataType,
    operator;
    T,
    quadorder,
    name,
    AT,
    factor,
    regions,
    time
) -> AssemblyPattern{ItemIntegrator, T, AT, Float64, Int32, DefaultUserAction{T1, Ti, dx, dt, di, dl, ndim, KernelType}} where {T<:Real, AT<:AssemblyType, T1, Ti, dx, dt, di, dl, ndim, KernelType}

Creates an ItemIntegrator that compares discrete FEVectorBlock operator-evaluations against the given comparedata and returns the L2-error || comparedata(x) - factor*discrete(x) ||. If quadorder is left on "auto" two times the quadorder of the data is used in the evaluation.

source
GradientRobustMultiPhysics.L2NormIntegratorMethod
L2NormIntegrator(
    ncomponents::Int64,
    operator;
    T,
    AT,
    name,
    quadorder,
    regions
) -> AssemblyPattern{ItemIntegrator, T, AT, Float64, Int32, DefaultUserAction{T1, Ti, dx, dt, di, dl, ndim, KernelType}} where {T<:Real, AT<:AssemblyType, T1, Ti, dx, dt, di, dl, ndim, KernelType}

Creates an ItemIntegrator that computes the L2 norm of an operator evaluation where ncomponents is the expected length of the operator evaluation.

source
GradientRobustMultiPhysics.evaluate!Method
function evaluate!(
    b::AbstractArray{T,2},
    AP::AssemblyPattern{APT,T,AT},
    FEB::Union{<:FEVector{T,Tv,Ti},<:FEVectorBlock{T,Tv,Ti},Array{<:FEVectorBlock{T,Tv,Ti},1}};
    skip_preps::Bool = false) where {APT <: APT_ItemIntegrator, T<: Real, AT <: AssemblyType, Tv, Ti}

Evaluation of an ItemIntegrator assembly pattern with given FEVectorBlock or FEVector FEB into given two-dimensional Array b.

source
GradientRobustMultiPhysics.evaluateMethod
function evaluate(
    AP::AssemblyPattern{APT,T,AT},
    FEB::Union{<:FEVector{T,Tv,Ti},<:FEVectorBlock{T,Tv,Ti},Array{<:FEVectorBlock{T,Tv,Ti},1}};
    skip_preps::Bool = false) where {APT <: APT_ItemIntegrator, T<: Real, AT <: AssemblyType, Tv, Ti}

Evaluation of an ItemIntegrator assembly pattern with given FEVectorBlock or FEVector FEB, only returns accumulation over all items.

source
GradientRobustMultiPhysics.DiscreteBilinearFormFunction
DiscreteBilinearForm(
    operators,
    FES;
    ...
) -> AssemblyPattern{BilinearForm, Float64, ON_CELLS, _A, _B, NoAction} where {_A<:Real, _B<:Integer}
DiscreteBilinearForm(
    operators,
    FES,
    action;
    T,
    AT,
    name,
    regions,
    apply_action_to
) -> AssemblyPattern{BilinearForm, Float64, ON_CELLS}

Creates a (discrete) BilinearForm assembly pattern based on:

  • operators : operators that should be evaluated for the coressponding FESpace (last two refer to ansatz and test function)
  • FES : FESpaces for each operator (last two refer to ansatz and test function)
  • action : an Action with kernel of interface (result, input, kwargs) that takes input (= all but last operator evaluations) and computes result to be dot-producted with test function evaluation (if no action is specified, the full input vector is dot-producted with the test function operator evaluation)

Optional arguments:

  • applyactionto : specifies which of the two linear arguments is part of the action input ([1] = ansatz, [2] = test)
  • regions : specifies in which regions the operator should assemble, default [0] means all regions
  • name : name for this LinearForm that is used in print messages
  • AT : specifies on which entities of the grid the LinearForm is assembled (default: ON_CELLS)
source
GradientRobustMultiPhysics.DiscreteLinearFormMethod
DiscreteLinearForm(
    operators,
    FES::Array{<:FESpace{Tv, Ti}, 1};
    ...
) -> AssemblyPattern{LinearForm, Float64, ON_CELLS, _A, _B, NoAction} where {_A<:Real, _B<:Integer}
DiscreteLinearForm(
    operators,
    FES::Array{<:FESpace{Tv, Ti}, 1},
    action;
    T,
    AT,
    regions,
    name
) -> AssemblyPattern{LinearForm, Float64, ON_CELLS}

Creates a (discrete) LinearForm assembly pattern based on:

  • operators : operators that should be evaluated for the coressponding FESpace (last one refers to test function)
  • FES : FESpaces for each operator (last one refers to test function)
  • action : an Action with kernel of interface (result, input, kwargs) that takes input (= all but last operator evaluations) and computes result to be dot-producted with test function evaluation (if no action is specified, the full input vector is dot-producted with the test function operator evaluation)

Optional arguments:

  • regions : specifies in which regions the operator should assemble, default [0] means all regions
  • name : name for this LinearForm that is used in print messages
  • AT : specifies on which entities of the grid the LinearForm is assembled (default: ON_CELLS)
source
GradientRobustMultiPhysics.DiscreteLumpedBilinearFormFunction
DiscreteLumpedBilinearForm(
    operators,
    FES;
    ...
) -> AssemblyPattern{LumpedBilinearForm, Float64, ON_CELLS, _A, _B, NoAction} where {_A<:Real, _B<:Integer}
DiscreteLumpedBilinearForm(
    operators,
    FES,
    action;
    T,
    AT,
    name,
    regions,
    apply_action_to
) -> AssemblyPattern{LumpedBilinearForm, Float64, ON_CELLS}

Creates a (discrete) LumpedBilinearForm assembly pattern. For more details see BilinearForm constructor.

source
GradientRobustMultiPhysics.DiscreteSymmetricBilinearFormFunction
DiscreteSymmetricBilinearForm(
    operators,
    FES;
    ...
) -> AssemblyPattern{SymmetricBilinearForm, Float64, ON_CELLS, _A, _B, NoAction} where {_A<:Real, _B<:Integer}
DiscreteSymmetricBilinearForm(
    operators,
    FES,
    action;
    T,
    AT,
    name,
    regions,
    apply_action_to
) -> AssemblyPattern{SymmetricBilinearForm, Float64, ON_CELLS}

Creates a (discrete) SymmetricBilinearForm assembly pattern. For more details see BilinearForm constructor.

source
GradientRobustMultiPhysics.assemble!Method
assemble!(
    b::Union{AbstractArray{T,1},AbstractArray{T,2}},    # target vector/matrix
    AP::AssemblyPattern{APT,T,AT};                      # LinearForm pattern
    factor = 1)                                         # factor that is multiplied
    where {APT <: APT_LinearForm, T, AT}

Assembly of a LinearForm pattern AP into a vector or matrix (if action is vetor-valued).

source
GradientRobustMultiPhysics.assemble!Method
assemble!(
    A::AbstractArray{T,2},                  # target matrix
    AP::AssemblyPattern{APT,T,AT};          # BilinearForm Pattern
    apply_action_to::Int = 1,               # action is applied to which argument?
    factor = 1,                             # factor that is multiplied
    transposed_assembly::Bool = false,      # transpose result?
    transpose_copy = Nothing)               # copy a transposed block to this matrix
    where {APT <: APT_BilinearForm, T, AT}

Assembly of a BilinearForm BLF into given two-dimensional AbstractArray (e.g. FEMatrixBlock or a ExtendableSparseMatrix).

source
GradientRobustMultiPhysics.full_assemble!Method
full_assemble!(
    A::AbstractArray{T,2},                 # target matrix
    b::AbstractArray{T,1},                 # target rhs
    AP::AssemblyPattern{APT,T,AT};         # NonlinearForm pattern
    FEB::Array{<:FEVectorBlock,1};         # coefficients of current solution for each operator
    factor = 1,                            # factor that is multiplied
    transposed_assembly::Bool = false)     # transpose result?
    where {APT <: APT_NonlinearForm, T, AT}

Assembly (of Newton terms) of a NonlinearForm assembly pattern (assembles both matrix and rhs!).

source