Gridap.ReferenceFEs

Gridap.ReferenceFEsModule

The exported names are

source
Gridap.ReferenceFEs.ExtrusionPolytopeMethod
ExtrusionPolytope(extrusion::Int...)

Generates an ExtrusionPolytope from the tuple extrusion. The values in extrusion are either equal to the constant HEX_AXIS or the constant TET_AXIS.

Examples

Creating a quadrilateral, a triangle, and a wedge

using Gridap.ReferenceFEs

quad = ExtrusionPolytope(HEX_AXIS,HEX_AXIS)

tri = ExtrusionPolytope(TET_AXIS,TET_AXIS)

wedge = ExtrusionPolytope(TET_AXIS,TET_AXIS,HEX_AXIS)

println(quad == QUAD)
println(tri == TRI)
println(wedge == WEDGE)

# output
true
true
true
source
Gridap.ReferenceFEs.GenericQuadratureType
struct GenericQuadrature{D,T,C <: AbstractVector{Point{D,T}},W <: AbstractVector{T}} <: Quadrature{D,T}
  coordinates::Vector{Point{D,T}}
  weights::Vector{T}
  name::String
end
source
Gridap.ReferenceFEs.GenericRefFEType
struct GenericRefFE{T,D} <: ReferenceFE{D}
  # + private fields
end

This type is a materialization of the ReferenceFE interface. That is, it is a struct that stores the values of all abstract methods in the ReferenceFE interface. This type is useful to build reference FEs from the underlying ingredients without the need to create a new type.

Note that some fields in this struct are type unstable deliberately in order to simplify the type signature. Don't access them in computationally expensive functions, instead extract the required fields before and pass them to the computationally expensive function.

source
Gridap.ReferenceFEs.GenericRefFEMethod
  GenericRefFE{T}(
  ndofs::Int,
  polytope::Polytope{D},
  prebasis::AbstractVector{<:Field},
  dofs::AbstractVector{<:Dof},
  conformity::Conformity,
  metadata,
  face_dofs::Vector{Vector{Int}},
  shapefuns::AbstractVector{<:Field}=compute_shapefuns(dofs,prebasis)) where {T,D}

Constructs a GenericRefFE object with the provided data.

source
Gridap.ReferenceFEs.GenericRefFEMethod
  GenericRefFE{T}(
  ndofs::Int,
  polytope::Polytope{D},
  prebasis::AbstractVector{<:Field},
  predofs::AbstractVector{<:Dof},
  conformity::Conformity,
  metadata,
  face_dofs::Vector{Vector{Int}},
  shapefuns::AbstractVector{<:Field},
  dofs::AbstractVector{<:Dof}=compute_dofs(predofs,shapefuns)) where {T,D}

Constructs a GenericRefFE object with the provided data.

source
Gridap.ReferenceFEs.LagrangianDofBasisType
struct LagrangianDofBasis{P,V} <: AbstractArray{<:Dof}
  nodes::Vector{P}
  dof_to_node::Vector{Int}
  dof_to_comp::Vector{Int}
  node_and_comp_to_dof::Vector{V}
end

Type that implements a Lagrangian dof basis.

Fields:

  • nodes::Vector{P} vector of points (P<:Point) storing the nodal coordinates
  • node_and_comp_to_dof::Vector{V} vector such that node_and_comp_to_dof[node][comp] returns the dof associated with node node and the component comp in the type V.
  • dof_to_node::Vector{Int} vector of integers such that dof_to_node[dof] returns the node id associated with dof id dof.
  • dof_to_comp::Vector{Int} vector of integers such that dof_to_comp[dof] returns the component id associated with dof id dof.
source
Gridap.ReferenceFEs.LagrangianDofBasisMethod
LagrangianDofBasis(::Type{T},nodes::Vector{<:Point}) where T

Creates a LagrangianDofBasis for fields of value type T associated with the vector of nodal coordinates nodes.

source
Gridap.ReferenceFEs.LagrangianRefFEType
abstract type LagrangianRefFE{D} <: ReferenceFE{D}

Abstract type representing a Lagrangian reference FE. Lagrangian in the sense that get_dof_basis returns an instance of LagrangianDofBasis. The interface for this type is defined with the methods of ReferenceFE plus the following ones

source
Gridap.ReferenceFEs.LagrangianRefFEMethod
LagrangianRefFE(::Type{T},p::Polytope,orders) where T
LagrangianRefFE(::Type{T},p::Polytope,order::Int) where T

Builds a LagrangianRefFE object on top of the given polytope. T is the type of the value of the approximation space (e.g., T=Float64 for scalar-valued problems, T=VectorValue{N,Float64} for vector-valued problems with N components). The arguments order or orders are for the polynomial order of the resulting space, which allows isotropic or anisotropic orders respectively (provided that the cell topology allows the given anisotropic order). The argument orders should be an indexable collection of D integers (e.g., a tuple or a vector), being D the number of space dimensions.

In order to be able to use this function, the type of the provided polytope p has to implement the following additional methods. They have been implemented for ExtrusionPolytope in the library. They need to be implemented for new polytope types in order to build Lagangian reference elements on top of them.

The following methods are also used in the construction of the LagrangianRefFE object. A default implementation of them is available in terms of the three previous methods. However, the user can also implement them for new polytope types increasing customization possibilities.

source
Gridap.ReferenceFEs.LinearCombinationDofVectorType
struct LinearCombinationDofVector{T} <: AbstractVector{Dof}
  change_of_basis::Matrix{T}
  dof_basis::AbstractVector{<:Dof}
end

Type that implements a dof basis (a) as the linear combination of a dof basis (b). The dofs are first evaluated at dof basis (b) (field dof_basis) and the dof values are next mapped to dof basis (a) applying a change of basis (field change_of_basis).

Fields:

  • change_of_basis::Matrix{T} the matrix of the change from dof basis (b) to (a)
  • dof_basis::AbstractVector{<:Dof} A type representing dof basis (b)
source
Gridap.ReferenceFEs.PolytopeType
abstract type Polytope{D} <: GridapType

Abstract type representing a polytope (i.e., a polyhedron in arbitrary dimensions). D is the environment dimension (typically, 0, 1, 2, or 3). This type parameter is needed since there are functions in the Polytope interface that return containers with Point{D} objects. We adopt the usual nomenclature for polytope-related objects. All objects in a polytope (from vertices to the polytope itself) are called n-faces or simply faces. The notation n-faces is used only when it is needed to refer to the object dimension n. Otherwise we simply use face. In addition, we say

  • vertex (pl. vertices): for 0-faces
  • edge: for 1-faces
  • facet: for (D-1)-faces

The Polytope interface is defined by overloading the following functions

And optionally these ones:

The interface can be tested with the function

source
Gridap.ReferenceFEs.PolytopeMethod
Polytope{N}(p::Polytope,faceid::Integer) where N

Returns a Polytope{N} object representing the "reference" polytope of the N-face with id faceid. The value faceid refers to the numeration restricted to the dimension N (it starts with 1 for the first N-face).

source
Gridap.ReferenceFEs.ReferenceFEType
abstract type ReferenceFE{D} <: GridapType

Abstract type representing a Reference finite element. D is the underlying coordinate space dimension. We follow the Ciarlet definition. A reference finite element is defined by a polytope (cell topology), a basis of an interpolation space of top of this polytope (denoted here as the prebasis), and a basis of the dual of this space (i.e. the degrees of freedom). From this information one can compute the shape functions (i.e, the canonical basis of w.r.t. the degrees of freedom) with a simple change of basis. In addition, we also encode in this type information about how the interpolation space in a reference finite element is "glued" with neighbors in order to build conforming cell-wise spaces.

The ReferenceFE interface is defined by overloading these methods:

The interface is tested with

source
Base.:==Method
(==)(a::LagrangianRefFE{D}, b::LagrangianRefFE{D}) where D
source
Base.:==Method
(==)(a::Polytope{D},b::Polytope{D}) where D

Returns true if the polytopes a and b are equivalent. Otherwise, it returns false. Note that the operator == returns false by default for polytopes of different dimensions. Thus, this function has to be overloaded only for the case of polytopes a and b of same dimension.

source
Gridap.ReferenceFEs.BDMRefFEMethod

BDMRefFE(::Type{et},p::Polytope,order::Integer) where et

The order argument has the following meaning: the divergence of the functions in this basis is in the P space of degree order-1.

source
Gridap.ReferenceFEs.ModalC0RefFEMethod

ModalC0RefFE(::Type{T},p::Polytope{D},orders) where {T,D}

Returns an instance of GenericRefFE{ModalC0} representing a ReferenceFE with Modal C0-continuous shape functions (multivariate scalar-valued, vector-valued, or tensor-valued, iso- or aniso-tropic).

For more details about the shape functions, see Section 1.1.5 in

Ern, A., & Guermond, J. L. (2013). Theory and practice of finite elements (Vol. 159). Springer Science & Business Media.

and references therein.

The constructor is only implemented for for n-cubes and the minimum order in any dimension must be greater than one. The DoFs are numbered by n-faces in the same way as with CLagrangianRefFEs.

source
Gridap.ReferenceFEs.NedelecRefFEMethod
NedelecRefFE(::Type{et},p::Polytope,order::Integer) where et

The order argument has the following meaning: the curl of the functions in this basis is in the Q space of degree order.

source
Gridap.ReferenceFEs.RaviartThomasRefFEMethod
RaviartThomasRefFE(::Type{et},p::Polytope,order::Integer) where et

The order argument has the following meaning: the divergence of the functions in this basis is in the Q space of degree order.

source
Gridap.ReferenceFEs.SerendipityRefFEMethod
SerendipityRefFE(::Type{T},p::Polytope,order::Int) where T
SerendipityRefFE(::Type{T},p::Polytope,orders::Tuple) where T

Returns an instance of LagrangianRefFE, whose underlying approximation space is the serendipity space of order order. Implemented for order from 1 to 4. The type of the polytope p has to implement all the queries detailed in the constructor LagrangianRefFE(::Type{T},p::Polytope{D},orders) where {T,D}.

Examples

using Gridap.ReferenceFEs

order = 2
reffe = SerendipityRefFE(Float64,QUAD,order)

println( num_dofs(reffe) )

# output
8
source
Gridap.ReferenceFEs.compute_dofsMethod
compute_dofs(predofs,shapefuns)

Helper function used to compute the dof basis associated with the dof basis predofs and the basis shapefuns.

It is equivalent to

change = inv(evaluate(predofs,shapefuns))
linear_combination(change,predofs) # i.e. transpose(change)*predofs
source
Gridap.ReferenceFEs.compute_face_ordersMethod
compute_face_orders(p::Polytope,face::Polytope,iface::Int,orders)

Returns a vector or a tuple with the order per direction at the face face of the polytope p when restricting the order per direction orders to this face. iface is the face id of face in the numeration restricted to the face dimension.

source
Gridap.ReferenceFEs.compute_lagrangian_reffacesMethod
compute_lagrangian_reffaces(::Type{T},p::Polytope,orders) where T

Returns a tuple of length D being the number of space dimensions. The entry d+1 of this tuple contains a vector of LagrangianRefFE one for each face of dimension d on the boundary of the polytope.

source
Gridap.ReferenceFEs.compute_monomial_basisMethod
compute_monomial_basis(::Type{T},p::Polytope,orders) where T -> MonomialBasis

Returns the monomial basis of value type T and order per direction described by orders on top of the polytope p.

source
Gridap.ReferenceFEs.compute_nodesMethod
compute_nodes(p::Polytope,orders)

When called

node_coords, face_own_nodes = compute_nodes(p,orders)

Returns node_coords, the nodal coordinates of all the Lagrangian nodes associated with the order per direction orders, and face_own_nodes, being a vector of vectors indicating which nodes are owned by each of the faces of the polytope p.

source
Gridap.ReferenceFEs.compute_own_nodesMethod
compute_own_nodes(p::Polytope{D},orders) where D -> Vector{Point{D,Float64}}

Returns the coordinates of the nodes owned by the interior of the polytope associated with a Lagrangian space with the order per direction described by orders.

source
Gridap.ReferenceFEs.compute_shapefunsMethod
compute_shapefuns(dofs,prebasis)

Helper function used to compute the shape function basis associated with the dof basis dofs and the basis prebasis.

It is equivalent to

change = inv(evaluate(dofs,prebasis))
linear_combination(change,prebasis) # i.e. transpose(change)*prebasis
source
Gridap.ReferenceFEs.get_dimrangesMethod
get_dimranges(p::Polytope) -> Vector{UnitRange{Int}}

Given a polytope p it returns a vector of ranges. The entry d+1 in this vector contains the range of face ids for the faces of dimension d.

Examples

using Gridap.ReferenceFEs

ranges = get_dimranges(SEGMENT)
println(ranges)

# output
UnitRange{Int}[1:2, 3:3]

Face ids for the vertices in the segment range from 1 to 2 (2 vertices), the face ids for edges in the segment range from 3 to 3 (only one edge with id 3).

source
Gridap.ReferenceFEs.get_edge_tangentMethod
get_edge_tangent(p::Polytope) -> Vector{VectorValue{D,Float64}}

Given a polytope p, returns a vector of VectorValue objects representing the unit tangent vectors to the polytope edges.

source
Gridap.ReferenceFEs.get_face_dofsMethod
get_face_dofs(reffe::ReferenceFE) -> Vector{Vector{Int}}

Returns a vector of vector that, for each face, stores the dofids in the closure of the face.

source
Gridap.ReferenceFEs.get_face_typeMethod
get_face_type(p::Polytope,d::Integer) -> Vector{Int}

Return a vector of integers denoting, for each face of dimension d, an index to the vector get_reffaces(Polytope{d},p)

Examples

Get the unique polytopes for the facets of a wedge and identify of which type each face is.

using Gridap.ReferenceFEs

reffaces = get_reffaces(Polytope{2},WEDGE)

face_types = get_face_type(WEDGE,2)

println(reffaces)
println(face_types)

# output
Gridap.ReferenceFEs.ExtrusionPolytope{2}[TRI, QUAD]
[1, 1, 2, 2, 2]

The three first facets are of type 1, i.e, QUAD, and the last ones of type 2, i.e., TRI.

source
Gridap.ReferenceFEs.get_facedimsMethod
get_facedims(p::Polytope) -> Vector{Int}

Given a polytope p, returns a vector indicating the dimension of each face in the polytope

Examples

using Gridap.ReferenceFEs

dims = get_facedims(SEGMENT)
println(dims)

# output
[0, 0, 1]

The first two faces in the segment (the two vertices) have dimension 0 and the third face (the segment itself) has dimension 1

source
Gridap.ReferenceFEs.get_facesMethod
get_faces(p::Polytope,dimfrom::Integer,dimto::Integer) -> Vector{Vector{Int}}

For dimfrom >= dimto returns a vector that for each face of dimension dimfrom stores a vector of the ids of faces of dimension dimto on its boundary.

For dimfrom < dimto returns a vector that for each face of dimfrom stores a vector of the face ids of faces of dimension dimto that touch it.

The numerations used in this function are the ones restricted to each dimension.

using Gridap.ReferenceFEs

edge_to_vertices = get_faces(QUAD,1,0)
println(edge_to_vertices)

vertex_to_edges_around = get_faces(QUAD,0,1)
println(vertex_to_edges_around)

# output
Array{Int,1}[[1, 2], [3, 4], [1, 3], [2, 4]]
Array{Int,1}[[1, 3], [1, 4], [2, 3], [2, 4]]
source
Gridap.ReferenceFEs.get_facesMethod
get_faces(p::Polytope) -> Vector{Vector{Int}}

Given a polytope p the function returns a vector of vectors defining the incidence relation of the faces in the polytope.

Each face in the polytope receives a unique integer id. The id 1 is assigned to the first 0-face. Consecutive increasing ids are assigned to the other 0-faces, then to 1-faces, and so on. The polytope itself receives the largest id which coincides with num_faces(p). For a face id iface, get_faces(p)[iface] is a vector of face ids, corresponding to the faces that are incident with the face labeled with iface. That is, faces that are either on its boundary or the face itself. In this vector of incident face ids, faces are ordered by dimension, starting with 0-faces. Within each dimension, the labels are ordered in a consistent way with the polyope object for the face iface itself.

Examples

using Gridap.ReferenceFEs

faces = get_faces(SEGMENT)
println(faces)

# output
Array{Int,1}[[1], [2], [1, 2, 3]]

The constant SEGMENT is bound to a predefined instance of polytope that represents a segment. The face labels associated with a segment are [1,2,3], being 1 and 2 for the vertices and 3 for the segment itself. In this case, this function returns the vector of vectors [[1],[2],[1,2,3]] meaning that vertex 1 is incident with vertex 1 (idem for vertex 2), and that the segment (id 3) is incident with the vertices 1 and 2 and the segment itself.

source
Gridap.ReferenceFEs.get_facet_normalMethod
get_facet_normal(p::Polytope) -> Vector{VectorValue{D,Float64}}

Given a polytope p, returns a vector of VectorValue objects representing the unit outward normal vectors to the polytope facets.

source
Gridap.ReferenceFEs.get_facet_orientationsMethod
get_facet_orientations(p::Polytope) -> Vector{Int}

Given a polytope p returns a vector of integers of length num_facets(p). Facets, whose vertices are ordered consistently with the outwards normal vector, receive value 1 in this vector. Otherwise, facets receive value -1.

source
Gridap.ReferenceFEs.get_offsetsMethod
get_offsets(p::Polytope) -> Vector{Int}

Given a polytope p, it returns a vector of integers. The position in the d+1 entry in this vector is the offset that transforms a face id in the global numeration in the polytope to the numeration restricted to faces to dimension d.

Examples

using Gridap.ReferenceFEs

offsets = get_offsets(SEGMENT)
println(offsets)

# output
[0, 2]
source
Gridap.ReferenceFEs.get_reffacesMethod
get_reffaces(::Type{Polytope{d}},p::Polytope) where d -> Vector{Polytope{d}}

Get a vector of the unique polytopes for the faces of dimension d.

Examples

Get the unique polytopes for the facets of a wedge.

using Gridap.ReferenceFEs

reffaces = get_reffaces(Polytope{2},WEDGE)

println(reffaces)

# output
Gridap.ReferenceFEs.ExtrusionPolytope{2}[TRI, QUAD]
source
Gridap.ReferenceFEs.get_reffacesMethod
get_reffaces(
  ::Type{ReferenceFE{d}},
  reffe::GenericLagrangianRefFE{GradConformity}) where d -> Vector{GenericLagrangianRefFE{GradConformity,M,d}}
source
Gridap.ReferenceFEs.get_shapefunsMethod
get_shapefuns(reffe::ReferenceFE) -> Field

Returns the basis of shape functions (i.e. the canonical basis) associated with the reference FE. The result is encoded as a Field object.

source
Gridap.ReferenceFEs.get_vertex_permutationsMethod
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.

source
Gridap.ReferenceFEs.next_corner!Method

Sweep through the dim-dimensional hypercube recursively, collecting all simplices.

We represent vertices as bit patterns. In dim dimensions, the lowermost dim bits are either zero or one. Interpreted as integer, this labels the vertices of the hypercube from the origin ("bottom") 0 to the diagonally opposite vertex ("top") 2^dim-1.

Each simplex contains both the bottom vertex 0 as well as the top vertex 2^dim-1. Its other vertices trace a path from the bottom to the top. The algorithm below finds all possible paths.

  • simplices is the accumulator where the simplices are collected.
  • vertices is the current set of vertices as we sweep from the origin to the diagonally opposide vertex.
  • corner is the current corner.
source
Gridap.ReferenceFEs.test_dofMethod
test_dof(dof,field,v;cmp::Function=(==))

Test that the Dof interface is properly implemented for object dof. It also checks if the object dof when evaluated at the field field returns the same value as v. Comparison is made with the comp function.

source
Gridap.ReferenceFEs.test_polytopeMethod
test_polytope(p::Polytope{D}; optional::Bool=false) where D

Function that stresses out the functions in the Polytope interface. It tests whether the function in the polytope interface are defined for the given object, and whether they return objects of the expected type. With optional=false (the default), only the mandatory functions are checked. With optional=true, the optional functions are also tested.

source