Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Update some missing Lagrange docs #889

Merged
merged 15 commits into from
Jul 12, 2024
Merged
Show file tree
Hide file tree
Changes from 6 commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
30 changes: 23 additions & 7 deletions docs/src/devdocs/elements.md
Original file line number Diff line number Diff line change
Expand Up @@ -2,25 +2,41 @@

## Type definitions

Elements or cells are subtypes of `AbstractCell{dim,N,M}`. They are parametrized by
the dimension of their nodes via `dim`, the number of nodes `N` and the number
of faces `M`.
Elements or cells are subtypes of `AbstractCell{<:AbstractRefShape}`. As shown, they are parametrized
by the associated reference element.

### Required methods to implement for all subtypes of `AbstractCell` to define a new element

```@docs
Ferrite.get_node_ids
```

Please note that if your cell holds a single tuple called `nodes` to carry the node numbers, then
this method will work automatically. Almost all functions work automatically with the information
provided by the reference element.
termi-official marked this conversation as resolved.
Show resolved Hide resolved

Please note that if you want to implement a custom element which does not have a reference shape,
then the relevant functions below need to be dispatched.
termi-official marked this conversation as resolved.
Show resolved Hide resolved

### Common utilities and definitions when working with grids internally.

First we have some topological queries on the element

```@docs
Ferrite.vertices(::Ferrite.AbstractCell)
Ferrite.edges(::Ferrite.AbstractCell)
Ferrite.reference_faces(::Ferrite.AbstractRefShape)
Ferrite.faces(::Ferrite.AbstractCell)
Ferrite.geometric_interpolation(::Ferrite.AbstractCell)
Ferrite.facets(::Ferrite.AbstractCell)
Ferrite.boundaryfunction(::Type{<:Ferrite.BoundaryIndex})
Ferrite.reference_vertices(::Ferrite.AbstractCell)
Ferrite.reference_edges(::Ferrite.AbstractCell)
Ferrite.reference_faces(::Ferrite.AbstractCell)
```

### Common utilities and definitions when working with grids internally.
and some generic utils which are commonly found in finite element codes

```@docs
Ferrite.BoundaryIndex
Ferrite.boundaryfunction(::Type{<:Ferrite.BoundaryIndex})
Ferrite.get_coordinate_eltype(::Ferrite.AbstractGrid)
Ferrite.get_coordinate_eltype(::Node)
Ferrite.toglobal
Expand Down
16 changes: 16 additions & 0 deletions docs/src/devdocs/reference_cells.md
Original file line number Diff line number Diff line change
Expand Up @@ -15,3 +15,19 @@ Ferrite.RefTetrahedron
Ferrite.RefHexahedron
Ferrite.RefPrism
```

### Required methods to implement for all subtypes of `AbstractRefShape` to define a new reference element
termi-official marked this conversation as resolved.
Show resolved Hide resolved

```@docs
Ferrite.geometric_interpolation(::Ferrite.AbstractRefShape)
termi-official marked this conversation as resolved.
Show resolved Hide resolved
Ferrite.reference_vertices(::Ferrite.AbstractRefShape)
Ferrite.reference_edges(::Ferrite.AbstractRefShape)
Ferrite.reference_faces(::Ferrite.AbstractRefShape)
termi-official marked this conversation as resolved.
Show resolved Hide resolved
```

which automatically defines


```@docs
Ferrite.reference_facets(::Ferrite.AbstractRefShape)
termi-official marked this conversation as resolved.
Show resolved Hide resolved
```
10 changes: 5 additions & 5 deletions docs/src/topics/degrees_of_freedom.md
Original file line number Diff line number Diff line change
Expand Up @@ -25,8 +25,8 @@ need to specify number of components for the field. Here we add a vector field `
(2 components for a 2D problem) and a scalar field `:p`.

```@example dofs
add!(dh, :u, Lagrange{2,RefTetrahedron,1}()^2)
add!(dh, :p, Lagrange{2,RefTetrahedron,1}())
add!(dh, :u, Lagrange{RefTriangle, 1}()^2)
add!(dh, :p, Lagrange{RefTriangle, 1}())
# hide
```

Expand All @@ -53,12 +53,12 @@ quadratic interpolation, and our `:p` field with a linear approximation.

```@example dofs
dh = DofHandler(grid) # hide
add!(dh, :u, Lagrange{2,RefTetrahedron,2}()^2)
add!(dh, :p, Lagrange{2,RefTetrahedron,1}())
add!(dh, :u, Lagrange{RefTriangle, 2}()^2)
add!(dh, :p, Lagrange{RefTriangle, 1}())
# hide
```

## Ordering of Dofs

ordered in the same order as we add to dofhandler
nodes -> (edges ->) faces -> cells
vertices -> edges -> faces -> volumes
91 changes: 59 additions & 32 deletions docs/src/topics/grid.md
Original file line number Diff line number Diff line change
Expand Up @@ -18,7 +18,7 @@ Either, a mesh is created on the fly with the gmsh API or a mesh in `.msh` or `.
```@docs
FerriteGmsh.togrid
```
`FerriteGmsh.jl` supports currently the translation of `cellsets` and `facesets`.
`FerriteGmsh.jl` supports currently the translation of `cellsets` and `facetsets`.
Such sets are defined in Gmsh as `PhysicalGroups` of dimension `dim` and `dim-1`, respectively.
In case only a part of the mesh is the domain, the domain can be specified by providing the keyword argument `domain` the name of the `PhysicalGroups` in the [`FerriteGmsh.togrid`](@ref) function.

Expand Down Expand Up @@ -71,67 +71,94 @@ julia> cells = [
```

where each Quadrilateral, which is a subtype of `AbstractCell` saves in the field `nodes` the tuple of node IDs.
Additionally, the data structure `Grid` can hold node-, face- and cellsets.
Additionally, the data structure `Grid` can hold node-, facet- and cellsets.
All of these three sets are defined by a dictionary that maps a string key to a `Set`.
For the special case of node- and cellsets the dictionary's value is of type `Set{Int}`, i.e. a keyword is mapped to a node or cell ID, respectively.

Facesets are a more elaborate construction. They map a `String` key to a `Set{FaceIndex}`, where each `FaceIndex` consists of `(global_cell_id, local_face_id)`.
In order to understand the `local_face_id` properly, one has to consider the reference space of the element, which typically is spanned by a product of the interval ``[-1, 1]`` and in this particular example ``[-1, 1] \times [-1, 1]``.
In this space a local numbering of nodes and faces exists, i.e.
Facesets are a more elaborate construction. They map a `String` key to a `Set{FaceIndex}`, where each `FaceIndex` consists of `(global_cell_id, local_facet_id)`.
In order to understand the `local_facet_id` properly, one has to consider the reference space of the element, which typically is spanned by a product of the interval ``[-1, 1]`` and in this particular example ``[-1, 1] \times [-1, 1]``.
In this space a local numbering of nodes and facets exists, i.e.


![local element](./assets/local_element.svg)


The example shows a local face ID ordering, defined as:
The example shows a local edge ordering, defined as:

```julia
faces(::Lagrange{2,RefCube,1}) = ((1,2), (2,3), (3,4), (4,1))
Ferrite.reference_facets(RefQuadrilateral) = ((1,2), (2,3), (3,4), (4,1))
```

Other face ID definitions [can be found in the src files](https://github.com/Ferrite-FEM/Ferrite.jl/blob/8224282ab4d67cb523ef342e4a6ceb1716764ada/src/interpolations.jl#L154) in the corresponding `faces` dispatch.
Other facet definitions can be found in the src file `src/Grid/grid.jl` in the corresponding dispatches for [`reference_facets`](@ref Ferrite.reference_facets). Furthermorem you can query specific information about subentities via [`reference_vertices`](@ref Ferrite.reference_vertices), [`reference_edges`](@ref Ferrite.reference_edges) and [`reference_faces`](@ref Ferrite.reference_faces).


The highlighted face, i.e. the two lines from node ID 3 to 6 and from 6 to 9, on the right hand side of our test mesh can now be described as
The highlighted facet, i.e. the two lines from node ID 3 to 6 and from 6 to 9, on the right hand side of our test mesh can now be described as

```julia
julia> faces = [
julia> facets = [
(3,6),
(6,9)
]
```

The local ID can be constructed based on elements, corresponding faces and chosen interpolation, since the face ordering is interpolation dependent.
Now to a full example. Let us reconstruct the local facets of a grid boundary with the decsribed machinery to see how everything
fits together:
```julia
julia> function compute_faceset(cells, global_faces, ip::Interpolation{dim}) where {dim}
local_faces = Ferrite.faces(ip)
nodes_per_face = length(local_faces[1])
d = Dict{NTuple{nodes_per_face, Int}, FaceIndex}()
julia> function recompute_facetset_on_boundary(cells, global_facets)
# Facets in 2D are edges and have two vertices
d = Dict{NTuple{2, Int}, FacetIndex}()
for (c, cell) in enumerate(cells) # c is global cell number
for (f, face) in enumerate(local_faces) # f is local face number
# store the global nodes for the particular element, local face combination
d[ntuple(i-> cell.nodes[face[i]], nodes_per_face)] = FaceIndex(c, f)
# Grab all local facets
local_facets = Ferrite.reference_facets(cell)
for (f, facet) in enumerate(local_facets) # f is local facet number
# Translate into global nodes
global_facet = ntuple(i-> cell.nodes[facet[i]], length(facet))
d[global_facet] = FacetIndex(c, f)
end
end

faces = Vector{FaceIndex}()
for face in global_faces
# lookup the element, local face combination for this face
push!(faces, d[face])
facets = Vector{FacetIndex}()
for facet in global_facets
# lookup the element, local facet combination for this facet
push!(facets, d[facet])
end

return faces
return facets
end

julia> interpolation = Lagrange{2, RefCube, 1}()

julia> compute_faceset(cells, faces, interpolation)
Vector{FaceIndex} with 2 elements:
FaceIndex((2, 2))
FaceIndex((4, 2))
julia> recompute_facetset_on_boundary(cells, facets)
Vector{FacetIndex} with 2 elements:
FacetIndex((2, 2))
FacetIndex((4, 2))
```

Ferrite considers edges only in the three dimensional space. However, they share the concepts of faces in terms of `(global_cell_id,local_edge_id)` identifier.
Note that this example can be simplified by directly using `Ferrite.facets`:

```julia
julia> function recompute_facetset_on_boundary_alternative(cells, global_facets)
# Facets in 2D are edges and have two vertices
d = Dict{NTuple{2, Int}, FacetIndex}()
for (c, cell) in enumerate(cells) # c is global cell number
# Grab all local facets
for (f, global_facet) in enumerate(Ferrite.facets(cell)) # f is local facet number
d[global_facet] = FacetIndex(c, f)
end
end

facets = Vector{FacetIndex}()
for facet in global_facets
# lookup the element, local facet combination for this facet
push!(facets, d[facet])
end

return facets
end

julia> recompute_facetset_on_boundary_alternative(cells, facets)
Vector{FacetIndex} with 2 elements:
FacetIndex((2, 2))
FacetIndex((4, 2))
```

## AbstractGrid

Expand Down Expand Up @@ -182,6 +209,6 @@ In order to use boundaries, e.g. for Dirichlet constraints in the ConstraintHand

## Topology

Ferrite.jl's `Grid` type offers experimental features w.r.t. topology information. The functions [`Ferrite.getneighborhood`](@ref) and [`Ferrite.faceskeleton`](@ref)
Ferrite.jl's `Grid` type offers experimental features w.r.t. topology information. The functions [`Ferrite.getneighborhood`](@ref) and [`Ferrite.facetskeleton`](@ref)
are the interface to obtain topological information. The [`Ferrite.getneighborhood`](@ref) can construct lists of directly connected entities based on a given entity (`CellIndex,FaceIndex,EdgeIndex,VertexIndex`).
The [`Ferrite.faceskeleton`](@ref) function can be used to evaluate integrals over material interfaces or computing element interface values such as jumps.
The [`Ferrite.facetskeleton`](@ref) function can be used to evaluate integrals over material interfaces or computing element interface values such as jumps.
36 changes: 34 additions & 2 deletions src/Grid/grid.jl
Original file line number Diff line number Diff line change
Expand Up @@ -55,6 +55,16 @@
nfaces( ::Type{T}) where {T <: AbstractRefShape} = length(reference_faces(T))
nfacets( ::Type{T}) where {T <: AbstractRefShape} = length(reference_facets(T))


"""
reference_vertices(::AbstractRefShape)
termi-official marked this conversation as resolved.
Show resolved Hide resolved
reference_vertices(::AbstractCell)

Returns a tuple of integers containing the ordered local vertex indices corresponding to
the corners or endpoints of an element.
termi-official marked this conversation as resolved.
Show resolved Hide resolved
"""
reference_vertices(::Union{AbstractRefShape, AbstractCell})
termi-official marked this conversation as resolved.
Show resolved Hide resolved

"""
Ferrite.vertices(::AbstractCell)

Expand All @@ -64,6 +74,20 @@
"""
vertices(::AbstractCell)

"""
reference_edges(::AbstractRefShape)
termi-official marked this conversation as resolved.
Show resolved Hide resolved
reference_edges(::AbstractCell)

Returns a tuple of 2-tuples containing the ordered local node indices corresponding to
the vertices that define an *oriented edge*.

An *oriented edge* is an edge with the first node having the lowest local index and the other
node the secon dindex.
KnutAM marked this conversation as resolved.
Show resolved Hide resolved

Note that the vertices are sufficient to define a face uniquely.
KnutAM marked this conversation as resolved.
Show resolved Hide resolved
"""
reference_edges(::Union{AbstractRefShape, AbstractCell})
termi-official marked this conversation as resolved.
Show resolved Hide resolved

"""
Ferrite.edges(::AbstractCell)

Expand All @@ -77,6 +101,7 @@

"""
reference_faces(::AbstractRefShape)
termi-official marked this conversation as resolved.
Show resolved Hide resolved
reference_faces(::AbstractCell)

Returns a tuple of n-tuples containing the ordered local node indices corresponding to
the vertices that define an *oriented face*.
Expand All @@ -86,7 +111,7 @@

Note that the vertices are sufficient to define a face uniquely.
"""
reference_faces(::AbstractRefShape)
reference_faces(::Union{AbstractRefShape, AbstractCell})
termi-official marked this conversation as resolved.
Show resolved Hide resolved

"""
Ferrite.faces(::AbstractCell)
Expand Down Expand Up @@ -119,6 +144,7 @@

"""
Ferrite.reference_facets(::Type{<:AbstractRefShape})
Ferrite.reference_facets(::AbstractCell)

Returns a tuple of n-tuples containing the ordered local node indices corresponding to
the vertices that define an oriented facet.
Expand All @@ -131,13 +157,19 @@
@inline reference_facets(refshape::Type{<:AbstractRefShape{2}}) = reference_edges(refshape)
@inline reference_facets(refshape::Type{<:AbstractRefShape{3}}) = reference_faces(refshape)

@inline reference_faces(::AbstractCell{refshape}) where refshape = reference_faces(refshape)
@inline reference_edges(::AbstractCell{refshape}) where refshape = reference_edges(refshape)
@inline reference_vertices(::AbstractCell{refshape}) where refshape = reference_vertices(refshape)
KnutAM marked this conversation as resolved.
Show resolved Hide resolved
@inline reference_facets(::AbstractCell{refshape}) where refshape = reference_facets(refshape)

Check warning on line 163 in src/Grid/grid.jl

View check run for this annotation

Codecov / codecov/patch

src/Grid/grid.jl#L160-L163

Added lines #L160 - L163 were not covered by tests
termi-official marked this conversation as resolved.
Show resolved Hide resolved

"""
geometric_interpolation(::AbstractCell)::Interpolation
geometric_interpolation(::Type{AbstractCell})::Interpolation

Each `AbstractCell` type has a unique geometric interpolation describing its geometry.
This function returns that interpolation.
"""
geometric_interpolation(::AbstractCell)
geometric_interpolation(::T) where T <: AbstractCell = geometric_interpolation(T)

Check warning on line 172 in src/Grid/grid.jl

View check run for this annotation

Codecov / codecov/patch

src/Grid/grid.jl#L172

Added line #L172 was not covered by tests
KnutAM marked this conversation as resolved.
Show resolved Hide resolved
KnutAM marked this conversation as resolved.
Show resolved Hide resolved

"""
Ferrite.get_node_ids(c::AbstractCell)
Expand Down
Loading