Gridap.Adaptivity
The adaptivity module provides a framework to work with adapted (refined/coarsened/mixed) meshes.
It provides
- A generic interface to represent adapted meshes and a set of tools to work with Finite Element spaces defined on them. In particular, moving
CellFields
between parent and child meshes. - Particular implementations for conformally refining/coarsening 2D/3D meshes using several well-known strategies. In particular, Red-Green refinement and longest-edge bisection.
Interface
The following types are defined in the module:
Gridap.Adaptivity.RefinementRule
— TypeStructure representing the map between a single parent cell and its children.
Contains:
- T ::
RefinementRuleType
, indicating the refinement method. - poly ::
Polytope
, representing the geometry of the parent cell. - ref_grid ::
DiscreteModel
defined onpoly
, giving the parent-to-children cell map.
Gridap.Adaptivity.AdaptivityGlue
— TypeGlue containing the map between two nested triangulations. The contained datastructures will depend on the type of glue. There are two types of AdaptivityGlue
:
RefinementGlue
:: All cells in the new mesh are children of cells in the old mesh. I.e given a new cell, it is possible to find a single old cell containing it (the new cell might be exactly the old cell if no refinement).MixedGlue
:: Some cells in the new mesh are children of cells in the old mesh, while others are parents of cells in the old mesh.
Contains:
n2o_faces_map
:: Given a new face gid, returns- if fine, the gid of the old face containing it.
- if coarse, the gids of its children (in child order)
n2o_cell_to_child_id
:: Given a new cell gid, returns- if fine, the local child id within the (old) coarse cell containing it.
- if coarse, a list of local child ids of the (old) cells containing it.
refinement_rules
:: Array conatining theRefinementRule
used for each coarse cell.
Gridap.Adaptivity.AdaptedDiscreteModel
— TypeDiscreteModel
created by refining/coarsening another DiscreteModel
.
The refinement/coarsening hierarchy can be traced backwards by following the parent
pointer chain. This allows the transfer of dofs between FESpaces
defined on this model and its ancestors.
Gridap.Adaptivity.AdaptedTriangulation
— TypeTriangulation produced from an AdaptedDiscreteModel.
Contains:
- adapted_model ::
AdaptedDiscreteModel
for the triangulation. - trian ::
Triangulation
extracted from the background model, i.eget_model(adapted_model)
.
The high-level interface is provided by the following methods:
Gridap.Adaptivity.refine
— Functionfunction refine(model::DiscreteModel,args...;kwargs...) :: AdaptedDiscreteModel
Returns an AdaptedDiscreteModel
that is the result of refining the given DiscreteModel
.
Gridap.Adaptivity.coarsen
— Functionfunction coarsen(model::DiscreteModel,args...;kwargs...) :: AdaptedDiscreteModel
Returns an AdaptedDiscreteModel
that is the result of coarsening the given DiscreteModel
.
Gridap.Adaptivity.adapt
— Functionfunction adapt(model::DiscreteModel,args...;kwargs...) :: AdaptedDiscreteModel
Returns an AdaptedDiscreteModel
that is the result of adapting (mixed coarsening and refining) the given DiscreteModel
.
Edge-Based refinement
The module provides a refine
method for UnstructuredDiscreteModel
. The method takes a string refinement_method
that determines the refinement strategy to be used. The following strategies are available:
"red_green"
:: Red-Green refinement, default."nvb"
:: Longest-edge bisection (only for meshes of TRIangles)"barycentric"
:: Barycentric refinement (only for meshes of TRIangles)"simplexify"
:: Simplexify refinement. Same resulting mesh as thesimplexify
method, but keeps track of the parent-child relationships.
Additionally, the method takes a kwarg cells_to_refine
that determines which cells will be refined. Possible input types are:
Nothing
:: All cells get refined.AbstractArray{<:Bool}
of sizenum_cells(model)
:: Only cells such thatcells_to_refine[iC] == true
get refined.AbstractArray{<:Integer}
:: Cells for whichgid ∈ cells_to_refine
get refined
The algorithms try to respect the cells_to_refine
input as much as possible, but some additional cells might get refined in order to guarantee that the mesh remains conforming.
function refine(model::UnstructuredDiscreteModel;refinement_method="red_green",kwargs...)
[...]
end
CartesianDiscreteModel refining
The module provides a refine
method for CartesianDiscreteModel
. The method takes a Tuple
of size Dc
(the dimension of the model cells) that will determine how many times cells will be refined in each direction. For example, for a 2D model, refine(model,(2,3))
will refine each QUAD cell into a 2x3 grid of cells.
function refine(model::CartesianDiscreteModel{Dc}, cell_partition::Tuple) where Dc
[...]
end
Macro Finite-Elements
The module also provides support for macro finite-elements. From an abstract point of view, a macro finite-element is a finite-element defined on a refined polytope, where polynomial basis are defined on each of the subcells (creating a broken piece-wise polynomial space on the original polytope). From Gridap's point of view, a macro finite-element is a ReferenceFE
defined on a RefinementRule
from an array of ReferenceFE
s defined on the subcells.
Although there are countless combinations, here are two possible applications:
- Linearized High-Order Lagrangian FESpaces: These are spaces which have the same DoFs as a high-order Lagrangian space, but where the basis functions are linear on each subcell.
- Barycentrically-refined elements for Stokes-like problems: These are spaces where the basis functions for velocity are defined on the barycentrically-refined mesh, whereas the basis functions for pressure are defined on the original cells. This allows for exact so-called Stokes sequences (see here).
The API is given by the following methods:
Gridap.Adaptivity.MacroReferenceFE
— FunctionMacroReferenceFE(rrule::RefinementRule,reffes::AbstractVector{<:ReferenceFE})
Constructs a ReferenceFE for a macro-element, given a RefinementRule and a set of ReferenceFEs for the subcells.
For performance, these should be paired with CompositeQuadratures.
Adaptive Mesh Refinement
One of the main uses of mesh refinement is Adaptive Mesh Refinement, where the mesh is refined only in regions of interest.
The typical AMR workflow is the so-called solve-estimate-mark-refine
loop. Since estimators will generally be problem-dependent, we only aim to provide some generic tools that can be combined by the user:
Gridap.Adaptivity.DorflerMarking
— Typestruct DorflerMarking
θ :: Float64
ν :: Float64
strategy :: Symbol
end
DorflerMarking(θ::Float64; ν::Float64 = 0.5, strategy::Symbol = :quickmark)
Implements the Dorfler marking strategy. Given a vector η
of real positive numbers, the marking strategy find a subset of indices I
such that
sum(η[I]) > θ * sum(η)
where 0 < θ < 1
is a threshold parameter.
For more details, see the following reference:
"Dörfler marking with minimal cardinality is a linear complexity problem", Pfeiler et al. (2020)
The marking algorithm is controlled by the strategy
parameter, which can take the following values:
:sort
: Optimal cardinality, O(N log N) complexity. See Algorithm 2 in the reference.:binsort
: Quasi-optimal cardinality, O(N) complexity. See Algorithm 7 in the reference.:quickmark
: Optimal cardinality, O(N) complexity. See Algorithm 10 in the reference.
Arguments
θ::Float64
: The threshold parameter. Between 0 and 1.ν::Float64
: Extra parameter for:binsort
. Default is 0.5.strategy::Symbol
: The marking strategy. Default is:quickmark
.
Usage
η = abs.(randn(1000))
m = DorflerMarking(0.5)
I = mark(m,η)
Gridap.Adaptivity.mark
— Functionmark(m::DorflerMarking, η::Vector{<:Real}) -> Vector{Int}
Given a vector η
of real positive numbers, returns a subset of indices I
such that satisfying the Dorfler marking condition.
Gridap.Adaptivity.estimate
— Functionestimate(f::Function, uh::Function) -> Vector{Float64}
Given a functional f
and a function uh
, such that f(uh)
produces a scalar-valued DomainContribution
, collects the estimator values for each cell in the background model.
Notes for users
Most of the tools provided by this module are showcased in the tests of the module itself, as well as the following tutorial (coming soon).
However, we want to stress a couple of key performance-critical points:
The refining/coarsening routines are not optimized for performance. In particular, they are not parallelized. If you require an optimized/parallel implementation, please consider leveraging specialised meshing libraries. For instance, we provide an implementation of
refine/coarsen
using P4est in the GridapP4est.jl library.Although the toolbox allows you to evaluate
CellFields
defined on both fine/coarse meshes on their parent/children mesh, both directions of evaluation are not equivalent. As a user, you should always try to evaluate/integrate on the finest mesh for maximal performance. Evaluating a fineCellField
on a coarse mesh relies on local tree searches, and is therefore a very expensive operation that should be avoided whenever possible.
Notes for developers
RefinementRule API
Given a RefinementRule
, the library provides a set of methods to compute the mappings between parent (coarse) face ids and child (fine) face ids (and vice-versa).
The most basic information (that can directly be hardcoded in the RefinementRule for performance) are the mappings between parent face ids and child face ids. These are provided by:
Gridap.Adaptivity.get_d_to_face_to_child_faces
— Functionget_d_to_face_to_child_faces(rr::RefinementRule)
Given a RefinementRule
, returns for each parent/coarse face the child/fine faces of the same dimension that it contains. Therefore, only fine faces at the coarse cell boundary are listed in the returned structure.
Returns: [Face dimension][Coarse Face id] -> [Fine faces]
Gridap.Adaptivity.get_d_to_face_to_parent_face
— Functionget_d_to_face_to_parent_face(rr::RefinementRule)
Given a RefinementRule
, returns for each fine/child face the parent/coarse face containing it. The parent face can have higher dimension.
Returns the tuple (A,B) with
- A = [Face dimension][Fine Face id] -> [Parent Face]
- B = [Face dimension][Fine Face id] -> [Parent Face Dimension]
On top of these two basic mappings, a whole plethora of additional topological mappings can be computed. These first set of routines extend the ReferenceFEs API to provide information on the face-to-node mappings and permutations:
Gridap.ReferenceFEs.get_face_vertices
— FunctionReferenceFEs.get_face_vertices(rr::RefinementRule)
Given a RefinementRule
, returns for each parent/coarse face the ids of the child/fine vertices it contains.
get_face_vertices(p::Polytope) -> Vector{Vector{Int}}
get_face_vertices(p::Polytope,dim::Integer) -> Vector{Vector{Int}}
get_face_vertices(g::GridTopology,d::Integer)
get_face_vertices(g::GridTopology)
Defaults to
compute_face_vertices(g)
Gridap.ReferenceFEs.get_face_coordinates
— FunctionReferenceFEs.get_face_coordinates(rr::RefinementRule)
Given a RefinementRule
, returns for each parent/coarse face the coordinates of the child/fine vertices it contains.
get_face_coordinates(p::Polytope)
get_face_coordinates(p::Polytope,d::Integer)
get_face_coordinates(g::GridTopology)
get_face_coordinates(g::GridTopology,d::Integer)
Gridap.ReferenceFEs.get_vertex_permutations
— FunctionReferenceFEs.get_vertex_permutations(rr::RefinementRule)
Given a RefinementRule
, returns all possible permutations of the child/fine vertices within the cell.
get_vertex_permutations(p::Polytope) -> Vector{Vector{Int}}
Given a polytope p
, returns a vector of vectors containing all admissible permutations of the polytope vertices. An admissible permutation is one such that, if the vertices of the polytope are re-labeled according to this permutation, the resulting polytope preserves the shape of the original one.
Examples
using Gridap.ReferenceFEs
perms = get_vertex_permutations(SEGMENT)
println(perms)
# output
Array{Int,1}[[1, 2], [2, 1]]
The first admissible permutation for a segment is [1,2]
,i.e., the identity. The second one is [2,1]
, i.e., the first vertex is relabeled as 2
and the second vertex is relabeled as 1
.
Gridap.ReferenceFEs.get_face_vertex_permutations
— FunctionReferenceFEs.get_face_vertex_permutations(rr::RefinementRule)
Given a RefinementRule
, returns for each parent/coarse face the possible permutations of the child/fine vertices it contains.
get_face_vertex_permutations(p::Polytope)
get_face_vertex_permutations(p::Polytope,d::Integer)
We also provide face-to-face maps:
Gridap.Adaptivity.get_cface_to_num_own_ffaces
— Functionget_cface_to_num_own_ffaces(rr::RefinementRule)
Given a RefinementRule
, returns for each parent/coarse face the number of child/fine faces of all dimensions that it owns.
Gridap.Adaptivity.get_cface_to_own_ffaces
— Functionget_cface_to_own_ffaces(rr::RefinementRule)
Given a RefinementRule
, returns for each parent/coarse face the child/fine faces of all dimensions that it owns.
Gridap.Adaptivity.get_cface_to_ffaces
— Functionget_cface_to_ffaces(rr::RefinementRule)
Given a RefinementRule
, returns for each parent/coarse face the child/fine faces of all dimensions that are on it (owned and not owned).
The implementation aggregates the results of get_cface_to_own_ffaces
.
Gridap.Adaptivity.get_cface_to_own_ffaces_to_lnodes
— Functionget_cface_to_own_ffaces_to_lnodes(rr::RefinementRule)
Given a RefinementRule
, returns
[coarse face][local owned child face] -> local fine node ids
where local refers to the local fine numbering within the coarse face.
Gridap.Adaptivity.get_cface_to_ffaces_to_lnodes
— Functionget_cface_to_ffaces_to_lnodes(rr::RefinementRule)
Given a RefinementRule
, returns
[coarse face][local child face] -> local fine node ids
where local refers to the local fine numbering within the coarse face.
Gridap.Adaptivity.get_cface_to_fface_permutations
— Functionget_cface_to_fface_permutations(rrule::RefinementRule)
get_cface_to_own_fface_permutations(rrule::RefinementRule)
Given a RefinementRule
, this function returns:
cface_to_cpindex_to_ffaces
: For each coarse face, for each coarse face permutation, the permuted ids of the fine faces.cface_to_cpindex_to_fpindex
: For each coarse face, for each coarse face permutation, the sub-permutation of the fine faces.
The idea is the following: A permutation on a coarse face induces a 2-level permutation for the fine faces, i.e
- First, the fine faces get shuffled amongs themselves.
- Second, each fine face has it's orientation changed (given by a sub-permutation).
For instance, let's consider a 1D example, where a SEGMENT is refined into 2 segments:
3 4 5
X-----X --> X-----X-----X
1 2 1 3 2
Then when aplying the coarse node permutation (1,2) -> (2,1), we get the following fine face permutation:
- Faces (1,2,3,4,5) get mapped to (2,1,3,5,4)
- Moreover, the orientation of faces 3 and 5 is changed, i.e we get the sub-permutation (1,1,1,2,2)
Gridap.Adaptivity.aggregate_cface_to_own_fface_data
— Functionaggregate_cface_to_own_fface_data(
rr::RefinementRule,
cface_to_own_fface_to_data :: AbstractVector{<:AbstractVector{T}}
) where T
Given a RefinementRule
, and a data structure cface_to_own_fface_to_data
that contains data for each child/fine face owned by each parent/coarse face, returns a data structure cface_to_fface_to_data
that contains the data for each child/fine face contained in the closure of each parent/coarse face (i.e the fine faces are owned and not owned).
The implementation makes sure that the resulting data structure is ordered according to the fine face numbering in get_cface_to_ffaces(rrule)
(which in turn is by increasing fine face id).
Gridap.Adaptivity.get_face_subface_ldof_to_cell_ldof
— FunctionGiven a RefinementRule
of dimension Dc and a Dc-Tuple fine_orders
of approximation orders, returns a map between the fine nodal dofs of order fine_orders
in the reference grid and the coarse nodal dofs of order 2⋅fine_orders
in the coarse parent cell.
The result is given for each coarse/parent face of dimension D
as a list of the corresponding fine dof lids, i.e
- [coarse face][coarse dof lid] -> fine dof lid
AdaptivityGlue API
Gridap.Adaptivity.get_n2o_reference_coordinate_map
— FunctionFor each fine cell, returns the map Φ st. xcoarse = ϕ(xfine)
Gridap.Adaptivity.get_old_cell_refinement_rules
— FunctionGiven a RefinementGlue
, returns an array containing the refinement rules for each cell in the old mesh.
Gridap.Adaptivity.get_new_cell_refinement_rules
— FunctionGiven a RefinementGlue
, returns an array containing the refinement rules for each cell in the new mesh.
Gridap.Adaptivity.get_d_to_fface_to_cface
— FunctionFor each child/fine face, returns the parent/coarse face containing it. The parent face might have higher dimension.
Returns two arrays:
- [dimension][fine face gid] -> coarse parent face gid
- [dimension][fine face gid] -> coarse parent face dimension
Gridap.Adaptivity.n2o_reindex
— Functionfunction n2o_reindex(new_data,g::AdaptivityGlue) -> old_data
Reindexes a cell-wise array from the new mesh to the old mesh.
Gridap.Adaptivity.o2n_reindex
— Functionfunction o2n_reindex(old_data,g::AdaptivityGlue) -> new_data
Reindexes a cell-wise array from the old mesh to the new mesh.
New-to-old field evaluations
When a cell is refined, we need to be able to evaluate the fields defined on the children cells on the parent cell. To do so, we bundle the fields defined on the children cells into a new type of Field
called FineToCoarseField
. When evaluated on a Point
, a FineToCoarseField
will select the child cell that contains the Point
and evaluate the mapped point on the corresponding child field.
Gridap.Adaptivity.FineToCoarseField
— Typestruct FineToCoarseField <: Field
fine_fields :: AbstractVector{<:Field}
rrule :: RefinementRule
id_map :: AbstractVector{<:Integer}
end
Given a domain and a non-overlapping refined cover, a FineToCoarseField
is a Field
defined in the domain and constructed by a set of fields defined on the subparts of the covering partition. The refined cover is represented by a RefinementRule
.
Parameters:
rrule
: Refinement rule representing the covering partition.fine_fields
: Fields defined on the subcells of the covering partition. To accomodate the case where not all subcells are present (e.g in distributed), we allow forlength(fine_fields) != num_subcells(rrule)
.id_map
: Mapping from the subcell ids to the indices of thefine_fields
array. Iflength(fine_fields) == num_subcells(rrule)
, this is the identity map. Otherwise, theid_map
is used to map the subcell ids to the indices of thefine_fields
array.
Gridap.Adaptivity.FineToCoarseDofBasis
— TypeGridap.Adaptivity.FineToCoarseRefFE
— TypeWrapper for a ReferenceFE which is specialised for efficiently evaluating FineToCoarseFields.