Gridap.FESpaces

Gridap.FESpacesModule

Exported names

AffineFEOperator, Assembler, AssemblyStrategy, BasisStyle, CLagrangianFESpace, CellConformity, CellFE, ConstantFESpace, Constrained, ConstraintStyle, DefaultAssemblyStrategy, DirichletFESpace, DiscreteModelWithFEMap, EvaluationFunction, FEBasis, FEFunction, FEOperator, FESolver, FESpace, FESpaceWithConstantFixed, FESpaceWithLinearConstraints, FiniteElements, GenericAssemblyStrategy, GenericSparseMatrixAssembler, GridWithFEMap, HomogeneousTrialFESpace, HomogeneousTrialFESpace!, LinearFESolver, LocalOperator, NonlinearFESolver, PatchAssembler, PolytopalFESpace, SingleFieldFEFunction, SingleFieldFESpace, SparseMatrixAssembler, TestFESpace, TrialFESpace, TrialFESpace!, UnConstrained, UnconstrainedFESpace, ZeroMeanFESpace, add_mesh_displacement!, allocate_matrix, allocate_matrix_and_vector, allocate_vector, assemble_matrix, assemble_matrix!, assemble_matrix_add!, assemble_matrix_and_vector, assemble_matrix_and_vector!, assemble_matrix_and_vector_add!, assemble_vector, assemble_vector!, assemble_vector_add!, col_map, col_mask, collect_and_merge_cell_matrix, collect_and_merge_cell_matrix_and_vector, collect_and_merge_cell_vector, collect_cell_matrix, collect_cell_matrix_and_vector, collect_cell_vector, compute_cell_bases_changes, compute_cell_space, compute_conforming_cell_dofs, compute_dirichlet_values_for_tags, compute_dirichlet_values_for_tags!, gather_dirichlet_values, gather_dirichlet_values!, gather_free_and_dirichlet_values!, gather_free_values, gather_free_values!, get_algebraic_operator, get_cell_constraints, get_cell_dof_ids, get_cell_dof_values, get_cell_is_dirichlet, get_cell_isconstrained, get_cell_shapefuns, get_cols, get_dirichlet_dof_ids, get_dirichlet_dof_tag, get_dirichlet_dof_values, get_dirichlet_values, get_dof_value_type, get_fe_basis, get_fe_dof_basis, get_fe_space, get_free_dof_ids, get_free_dof_values, get_free_values, get_matrix_builder, get_matrix_type, get_rows, get_test, get_trial, get_trial_fe_basis, get_vector_builder, get_vector_type, has_constraints, interpolate, interpolate!, interpolate_dirichlet, interpolate_dirichlet!, interpolate_everywhere, interpolate_everywhere!, num_dirichlet_dofs, num_dirichlet_tags, num_free_dofs, numeric_loop_matrix!, numeric_loop_matrix_and_vector!, numeric_loop_vector!, row_map, row_mask, scatter_free_and_dirichlet_values, symbolic_loop_matrix!, symbolic_loop_matrix_and_vector!, symbolic_loop_vector!, test_assembler, test_fe_function, test_fe_operator, test_fe_solver, test_fe_space, test_single_field_fe_space, test_sparse_matrix_assembler, update_coordinates!, zero_dirichlet_values, zero_free_values,

source
Gridap.FESpaces.CellConformityMethod
CellConformity(cell_reffe::AbstractArray{<:ReferenceFE},conf::Conformity)
CellConformity(cell_polys::AbstractArray{<:ReferenceFE},cell_basis::AbstractArray{<:AbstractArray{<:Field}},conf::Conformity)
source
Gridap.FESpaces.CellFEType
struct CellFE <: GridapType

Minimum data required to build a conforming FE space. At this moment, the some cell-wise info is compressed on cell types. This can be relaxed in the future, and have an arbitrary cell-wise data.

source
Gridap.FESpaces.CellFEMethod
CellFE(model::DiscreteModel,cell_reffe::AbstractArray{<:ReferenceFE},conformity::Conformity, args...)

Generate a CellFE from a vector of reference fes

source
Gridap.FESpaces.ConstantFESpaceType
struct ConstantFESpace <: SingleFieldFESpace
  # private fields
end

ConstantFESpace(model::DiscreteModel; vector_type=Vector{Float64}, field_type=Float64)
ConstantFESpace(trian::Triangulation; vector_type=Vector{Float64}, field_type=Float64)

FESpace that is constant over the provided model/triangulation. Typically used as lagrange multipliers. The kwargs vector_type and field_type are used to specify the types of the dof-vector and dof-value respectively.

source
Gridap.FESpaces.DOFScalingMapMethod
DOFScalingMap(model, cell_reffe, pushforward, nothing)
DOFScalingMap(model, cell_reffe, pushforward, global_meshsize::Real)

Evaluated at a cell id of model to return a Diagonal matrix that rescales the cell physical shape-function basis, which is the result of mapping the cell_reffe[cell] basis using pushforward.

global_meshsize is either nothing or a Real number. If nothing, the local meshsize is estimated on each d dimensional face using the d-root of its d-volume.

DOFScalingMap is designed to be robust to heterogeneous cell reffes, and to reffes having heterogeneous DOF scaling (even within a face) like Mardal-Tai-Winter or C1 reffes.

The method to overload in order to implement nontrivial DOFScalingMap to a new reffe is get_dofscale_setter_function. See also "FE basis transformations".

source
Gridap.FESpaces.DiscontinuousCellConformityType
struct DiscontinuousCellConformity <: CellConformity

A CellConformity implementation for discontinuous FE spaces. It assumes that all the local DOFs belong to the cell interior, and allows for an arbitrary number of DOFs per cell (no compression).

source
Gridap.FESpaces.FEFunctionMethod
FEFunction(
  fs::SingleFieldFESpace, free_values::AbstractVector, dirichlet_values::AbstractVector)

The resulting FEFunction will be in the space if and only if dirichlet_values are the ones provided by get_dirichlet_dof_values(fs) Return a SingleFieldFEFunction.

source
Gridap.FESpaces.FEOperatorType
abstract type FEOperator <: GridapType

A FEOperator contains finite element problem, that is assembled as far as possible and ready to be solved.

source
Gridap.FESpaces.FESpaceMethod
FESpace(model::DiscreteModel, reffe; kwargs...)
FESpace(trian::Triangulation, reffe; kwargs...)

High level finite element space constructors. The reffe can be anything returned by a ReferenceFE constructor, but using one without the Polytope argument is safer, because the polytope(s) is(are) then chosen consistently with model or trian.

Keyword Arguments

  • conformity=nothing: prescribes the space Conformity, or use reffe's default if nothing.

  • dirichlet_tags=Int[]: tags of entities where Dirichlet boundary conditions are prescribed.

  • dirichlet_masks=nothing: for Cartesian product finite elements of value type V::MultiValue, specifies to which components of V do the Dirichlet BC apply, for each tag in dirichlet_tags. The type of dirichlet_masks is Nothing or Vector{NTuple{N,Bool}} where N = num_indep_components(V).

  • vector_type=nothing: sets the type of the vector storing the DOF values, e.g. Vector{Float64}. Useful to create complex valued function space from a real finite element basis.

  • constraint=nothing: if set to :zeromean, adds a constraint to the space to ensure its functions have zero mean-value.

  • scale_dof=false: if true, rescales the DOF to make them independent from the mesh size. Usefull for mixed-elements using different Piola/physical mappings, that would introduce heterogeneous scaling with meshsize in the systems. By default, the local mesh size is estimated on each physical d-faces/cell using its d-volume.

  • global_meshsize=nothing: if scale_dof, can be used to prescribe a constant homogeneous mesh size estimate, to avoid computing the volume of each physical face owning DOF. This should only be used on very shape-regular and quasi-uniform meshes (e.g. unstretched CartesianDiscreteModel).

  • labels=get_face_labeling(model): FaceLabeling on top of model that dirichlet_tags refer to.

  • trian=Triangulation(model).

Extended help

Internal constructors:

FESpace(model::DiscreteModel, cell_reffe::AbstractArray{<:ReferenceFE}; kwargs...)

FESpace(model::DiscreteModel, cell_fe::CellFE;
  trian, labels, dirichlet_tags, dirichlet_masks, constraint, vector_type)
source
Gridap.FESpaces.GridWithFEMapType

Given a discrete model and a reffe, builds a new grid in which the geometrical map is a FEFunction. This is useful when considering geometrical maps that are the result of a FE problem (mesh displacement).

source
Gridap.FESpaces.MultiConstantFESpaceType
struct MultiConstantFESpace{N,V,T} <: SingleFieldFESpace end

MultiConstantFESpace(model::DiscreteModel,tags::Vector,D::Integer)
MultiConstantFESpace(trians::Vector{<:BoundaryTriangulation{Df}})

Extension of ConstantFESpace, representing a FESpace which is constant on each of it's N triangulations.

source
Gridap.FESpaces.NodeToDofGlueType
struct NodeToDofGlue{T}
  free_dof_to_node::Vector{Int32}
  free_dof_to_comp::Vector{Int16}
  dirichlet_dof_to_node::Vector{Int32}
  dirichlet_dof_to_comp::Vector{Int16}
  node_and_comp_to_dof::Vector{T}
end
source
Gridap.FESpaces.NormalSignMapType
struct NormalSignMap <: Map
  ...
end

The NormalSignMap compute the signs to apply to the mapped reference normals, for each facet of a physical cell.

Each physical facet $f$ is shared by up to two cells $K$ and $K'$. The orientation of the physical/global normal to $f$ is chosen by the main cell $K$, the first one in the list of adjascent cells to $f$ in the grid topology.

The physical/global normal is the (normalized) Piola mapped reference normal to $f̂ = F⁻¹(f)$ where $F$ is the geometrical map $F:K̂->K$. It is also minus the (normalized) Piola mapped reference normal to $f̂' = F'⁻¹(f)$ where $F'$ is the geometrical map $F':K̂->K'$.

source
Gridap.Algebra.solve!Method
uh, cache = solve!(uh,solver,op,cache)

This function changes the state of the input and can render it in a corrupted state. It is recommended to rewrite the input uh with the output as illustrated to prevent any issue. If cache===nothing, then it creates a new cache object.

source
Gridap.Algebra.solve!Method
uh, cache = solve!(uh,solver,op)

This function changes the state of the input and can render it in a corrupted state. It is recommended to rewrite the input uh with the output as illustrated to prevent any issue.

source
Gridap.Algebra.solveMethod
solve(nls::FESolver, op::FEOperator)

Solve that allocates, and sets initial guess to zero and returns the solution.

source
Gridap.FESpaces.HomogeneousTrialFESpace!Method
HomogeneousTrialFESpace!(dirichlet_values::AbstractVector,U::SingleFieldFESpace)

Return a SingleFieldFESpace that is U but with Dirichlet values set to zero. Uses dirichlet_values with zeros set in place for the container of the Dirichlet values of the returned space.

source
Gridap.FESpaces.compute_cell_bases_changesMethod
compute_cell_bases_changes(name::ReferenceFEName, push::Pushforward, model, cell_reffe, cell_Jt)

Computes, in each cell, the change of basis $M$ between the pushforwarded reference shape-function basis and the expected physical shape-functions, as well as it transposed inverse $M⁻ᵀ$, the change of basis between the inverse pullback of the reference cell DOF and the expected physical cell DOF.

For that, the model is provided. As well as the reference FE + lazy gradient (transposed Jacobian) of the geometrical map, in each cell.

See also the manual page on "FE basis transformations". Return either nothing (no change required) or a couple of cell arrays (cell_M, cell_M⁻ᵀ)

The dof_scale and global_meshsize kwargs are not handeled by this function. For them to work properly, it might be necessary to also overload get_dofscale_setter_function.

source
Gridap.FESpaces.compute_conforming_cell_dofsFunction
compute_conforming_cell_dofs(
  cell_conformity::CompressedCellConformity,
  grid_topology, face_labeling, dirichlet_tags, dirichlet_components=nothing
)

The result is the tuple

(cell_dofs, nfree, ndiri, dirichlet_dof_tag, dirichlet_cells)

If dirichlet_components is given, then get_dof_to_comp has to be defined for the reference elements in reffes.

source
Gridap.FESpaces.generate_cell_dof_maskMethod
generate_cell_dof_mask(
  scell_conformity::CellConformity,
  tcell_to_scell::AbstractVector{<:Integer},
  d_to_tcell_to_tdface,
  d_to_tdface_to_mask;
  reverse::Bool=false
)

Given a CellConformity object, defined on a set of source cells (scell), and given a cell/face/edge/node masks on a set of target cells (tcell), this function generates a mask for the degrees of freedom (dofs) on the target cells.

Parameters:

  • scell_conformity: The CellConformity object on the source cells (scell).
  • tcell_to_scell: A vector mapping target cells (tcell) to source cells (scell).
  • d_to_tcell_to_tdface: For each dimension d, a Table mapping target cells to its d-faces (edges, faces, ...)
  • d_to_tdface_to_mask: For each dimension d, an array of booleans indicating whether the d-face is masked or not.

Returns:

  • tcell_dof_mask: A vector that for each target cell contains a boolean vector of size equal to the number of dofs in that cell, containing the mask.

Modes of operation:

  • If reverse = false, the function generates a mask where the dofs are true if the corresponding d-face is masked (true).
  • If reverse = true, the mask is reversed, meaning that the dofs are true if the corresponding d-face is not masked (false).
source
Gridap.FESpaces.generate_dof_maskMethod
generate_dof_mask(
  space::FESpace, labels::FaceLabeling, tags;
  cell_conformity = get_cell_conformity(space), reverse::Bool=false
)

Generate a mask for the dofs in the FESpace space, based on the face labeling labels and the tags provided. If reverse is false, the dofs are masked (true) if they belong to a tagged face. If reverse is true, the mask is reversed.

source
Gridap.FESpaces.get_algebraic_operatorMethod
get_algebraic_operator(feop)

Return an "algebraic view" of an operator. Algebraic means, that the resulting operator acts on plain arrays, instead of FEFunctions. This can be useful for solving with external tools like NLsolve.jl. See also FEOperator.

source
Gridap.FESpaces.get_cell_shapefuns_and_dof_basisMethod
get_cell_shapefuns_and_dof_basis(
    model::DiscreteModel, cell_reffe::AbstractArray, conf::Conformity; kwargs...)

Return (shapefuns, dof_basis), the physical cell arrays of shape-function and DOF.

These arrays are computed by mapping the reference shape-function and DOF bases to the physical elements cell-wise, taking into account possible non-trivial Piola maps, DOF change of basis, sign change and rescaling to build conf-conforming physical FEs.

cell_reffe is an array of reference FE per cell of model, all reffes must share the same ReferenceFEName and pushforward p when used with the given conformity.

The kwargs are scale_dof and global_meshsize, c.f. FESpace.

For performances, cell_reffe should efficiently be compressed by [compress_cell_data].

source
Gridap.FESpaces.interpolateMethod

The resulting FE function is in the space (in particular it fulfills Dirichlet BCs even in the case that the given cell field does not fulfill them)

source