ExtendableFEM.BilinearOperatorType
function BilinearOperator(
	A::AbstractMatrix,
	u_test,
	u_ansatz = u_test;
	kwargs...)

Generates a bilinear form from a user-provided matrix, which can be a sparse matrix or a FEMatrix with multiple blocks. The arguments utest and uansatz specify where to put the (blocks of the) matrix in the system.

ExtendableFEM.BilinearOperatorType
function BilinearOperator(
	[kernel!::Function],
	oa_test::Array{<:Tuple{Union{Unknown,Int}, DataType},1},
	oa_ansatz::Array{<:Tuple{Union{Unknown,Int}, DataType},1} = oa_test;
	kwargs...)

Generates a bilinear form that evaluates the vector product of the operator evaluation(s) of the test function(s) with the operator evaluation(s) of the ansatz function(s). If a function is provided in the first argument, the ansatz function evaluations can be customized by the kernel function and its result vector is then used in a dot product with the test function evaluations. In this case the header of the kernel functions needs to be conform to the interface

kernel!(result, eval_ansatz, qpinfo)

where qpinfo allows to access information at the current quadrature point.

Operator evaluations are tuples that pair an unknown identifier or integer with a Function operator.

Example: BilinearOperator([grad(1)], [grad(1)]; kwargs...) generates a weak Laplace operator.

Keyword arguments:

  • lump: diagonal lumping (=0 no lumping, =1 only keep diagonal entry, =2 accumulate full row to diagonal). Default: 0

  • factor: factor that should be multiplied during assembly. Default: 1

  • verbosity: verbosity level. Default: 0

  • time_dependent: operator is time-dependent ?. Default: false

  • name: name for operator used in printouts. Default: ''BilinearOperator''

  • entities: assemble operator on these grid entities (default = ONCELLS). Default: ONCELLS

  • quadorder: quadrature order. Default: ''auto''

  • bonus_quadorder: additional quadrature order added to quadorder. Default: 0

  • entry_tolerance: threshold to add entry to sparse matrix. Default: 0

  • regions: subset of regions where operator should be assembly only. Default: Any[]

  • params: array of parameters that should be made available in qpinfo argument of kernel function. Default: nothing

  • parallel_groups: assemble operator in parallel using CellAssemblyGroups. Default: false

  • transposed_copy: assemble a transposed copy of that operator into the transposed matrix block(s), 0 = no, 1 = symmetric, -1 = skew-symmetric. Default: 0

  • usesparsitypattern: read sparsity pattern of jacobian of kernel to find out which components couple. Default: ''auto''

  • store: store matrix separately (and copy from there when reassembly is triggered). Default: false

ExtendableFEM.BilinearOperatorMethod
function BilinearOperator(
	kernel::Function,
	oa_test::Array{<:Tuple{Union{Unknown,Int}, DataType},1},
	oa_ansatz::Array{<:Tuple{Union{Unknown,Int}, DataType},1},
	oa_args::Array{<:Tuple{Union{Unknown,Int}, DataType},1};
	kwargs...)

Generates a nonlinear bilinear form that evaluates a kernel function that depends on the operator evaluation(s) of the ansatz function(s) and the operator evaluations of the current solution. The result of the kernel function is used in a vector product with the operator evaluation(s) of the test function(s). Hence, this can be used as a linearization of a nonlinear operator. The header of the kernel functions needs to be conform to the interface

kernel!(result, eval_ansatz, eval_args, qpinfo)

where qpinfo allows to access information at the current quadrature point.

Operator evaluations are tuples that pair an unknown identifier or integer with a Function operator.

Example: BilinearOperator([grad(1)], [grad(1)]; kwargs...) generates a weak Laplace operator.

Keyword arguments:

  • lump: diagonal lumping (=0 no lumping, =1 only keep diagonal entry, =2 accumulate full row to diagonal). Default: 0

  • factor: factor that should be multiplied during assembly. Default: 1

  • verbosity: verbosity level. Default: 0

  • time_dependent: operator is time-dependent ?. Default: false

  • name: name for operator used in printouts. Default: ''BilinearOperator''

  • entities: assemble operator on these grid entities (default = ONCELLS). Default: ONCELLS

  • quadorder: quadrature order. Default: ''auto''

  • bonus_quadorder: additional quadrature order added to quadorder. Default: 0

  • entry_tolerance: threshold to add entry to sparse matrix. Default: 0

  • regions: subset of regions where operator should be assembly only. Default: Any[]

  • params: array of parameters that should be made available in qpinfo argument of kernel function. Default: nothing

  • parallel_groups: assemble operator in parallel using CellAssemblyGroups. Default: false

  • transposed_copy: assemble a transposed copy of that operator into the transposed matrix block(s), 0 = no, 1 = symmetric, -1 = skew-symmetric. Default: 0

  • usesparsitypattern: read sparsity pattern of jacobian of kernel to find out which components couple. Default: ''auto''

  • store: store matrix separately (and copy from there when reassembly is triggered). Default: false

ExtendableFEM.BilinearOperatorDGType
function BilinearOperatorDG(
	[kernel!::Function],
	oa_test::Array{<:Tuple{Union{Unknown,Int}, DataType},1},
	oa_ansatz::Array{<:Tuple{Union{Unknown,Int}, DataType},1} = oa_test;
	kwargs...)

Generates a bilinear form that evaluates the vector product of the (discontinuous) operator evaluation(s) of the test function(s) with the (discontinuous) operator evaluation(s) of the ansatz function(s). If a function is provided in the first argument, the ansatz function evaluations can be customized by the kernel function and its result vector is then used in a dot product with the test function evaluations. In this case the header of the kernel functions needs to be conform to the interface

kernel!(result, eval_ansatz, qpinfo)

where qpinfo allows to access information at the current quadrature point.

Operator evaluations are tuples that pair an unknown identifier or integer with a Function operator.

Example: BilinearOperatorDG([jump(grad(1))], [jump(grad(1))]; kwargs...) generates an interior penalty stabilisation.

Keyword arguments:

  • lump: lump the operator (= only assemble the diagonal). Default: false

  • factor: factor that should be multiplied during assembly. Default: 1

  • callback!: function with interface (A, b, sol) that is called in each assembly step. Default: nothing

  • time_dependent: operator is time-dependent ?. Default: false

  • name: name for operator used in printouts. Default: ''BilinearOperatorDG''

  • entities: assemble operator on these grid entities (default = ONFACES). Default: ONFACES

  • quadorder: quadrature order. Default: ''auto''

  • bonus_quadorder: additional quadrature order added to quadorder. Default: 0

  • entry_tolerance: threshold to add entry to sparse matrix. Default: 0

  • verbosity: verbosity level. Default: 0

  • regions: subset of regions where operator should be assembly only. Default: Any[]

  • params: array of parameters that should be made available in qpinfo argument of kernel function. Default: nothing

  • parallel_groups: assemble operator in parallel using CellAssemblyGroups. Default: false

  • transposed_copy: assemble a transposed copy of that operator into the transposed matrix block(s), 0 = no, 1 = symmetric, -1 = skew-symmetric. Default: 0

  • usesparsitypattern: read sparsity pattern of jacobian of kernel to find out which components couple. Default: ''auto''

  • store: store matrix separately (and copy from there when reassembly is triggered). Default: false

ExtendableFEM.BilinearOperatorDGMethod
function BilinearOperatorDG(
	kernel::Function,
	oa_test::Array{<:Tuple{Union{Unknown,Int}, DataType},1},
	oa_ansatz::Array{<:Tuple{Union{Unknown,Int}, DataType},1},
	oa_args::Array{<:Tuple{Union{Unknown,Int}, DataType},1};
	kwargs...)

Generates a nonlinear bilinear form that evaluates a kernel function that depends on the (discontinuou) operator evaluation(s) of the ansatz function(s) and the (discontinuous) operator evaluations of the current solution. The result of the kernel function is used in a vector product with the operator evaluation(s) of the test function(s). Hence, this can be used as a linearization of a nonlinear operator. The header of the kernel functions needs to be conform to the interface

kernel!(result, eval_ansatz, eval_args, qpinfo)

where qpinfo allows to access information at the current quadrature point.

Operator evaluations are tuples that pair an unknown identifier or integer with a Function operator.

Keyword arguments:

  • lump: diagonal lumping (=0 no lumping, =1 only keep diagonal entry, =2 accumulate full row to diagonal). Default: 0

  • factor: factor that should be multiplied during assembly. Default: 1

  • verbosity: verbosity level. Default: 0

  • time_dependent: operator is time-dependent ?. Default: false

  • name: name for operator used in printouts. Default: ''BilinearOperator''

  • entities: assemble operator on these grid entities (default = ONCELLS). Default: ONCELLS

  • quadorder: quadrature order. Default: ''auto''

  • bonus_quadorder: additional quadrature order added to quadorder. Default: 0

  • entry_tolerance: threshold to add entry to sparse matrix. Default: 0

  • regions: subset of regions where operator should be assembly only. Default: Any[]

  • params: array of parameters that should be made available in qpinfo argument of kernel function. Default: nothing

  • parallel_groups: assemble operator in parallel using CellAssemblyGroups. Default: false

  • transposed_copy: assemble a transposed copy of that operator into the transposed matrix block(s), 0 = no, 1 = symmetric, -1 = skew-symmetric. Default: 0

  • usesparsitypattern: read sparsity pattern of jacobian of kernel to find out which components couple. Default: ''auto''

  • store: store matrix separately (and copy from there when reassembly is triggered). Default: false

ExtendableFEM.CallbackOperatorType
function CallbackOperator(
	callback!::Function,
	u_args = [];
	kwargs...)

Generates an operator that simply passes the matrix and rhs to a user-specified call back function. The callback function needs to be conform to the interface

callback!(A, b, args; assemble_matrix = true, assemble_rhs = true, time = 0, kwargs...)

The u_args argument can be used to specify the arguments of the solution that should be passed as args (a vector of FEVectorBlocks) to the callback.

Keyword arguments:

  • modifies_matrix: callback function modifies the matrix?. Default: true

  • time_dependent: operator is time-dependent ?. Default: false

  • name: name for operator used in printouts. Default: ''CallbackOperator''

  • linearizeddependencies: [uansatz, u_test] when linearized. Default: auto

  • modifies_rhs: callback function modifies the rhs?. Default: true

  • verbosity: verbosity level. Default: 0

  • store: store matrix and rhs separately (and copy from there when reassembly is triggered). Default: false

ExtendableFEM.CombineDofsType
function CombineDofs(uX, uY, dofsX, dofsY, factors; kwargs...)

When assembled, the dofsX of the unknown uX will be coupled with the dofsY of uY, e.g., for periodic boundary conditions.

Keyword arguments:

  • verbosity: verbosity level. Default: 0

  • penalty: penalty for fixed degrees of freedom. Default: 1.0e30

ExtendableFEM.FixDofsMethod
function FixDofs(u; vals = [], dofs = [], kwargs...)

When assembled, all specified dofs of the unknown u will be penalized to the specified values.

Keyword arguments:

  • verbosity: verbosity level. Default: 0

  • name: name for operator used in printouts. Default: ''FixDofs''

  • penalty: penalty for fixed degrees of freedom. Default: 1.0e30

ExtendableFEM.HomogeneousDataMethod
function HomogeneousData(u; entities = ON_CELLS, kwargs...)

When assembled, the unknown u of the Problem will be penalized to zero on the specifies entities and entity regions (via kwargs).

Keyword arguments:

  • value: constant value of the data. Default: 0

  • regions: subset of regions where operator should be assembly only. Default: Any[]

  • verbosity: verbosity level. Default: 0

  • name: name for operator used in printouts. Default: ''HomogeneousData''

  • penalty: penalty for fixed degrees of freedom. Default: 1.0e30

  • mask: array of zeros/ones to set which components should be set by the operator (only works with componentwise dofs, add a 1 or 0 to mask additional dofs). Default: Any[]

ExtendableFEM.InterpolateBoundaryDataType
function InterpolateBoundaryData(u, data!::Function; kwargs...)

When assembled, the unknown u of the Problem will be penalized to match the standard interpolation of the provided data! function. The header of this function needs to be conform to the interface

data!(result, qpinfo)

where qpinfo allows to access information at the current quadrature point, e.g. qpinfo.x provides the global coordinates of the quadrature/evaluation point.

Keyword arguments:

  • params: array of parameters that should be made available in qpinfo argument of kernel function. Default: nothing

  • regions: subset of regions where operator should be assembly only. Default: Any[]

  • name: name for operator used in printouts. Default: ''BoundaryData''

  • bonus_quadorder: additional quadrature order added to the quadorder chosen by the interpolator. Default: 0

  • penalty: penalty for fixed degrees of freedom. Default: 1.0e30

  • verbosity: verbosity level. Default: 0

  • plot: plot unicode plot of boundary data into terminal when assembled. Default: false

ExtendableFEM.ItemIntegratorMethod
function ItemIntegrator(
	[kernel!::Function],
	oa_args::Array{<:Tuple{Union{Unknown,Int}, DataType},1};
	kwargs...)

Generates an ItemIntegrator that evaluates the specified operator evaluations, puts it into the kernel function and integrates the results over the entities (see kwargs). If no kernel is given, the arguments are integrated directly. If a kernel is provided it has be conform to the interface

kernel!(result, eval_args, qpinfo)

where qpinfo allows to access information at the current quadrature point. Additionally the length of the result needs to be specified via the kwargs.

Evaluation can be triggered via the evaluate function.

Operator evaluations are tuples that pair an unknown identifier or integer with a Function operator.

Keyword arguments:

  • piecewise: returns piecewise integrations, otherwise a global integration. Default: true

  • factor: factor that should be multiplied during assembly. Default: 1

  • resultdim: dimension of result field (default = length of arguments). Default: 0

  • params: array of parameters that should be made available in qpinfo argument of kernel function. Default: nothing

  • name: name for operator used in printouts. Default: ''ItemIntegrator''

  • entities: assemble operator on these grid entities (default = ONCELLS). Default: ONCELLS

  • quadorder: quadrature order. Default: ''auto''

  • bonus_quadorder: additional quadrature order added to quadorder. Default: 0

  • verbosity: verbosity level. Default: 0

  • regions: subset of regions where the item integrator should be evaluated. Default: Any[]

ExtendableFEM.ItemIntegratorDGMethod
function ItemIntegratorDG(
	kernel::Function,
	oa_args::Array{<:Tuple{Union{Unknown,Int}, DataType},1};
	kwargs...)

Generates an ItemIntegrator that evaluates the specified (discontinuous) operator evaluations, puts it into the kernel function and integrates the results over the entities (see kwargs) along cell boundaries. If no kernel is given, the arguments are integrated directly. If a kernel is provided it has be conform to the interface

kernel!(result, eval_args, qpinfo)

where qpinfo allows to access information at the current quadrature point. Additionally the length of the result needs to be specified via the kwargs.

Evaluation can be triggered via the evaluate function.

Operator evaluations are tuples that pair an unknown identifier or integer with a Function operator.

Keyword arguments:

  • piecewise: returns piecewise integrations, otherwise a global integration. Default: true

  • factor: factor that should be multiplied during assembly. Default: 1

  • resultdim: dimension of result field (default = length of arguments). Default: 0

  • params: array of parameters that should be made available in qpinfo argument of kernel function. Default: nothing

  • name: name for operator used in printouts. Default: ''ItemIntegratorDG''

  • entities: assemble operator on these grid entities (default = ONCELLS). Default: ONFACES

  • quadorder: quadrature order. Default: ''auto''

  • bonus_quadorder: additional quadrature order added to quadorder. Default: 0

  • verbosity: verbosity level. Default: 0

  • regions: subset of regions where the item integrator should be evaluated. Default: Any[]

ExtendableFEM.LinearOperatorMethod
function LinearOperator(
	A,
	u_test,
	u_args;
	kwargs...)

Generates a linear form from a user-provided matrix A, which can be an AbstractMatrix or a FEMatrix with multiple blocks. The arguments uargs specify which coefficients of the current solution should be multiplied with the matrix and utest specifies where to put the (blocks of the) resulting vector in the system right-hand side.

ExtendableFEM.LinearOperatorMethod
function LinearOperator(
	b,
	u_test;
	kwargs...)

Generates a linear form from a user-provided vector b, which can be an AbstractVector or a FEVector with multiple blocks. The argument u_test specifies where to put the (blocks of the) vector in the system right-hand side.

ExtendableFEM.LinearOperatorMethod
function LinearOperator(
	kernel!::Function,
	oa_test::Array{<:Tuple{Union{Unknown,Int}, DataType},1},
	oa_args::Array{<:Tuple{Union{Unknown,Int}, DataType},1};
	kwargs...)

Generates a nonlinear linear form that evaluates a kernel function that depends on the operator evaluations of the current solution. The result of the kernel function is used in a vector product with the operator evaluation(s) of the test function(s). Hence, this can be used as a linearization of a nonlinear operator. The header of the kernel functions needs to be conform to the interface

kernel!(result, eval_args, qpinfo)

where qpinfo allows to access information at the current quadrature point.

Operator evaluations are tuples that pair an unknown identifier or integer with a Function operator.

Keyword arguments:

  • factor: factor that should be multiplied during assembly. Default: 1

  • verbosity: verbosity level. Default: 0

  • time_dependent: operator is time-dependent ?. Default: false

  • name: name for operator used in printouts. Default: ''LinearOperator''

  • entities: assemble operator on these grid entities (default = ONCELLS). Default: ONCELLS

  • quadorder: quadrature order. Default: ''auto''

  • bonus_quadorder: additional quadrature order added to quadorder. Default: 0

  • regions: subset of regions where operator should be assembly only. Default: Any[]

  • params: array of parameters that should be made available in qpinfo argument of kernel function. Default: nothing

  • parallel_groups: assemble operator in parallel using CellAssemblyGroups. Default: false

  • store: store matrix separately (and copy from there when reassembly is triggered). Default: false

ExtendableFEM.LinearOperatorMethod
function LinearOperator(
	[kernel!::Function],
	oa_test::Array{<:Tuple{Union{Unknown,Int}, DataType},1};
	kwargs...)

Generates a linear form that evaluates, in each quadrature point, the kernel function (if non is provided, a constant function one is used) and computes the vector product of the result with with the operator evaluation(s) of the test function(s). The header of the kernel functions needs to be conform to the interface

kernel!(result, qpinfo)

where qpinfo allows to access information at the current quadrature point, e.g. qpinfo.x are the global coordinates of the quadrature point.

Operator evaluations are tuples that pair an unknown identifier or integer with a Function operator.

Example: LinearOperator(kernel!, [id(1)]; kwargs...) generates the right-hand side for a Poisson problem, where kernel! evaluates the right-hand side.

Keyword arguments:

  • factor: factor that should be multiplied during assembly. Default: 1

  • verbosity: verbosity level. Default: 0

  • time_dependent: operator is time-dependent ?. Default: false

  • name: name for operator used in printouts. Default: ''LinearOperator''

  • entities: assemble operator on these grid entities (default = ONCELLS). Default: ONCELLS

  • quadorder: quadrature order. Default: ''auto''

  • bonus_quadorder: additional quadrature order added to quadorder. Default: 0

  • regions: subset of regions where operator should be assembly only. Default: Any[]

  • params: array of parameters that should be made available in qpinfo argument of kernel function. Default: nothing

  • parallel_groups: assemble operator in parallel using CellAssemblyGroups. Default: false

  • store: store matrix separately (and copy from there when reassembly is triggered). Default: false

ExtendableFEM.LinearOperatorDGMethod
function LinearOperatorDG(
	kernel::Function,
	oa_test::Array{<:Tuple{Union{Unknown,Int}, DataType},1},
	oa_args::Array{<:Tuple{Union{Unknown,Int}, DataType},1};
	kwargs...)

Generates a nonlinear linear form that evaluates a kernel function that depends on the (discontinous) operator evaluations of the current solution. The result of the kernel function is used in a vector product with the operator evaluation(s) of the test function(s). Hence, this can be used as a linearization of a nonlinear operator. The header of the kernel functions needs to be conform to the interface

kernel!(result, eval_args, qpinfo)

where qpinfo allows to access information at the current quadrature point.

Operator evaluations are tuples that pair an unknown identifier or integer with a Function operator.

Keyword arguments:

  • factor: factor that should be multiplied during assembly. Default: 1

  • verbosity: verbosity level. Default: 0

  • time_dependent: operator is time-dependent ?. Default: false

  • name: name for operator used in printouts. Default: ''LinearOperatorDG''

  • entities: assemble operator on these grid entities (default = ONFACES). Default: ONFACES

  • quadorder: quadrature order. Default: ''auto''

  • bonus_quadorder: additional quadrature order added to quadorder. Default: 0

  • entry_tolerance: threshold to add entry to sparse matrix. Default: 0

  • regions: subset of regions where operator should be assembly only. Default: Any[]

  • params: array of parameters that should be made available in qpinfo argument of kernel function. Default: nothing

  • parallel_groups: assemble operator in parallel using CellAssemblyGroups. Default: false

  • store: store matrix separately (and copy from there when reassembly is triggered). Default: false

ExtendableFEM.LinearOperatorDGMethod
function LinearOperatorDG(
	[kernel!::Function],
	oa_test::Array{<:Tuple{Union{Unknown,Int}, DataType},1};
	kwargs...)

Generates a linear form that evaluates, in each quadrature point, the kernel function (if non is provided, a constant function one is used) and computes the vector product of the result with with the (discontinuous) operator evaluation(s) of the test function(s). The header of the kernel functions needs to be conform to the interface

kernel!(result, qpinfo)

where qpinfo allows to access information at the current quadrature point, e.g. qpinfo.x are the global coordinates of the quadrature point.

Keyword arguments:

  • factor: factor that should be multiplied during assembly. Default: 1

  • verbosity: verbosity level. Default: 0

  • time_dependent: operator is time-dependent ?. Default: false

  • name: name for operator used in printouts. Default: ''LinearOperatorDG''

  • entities: assemble operator on these grid entities (default = ONFACES). Default: ONFACES

  • quadorder: quadrature order. Default: ''auto''

  • bonus_quadorder: additional quadrature order added to quadorder. Default: 0

  • entry_tolerance: threshold to add entry to sparse matrix. Default: 0

  • regions: subset of regions where operator should be assembly only. Default: Any[]

  • params: array of parameters that should be made available in qpinfo argument of kernel function. Default: nothing

  • parallel_groups: assemble operator in parallel using CellAssemblyGroups. Default: false

  • store: store matrix separately (and copy from there when reassembly is triggered). Default: false

ExtendableFEM.NonlinearOperatorType
function NonlinearOperator(
	[kernel!::Function],
	oa_test::Array{<:Tuple{Union{Unknown,Int}, DataType},1},
	oa_args::Array{<:Tuple{Union{Unknown,Int}, DataType},1} = oa_test;
	jacobian = nothing,
	kwargs...)

Generates a nonlinear form for the specified kernel function, test function operators, and argument operators evaluations. Operator evaluations are tuples that pair an unknown identifier or integer with a FunctionOperator. The header of the kernel functions needs to be conform to the interface

kernel!(result, input, qpinfo)

where qpinfo allows to access information at the current quadrature point.

During assembly the Newton update is computed via local jacobians of the kernel which are calculated by automatic differentiation or by the user-provided jacobian function with interface

jacobian!(jac, input_args, params)

Keyword arguments:

  • extra_inputsize: additional entries in input vector (e.g. for type-stable storage for intermediate resutls). Default: 0

  • factor: factor that should be multiplied during assembly. Default: 1

  • regions: subset of regions where operator should be assembly only. Default: Any[]

  • sparse_jacobians: use sparse jacobians. Default: true

  • name: name for operator used in printouts. Default: ''NonlinearOperator''

  • bonus_quadorder: additional quadrature order added to quadorder. Default: 0

  • entities: assemble operator on these grid entities (default = ONCELLS). Default: ONCELLS

  • quadorder: quadrature order. Default: ''auto''

  • entry_tolerance: threshold to add entry to sparse matrix. Default: 0

  • time_dependent: operator is time-dependent ?. Default: false

  • verbosity: verbosity level. Default: 0

  • params: array of parameters that should be made available in qpinfo argument of kernel function. Default: nothing

  • parallel_groups: assemble operator in parallel using CellAssemblyGroups. Default: false

  • sparsejacobianspattern: user provided sparsity pattern for the sparse jacobians (in case automatic detection fails). Default: nothing

ExtendableFEM.ProblemDescriptionType
struct ProblemDescription

Structure holding data for a problem description with the following fields:

  • name::String: The name of the problem used for printout messages. Default: "My Problem"
  • unknowns::Vector{Unknown}: A vector of Unknowns that are involved in the problem.
  • operators::Vector{AbstractOperator}: A vector of operators that are involved in the problem.
ExtendableFEM.SolverConfigurationMethod
function iterate_until_stationarity(
	SolverConfiguration(Problem::ProblemDescription
	[FES::Union{<:FESpace, Vector{<:FESpace}}];
	init = nothing,
	unknowns = Problem.unknowns,
	kwargs...)

Returns a solver configuration for the ProblemDescription that can be passed to the solve function. Here, FES are the FESpaces that should be used to discretize the selected unknowns. If no FES is provided an initial FEVector (see keyword init) must be provided (which is used to built the FES).

Keyword arguments:

  • verbosity: verbosity level. Default: 0

  • method_linear: any solver or custom LinearSolveFunction compatible with LinearSolve.jl (default = UMFPACKFactorization()). Default: LinearSolve.UMFPACKFactorization(true, true)

  • time: current time to be used in all time-dependent operators. Default: 0.0

  • spy: show unicode spy plot of system matrix during solve. Default: false

  • show_config: show configuration at the beginning of solve. Default: false

  • preconlinear: function that computes preconditioner for methodlinear incase an iterative solver is chosen. Default: nothing

  • init: initial solution (also used to save the new solution). Default: nothing

  • initialized: linear system in solver configuration is already assembled (turns true after first solve). Default: false

  • return_config: solver returns solver configuration (including A and b of last iteration). Default: false

  • show_matrix: show system matrix after assembly. Default: false

  • reltol: reltol for linear solver (if iterative). Default: 1.0e-11

  • is_linear: linear problem (avoid reassembly of nonlinear operators to check residual). Default: ''auto''

  • symmetrize: make system matrix symmetric (replace by (A+A')/2). Default: false

  • inactive: inactive unknowns (are made available in assembly, but not updated in solve). Default: ExtendableFEM.Unknown[]

  • plot: plot all solved unknowns with a (very rough but fast) unicode plot. Default: false

  • damping: amount of damping, value should be between in (0,1). Default: 0

  • check_matrix: check matrix for symmetry and positive definiteness and largest/smallest eigenvalues. Default: false

  • restrict_dofs: array of dofs for each unknown that should be solved (default: all dofs). Default: Any[]

  • constant_rhs: right-hand side is constant (skips reassembly). Default: false

  • abstol: abstol for linear solver (if iterative). Default: 1.0e-11

  • symmetrize_structure: make the system sparse matrix structurally symmetric (e.g. if [j,k] is also [k,j] must be set, all diagonal entries must be set). Default: false

  • target_residual: stop if the absolute (nonlinear) residual is smaller than this number. Default: 1.0e-10

  • maxiterations: maximal number of nonlinear iterations/linear solves. Default: 10

  • constant_matrix: matrix is constant (skips reassembly and refactorization in solver). Default: false

ExtendableFEM.UnknownType
struct Unknown

Structure holding information for an unknwon with the following fields:

  • name::String: The name of the unknown used for printout messages.
  • identifier::Any: The identifier of the unknown used for assignments to operators.
  • parameters::Dict{Symbol, Any}: Further properties of the unknown can be stored in a Dict, see constructor.
ExtendableFEM.UnknownMethod
function Unknown(
	u::String;
	identifier = Symbol(u),
	name = u,
	kwargs...)

Generates and returns an Unknown with the specified name, identifier and other traits.

Example: BilinearOperator([grad(1)], [grad(1)]; kwargs...) generates a weak Laplace operator.

Keyword arguments:

  • algebraic_constraint: is this unknown an algebraic constraint?. Default: nothing

  • symbol_ansatz: symbol for ansatz functions of this unknown in printouts. Default: nothing

  • symbol_test: symbol for test functions of this unknown in printouts. Default: nothing

  • dimension: dimension of the unknown. Default: nothing

CommonSolve.solveFunction
function solve(
	PD::ProblemDescription,
	[FES::Union{<:FESpace,Vector{<:FESpace}}],
	SC = nothing;
	unknowns = PD.unknowns,
	kwargs...)

Returns a solution of the PDE as an FEVector for the provided FESpace(s) FES (to be used to discretised the unknowns of the PDEs). If no FES is provided an initial FEVector (see keyword init) must be provided (which is used to built the FES).

This function extends the CommonSolve.solve interface and the PDEDescription takes the role of the ProblemType and FES takes the role of the SolverType.

Keyword arguments:

  • verbosity: verbosity level. Default: 0

  • method_linear: any solver or custom LinearSolveFunction compatible with LinearSolve.jl (default = UMFPACKFactorization()). Default: LinearSolve.UMFPACKFactorization(true, true)

  • time: current time to be used in all time-dependent operators. Default: 0.0

  • spy: show unicode spy plot of system matrix during solve. Default: false

  • show_config: show configuration at the beginning of solve. Default: false

  • preconlinear: function that computes preconditioner for methodlinear incase an iterative solver is chosen. Default: nothing

  • init: initial solution (also used to save the new solution). Default: nothing

  • initialized: linear system in solver configuration is already assembled (turns true after first solve). Default: false

  • return_config: solver returns solver configuration (including A and b of last iteration). Default: false

  • show_matrix: show system matrix after assembly. Default: false

  • reltol: reltol for linear solver (if iterative). Default: 1.0e-11

  • is_linear: linear problem (avoid reassembly of nonlinear operators to check residual). Default: ''auto''

  • symmetrize: make system matrix symmetric (replace by (A+A')/2). Default: false

  • inactive: inactive unknowns (are made available in assembly, but not updated in solve). Default: ExtendableFEM.Unknown[]

  • plot: plot all solved unknowns with a (very rough but fast) unicode plot. Default: false

  • damping: amount of damping, value should be between in (0,1). Default: 0

  • check_matrix: check matrix for symmetry and positive definiteness and largest/smallest eigenvalues. Default: false

  • restrict_dofs: array of dofs for each unknown that should be solved (default: all dofs). Default: Any[]

  • constant_rhs: right-hand side is constant (skips reassembly). Default: false

  • abstol: abstol for linear solver (if iterative). Default: 1.0e-11

  • symmetrize_structure: make the system sparse matrix structurally symmetric (e.g. if [j,k] is also [k,j] must be set, all diagonal entries must be set). Default: false

  • target_residual: stop if the absolute (nonlinear) residual is smaller than this number. Default: 1.0e-10

  • maxiterations: maximal number of nonlinear iterations/linear solves. Default: 10

  • constant_matrix: matrix is constant (skips reassembly and refactorization in solver). Default: false

Depending on the detected/configured nonlinearities the whole system is either solved directly in one step or via a fixed-point iteration.

ExtendableFEM.HomogeneousBoundaryDataMethod
function HomogeneousBoundaryData(u; entities = ON_BFACES, kwargs...)

When assembled, the unknown u of the Problem will be penalized to zero on the boundary faces and boundary regions (via kwargs).

Keyword arguments:

  • value: constant value of the data. Default: 0

  • regions: subset of regions where operator should be assembly only. Default: Any[]

  • verbosity: verbosity level. Default: 0

  • name: name for operator used in printouts. Default: ''HomogeneousData''

  • penalty: penalty for fixed degrees of freedom. Default: 1.0e30

  • mask: array of zeros/ones to set which components should be set by the operator (only works with componentwise dofs, add a 1 or 0 to mask additional dofs). Default: Any[]

ExtendableFEM.assign_operator!Method
assign_operator!(PD::ProblemDescription, o::AbstractOperator)

Assigns the AbstractOperator o to the ProblemDescription PD and returns its position in the operators array of the ProblemDescription.

ExtendableFEM.assign_unknown!Method
assign_unknown!(PD::ProblemDescription, u::Unknown)

Assigns the Unknown u to the ProblemDescription PD and returns its position in the unknowns array of the ProblemDescription.

ExtendableFEM.evaluateMethod
function evaluate(
	O::ItemIntegratorDG,
	sol;
	time = 0,
	kwargs...)

Evaluates the ItemIntegratorDG for the specified solution and returns an matrix of size resultdim x num_items.

ExtendableFEM.evaluateMethod
function evaluate(
	O::ItemIntegrator,
	sol;
	time = 0,
	kwargs...)

Evaluates the ItemIntegrator for the specified solution and returns an matrix of size resultdim x num_items.

ExtendableFEM.generate_ODEProblemFunction
function generate_ODEProblem(
	PD::ProblemDescription,
	FES,
	tspan;
	mass_matrix = nothing)
	kwargs...)

Reframes the ProblemDescription inside the SolverConfiguration into an ODEProblem, for DifferentialEquations.jl where tspan is the desired time interval.

If no mass matrix is provided the standard mass matrix for the respective finite element space(s) for all unknowns is assembled.

Additional keyword arguments:

  • verbosity: verbosity level. Default: 0

  • sametol: tolerance to identify two solution vectors to be identical (and to skip reassemblies called by DifferentialEquations.jl). Default: 1.0e-15

  • constant_rhs: right-hand side is constant (skips reassembly). Default: false

  • init: initial solution (otherwise starts with a zero vector). Default: nothing

  • initialized: linear system in solver configuration is already assembled (turns true after first assembly). Default: false

  • constant_matrix: matrix is constant (skips reassembly and refactorization in solver). Default: false

ExtendableFEM.get_periodic_coupling_infoMethod
function get_periodic_coupling_info(FES, xgrid, b1, b2, is_opposite::Function; factor_vectordofs = "auto")

computes the dofs that have to be coupled for periodic boundary conditions on the given xgrid for boundary regions b1, b2. The isopposite function evaluates if two provided face midpoints are on opposite sides to each other (the mesh xgrid should be appropriate). For vector-valued FETypes the user can provide factorvectordofs to incorporate a sign change if needed. This is automatically done for all Hdiv-conforming elements and (for the normal-weighted face bubbles of) the Bernardi-Raugel element H1BR.

ExtendableFEM.iterate_until_stationarityFunction
function iterate_until_stationarity(
	SCs::Array{<:SolverConfiguration, 1},
	FES = nothing;
	maxsteps = 1000,
	init = nothing,
	unknowns = [SC.PD.unknowns for SC in SCs],
	kwargs...)

Iterates consecutively over all SolverConfigurations (each contains the ProblemDescription of the corressponding subproblem) until the residuals of all subproblems are below their tolerances and returns the solution of the combined unknowns of all subproblems. The additional argument maxsteps limits the number of these iterations If an initial vector init is provided it should contain all unknowns of the subproblems.

Using the SolverConfiguration instead of the ProblemDescription in the first argument allows to use different kwargs for each subproblem. The SolverConfiguration for each subproblem can be generated by

SolverConfiguration(PD::ProblemDescription; init = sol, kwargs...)

with the usual keyword arguments.

ExtendableFEM.replace_operator!Method
replace_operator!(PD::ProblemDescription, j::Int, o::AbstractOperator)

Replaces the j-th operator of the ProblemDescription PD by the new operator o. Here, j is the position in operator array returned by the assign_operator! function. Nothing is returned (as the new operator gets position j).

ExtendableFEMBase.evaluate!Method
function evaluate(
	b::AbstractMatrix,
	O::ItemIntegrator,
	sol::FEVector;
	time = 0,
	kwargs...)

Evaluates the ItemIntegrator for the specified solution into the matrix b.

ExtendableFEMBase.evaluate!Method
function evaluate(
	b::AbstractMatrix,
	O::ItemIntegrator,
	sol::Array{FEVEctorBlock};
	time = 0,
	kwargs...)

Evaluates the ItemIntegrator for the specified solution into the matrix b.

ExtendableFEMBase.evaluate!Method
function evaluate(
	b::AbstractMatrix,
	O::ItemIntegratorDG,
	sol::FEVector;
	time = 0,
	kwargs...)

Evaluates the ItemIntegratorDG for the specified solution into the matrix b.

ExtendableFEMBase.evaluate!Method
function evaluate(
	b::AbstractMatrix,
	O::ItemIntegratorDG,
	sol::Array{FEVEctorBlock};
	time = 0,
	kwargs...)

Evaluates the ItemIntegratorDG for the specified solution into the matrix b.