Gridap.ReferenceFEs
Gridap.ReferenceFEs — ModuleThe exported names are
DofExtrusionPolytopeGenericNodalRefFEGenericRefFEHEXHEX8HEX_AXISINVALID_PERMLagrangianDofBasisLagrangianRefFENodalReferenceFEPDiscRefFEPYRAMIDPolytopeQUADQUAD4RaviartThomasRefFEReferenceFESEG2SEGMENTSerendipityRefFETETTET4TET_AXISTRITRI3VERTEXVERTEX1WEDGEcompute_face_orderscompute_lagrangian_reffacescompute_monomial_basiscompute_nodescompute_own_nodescompute_own_nodes_permutationscompute_shapefunsdof_cachedof_return_typeevaluate_dofevaluate_dof!evaluate_dof_arrayget_bounding_boxget_dimrangeget_dimrangesget_dof_basisget_dof_to_compget_dof_to_nodeget_edge_tangentsget_extrusionget_face_coordinatesget_face_dimrangesget_face_dofsget_face_nodesget_face_own_dofsget_face_own_dofs_permutationsget_face_own_nodesget_face_own_nodes_permutationsget_face_typeget_face_vertex_permutationsget_face_verticesget_facedimsget_facesget_facet_normalsget_facet_orientationsget_node_and_comp_to_dofget_node_coordinatesget_offsetget_offsetsget_own_dofs_permutationsget_own_nodes_permutationsget_polytopeget_prebasisget_reffacesget_shapefunsget_vertex_coordinatesget_vertex_nodeget_vertex_permutationsis_Pis_Qis_Sis_affineis_first_orderis_n_cubeis_simplexnum_cell_dimsnum_dimsnum_dofsnum_edgesnum_facesnum_facetsnum_nodesnum_point_dimsnum_verticessimplexifytest_doftest_nodal_reference_fetest_polytopetest_reference_fe
Polytopes
Interface
Gridap.ReferenceFEs.Polytope — Typeabstract type Polytope{D} <: GridapTypeAbstract 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
get_faces(p::Polytope)get_dimranges(p::Polytope)Polytope{N}(p::Polytope,faceid::Integer) where Nget_vertex_coordinates(p::Polytope)(==)(a::Polytope{D},b::Polytope{D}) where D
And optionally these ones:
get_edge_tangents(p::Polytope)get_facet_normals(p::Polytope)get_facet_orientations(p::Polytope)get_vertex_permutations(p::Polytope)is_n_cube(p::Polytope)is_simplex(p::Polytope)simplexify(p::Polytope)
The interface can be tested with the function
Gridap.ReferenceFEs.get_faces — Methodget_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{Int64,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.
Gridap.ReferenceFEs.get_dimranges — Methodget_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{Int64}[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).
Gridap.ReferenceFEs.get_dimrange — Methodget_dimrange(p::Polytope,d::Integer)Equivalent to
get_dimranges(p)[d+1]Gridap.ReferenceFEs.Polytope — MethodPolytope{N}(p::Polytope,faceid::Integer) where NReturns 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).
Gridap.ReferenceFEs.get_vertex_coordinates — Methodget_vertex_coordinates(p::Polytope) -> Vector{Point{D,Float64}}Given a polytope p return a vector of points representing containing the coordinates of the vertices.
Base.:== — Method(==)(a::Polytope{D},b::Polytope{D}) where DReturns 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.
Gridap.ReferenceFEs.get_edge_tangents — Methodget_edge_tangents(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.
Gridap.ReferenceFEs.get_facet_normals — Methodget_facet_normals(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.
Gridap.ReferenceFEs.get_facet_orientations — Methodget_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.
Gridap.ReferenceFEs.get_vertex_permutations — Methodget_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{Int64,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.is_simplex — Methodis_simplex(p::Polytope) -> BoolGridap.ReferenceFEs.is_n_cube — Methodis_n_cube(p::Polytope) -> BoolGridap.ReferenceFEs.simplexify — Methodsimplexify(p::Polytope) -> Tuple{Vector{Vector{Int}},Polytope}Gridap.ReferenceFEs.test_polytope — Methodtest_polytope(p::Polytope{D}; optional::Bool=false) where DFunction 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.
Gridap.Integration.num_dims — Methodnum_dims(::Type{<:Polytope{D}}) where D
num_dims(p::Polytope{D}) where DReturns D.
Gridap.ReferenceFEs.num_faces — Methodnum_faces(p::Polytope)Returns the total number of faces in polytope p (from vertices to the polytope itself).
Gridap.ReferenceFEs.num_faces — Methodnum_faces(p::Polytope,dim::Integer)Returns the number of faces of dimension dim in polytope p.
Gridap.ReferenceFEs.num_facets — Methodnum_facets(p::Polytope)Returns the number of facets in the polytope p.
Gridap.ReferenceFEs.num_edges — Methodnum_edges(p::Polytope)Returns the number of edges in the polytope p.
Gridap.ReferenceFEs.num_vertices — Methodnum_vertices(p::Polytope)Returns the number of vertices in the polytope p.
Gridap.ReferenceFEs.get_facedims — Methodget_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
Gridap.ReferenceFEs.get_offsets — Methodget_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]
Gridap.ReferenceFEs.get_offset — Methodget_offset(p::Polytope,d::Integer)Equivalent to get_offsets(p)[d+1].
Gridap.ReferenceFEs.get_faces — Methodget_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 funcitons 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{Int64,1}[[1, 2], [3, 4], [1, 3], [2, 4]]
Array{Int64,1}[[1, 3], [1, 4], [2, 3], [2, 4]]Gridap.ReferenceFEs.get_face_vertices — Methodget_face_vertices(p::Polytope) -> Vector{Vector{Int}}
get_face_vertices(p::Polytope,dim::Integer) -> Vector{Vector{Int}}Gridap.ReferenceFEs.get_face_coordinates — Methodget_face_coordinates(p::Polytope,d::Integer)Gridap.ReferenceFEs.get_face_dimranges — Methodget_face_dimranges(p::Polytope,d::Integer)Gridap.ReferenceFEs.get_reffaces — Methodget_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]
Gridap.ReferenceFEs.get_face_type — Methodget_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.
Gridap.ReferenceFEs.get_bounding_box — Methodget_bounding_box(p::Polytope{D}) where DGridap.ReferenceFEs.get_face_vertex_permutations — Methodget_face_vertex_permutations(p::Polytope)
get_face_vertex_permutations(p::Polytope,d::Integer)Extrusion polytopes
Gridap.ReferenceFEs.ExtrusionPolytope — Typestruct ExtrusionPolytope{D} <: Polytope{D}
extrusion::Point{D,Int}
# + private fields
endConcrete type for polytopes that can be represented with an "extrusion" tuple. The underlying extrusion is available in the field extrusion. Instances of this type can be obtained with the constructors
Gridap.ReferenceFEs.ExtrusionPolytope — MethodExtrusionPolytope(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
Gridap.ReferenceFEs.Polytope — MethodPolytope(extrusion::Int...)Equivalent to ExtrusionPolytope(extrusion...)
Gridap.ReferenceFEs.HEX_AXIS — ConstantConstant to be used in order to indicate a hex-like extrusion axis.
Gridap.ReferenceFEs.TET_AXIS — ConstantConstant to be used in order to indicate a tet-like extrusion axis.
Gridap.ReferenceFEs.get_extrusion — Methodget_extrusion(p::ExtrusionPolytope)Equivalent to p.extrusion.
Pre-defined polytope instances
Gridap.ReferenceFEs.VERTEX — Constantconst VERTEX = Polytope()Gridap.ReferenceFEs.SEGMENT — Constantconst SEGMENT = Polytope(HEX_AXIS)Gridap.ReferenceFEs.TRI — Constantconst TRI = Polytope(TET_AXIS,TET_AXIS)Gridap.ReferenceFEs.QUAD — Constantconst QUAD = Polytope(HEX_AXIS,HEX_AXIS)Gridap.ReferenceFEs.TET — Constantconst TET = Polytope(TET_AXIS,TET_AXIS,TET_AXIS)Gridap.ReferenceFEs.HEX — Constantconst HEX = Polytope(HEX_AXIS,HEX_AXIS,HEX_AXIS)Gridap.ReferenceFEs.WEDGE — Constantconst WEDGE = Polytope(TET_AXIS,TET_AXIS,HEX_AXIS)Gridap.ReferenceFEs.PYRAMID — Constantconst PYRAMID = Polytope(HEX_AXIS,HEX_AXIS,TET_AXIS)Degrees of freedom
Interface
Gridap.ReferenceFEs.Dof — Typeabstract type Dof <: KernelAbstract type representing a degree of freedom (DOF), a basis of DOFs, and related objects. These different cases are distinguished by the return type obtained when evaluating the Dof object on a Field object. See function evaluate_dof! for more details.
The following functions needs to be overloaded
The following functions can be overloaded optionally
The interface is tested with
In most of the cases it is not strictly needed that types that implement this interface inherit from Dof. However, we recommend to inherit from Dof, when possible.
Gridap.ReferenceFEs.evaluate_dof! — Methodevaluate_dof!(cache,dof,field)Evaluates the dof dof with the field field. It can return either an scalar value or an array of scalar values depending the case. The cache object is computed with function dof_cache.
When a mathematical dof is evaluated on a physical field, a scalar number is returned. If either the Dof object is a basis of DOFs, or the Field object is a basis of fields, or both objects are bases, then the returned object is an array of scalar numbers. The first dimensions in the resulting array are for the Dof object and the last ones for the Field object. E.g, a basis of nd DOFs evaluated at physical field returns a vector of nd entries. A basis of nd DOFs evaluated at a basis of nf fields returns a matrix of size (nd,nf).
Gridap.ReferenceFEs.dof_cache — Methoddof_cache(dof,field)Returns the cache needed to call evaluate_dof!(cache,dof,field)
Gridap.ReferenceFEs.dof_return_type — Methoddof_return_type(dof,field)Returns the type for the value obtained with evaluating dof with field.
It defaults to
typeof(evaluate_dof(dof,field))Gridap.ReferenceFEs.test_dof — Methodtest_dof(dof,field,v,comp::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.
Gridap.ReferenceFEs.evaluate_dof — Methodevaluate_dof(dof,field)Equivalent to
cache = dof_cache(dof,field)
evaluate_dof!(cache,dof,field)Gridap.Fields.evaluate — Methodevaluate(dof::Dof,field)Equivalent to evaluate_dof(dof,field).
Working with arrays of DOFs
Gridap.ReferenceFEs.evaluate_dof_array — Methodevaluate_dof_array(dof::AbstractArray,field::AbstractArray)Evaluates the Dof objects in the array dof at the Field objects at the array field element by element.
The result is numerically equivalent to
map(evaluate_dof, dof, field)but it is described with a more memory-friendly lazy type.
Gridap.Fields.evaluate — Methodevaluate(dof::AbstractArray{<:Dof},field::AbstractArray)Equivalent to evaluate_dof_array(dof,field)
Lagrangian dof bases
Gridap.ReferenceFEs.LagrangianDofBasis — Typestruct LagrangianDofBasis{P,V} <: Dof
nodes::Vector{P}
dof_to_node::Vector{Int}
dof_to_comp::Vector{Int}
node_and_comp_to_dof::Vector{V}
endType that implements a Lagrangian dof basis.
Fields:
nodes::Vector{P}vector of points (P<:Point) storing the nodal coordinatesnode_and_comp_to_dof::Vector{V}vector such thatnode_and_comp_to_dof[node][comp]returns the dof associated with nodenodeand the componentcompin the typeV.dof_to_node::Vector{Int}vector of integers such thatdof_to_node[dof]returns the node id associated with dof iddof.dof_to_comp::Vector{Int}vector of integers such thatdof_to_comp[dof]returns the component id associated with dof iddof.
Gridap.ReferenceFEs.LagrangianDofBasis — MethodLagrangianDofBasis(::Type{T},nodes::Vector{<:Point}) where TCreates a LagrangianDofBasis for fields of value type T associated with the vector of nodal coordinates nodes.
Reference Finite Elements
Interface
Gridap.ReferenceFEs.ReferenceFE — Typeabstract type ReferenceFE{D} <: GridapTypeAbstract 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:
num_dofs(reffe::ReferenceFE)get_polytope(reffe::ReferenceFE)get_prebasis(reffe::ReferenceFE)get_dof_basis(reffe::ReferenceFE)get_face_own_dofs(reffe::ReferenceFE)get_face_own_dofs_permutations(reffe::ReferenceFE)get_face_dofs(reffe::ReferenceFE)
The interface is tested with
Gridap.ReferenceFEs.num_dofs — Methodnum_dofs(reffe::ReferenceFE) -> IntReturns the number of DOFs.
Gridap.ReferenceFEs.get_polytope — Methodget_polytope(reffe::ReferenceFE) -> PolytopeReturns the underlying polytope object.
Gridap.ReferenceFEs.get_prebasis — Methodget_prebasis(reffe::ReferenceFE) -> FieldReturns the underlying prebasis encoded as a Field object.
Gridap.ReferenceFEs.get_dof_basis — Methodget_dof_basis(reffe::ReferenceFE) -> DofReturns the underlying dof basis encoded in a Dof object.
Gridap.ReferenceFEs.get_face_own_dofs — Methodget_face_own_dofs(reffe::ReferenceFE) -> Vector{Vector{Int}}Gridap.ReferenceFEs.get_face_own_dofs_permutations — Methodget_face_own_dofs_permutations(reffe::ReferenceFE) -> Vector{Vector{Vector{Int}}}Gridap.ReferenceFEs.get_face_dofs — Methodget_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.
Gridap.ReferenceFEs.INVALID_PERM — ConstantConstant of type Int used to signal that a permutation is not valid.
Gridap.ReferenceFEs.test_reference_fe — Methodtest_reference_fe(reffe::ReferenceFE{D}) where DTest if the methods in the ReferenceFE interface are defined for the object reffe.
Gridap.Integration.num_dims — Methodnum_dims(::Type{<:ReferenceFE{D}}) where D
num_dims(reffe::ReferenceFE{D}) where DReturns D.
Gridap.ReferenceFEs.num_cell_dims — Methodnum_cell_dims(::Type{<:ReferenceFE{D}}) where D
num_cell_dims(reffe::ReferenceFE{D}) where DReturns D.
Gridap.Integration.num_point_dims — Methodnum_point_dims(::Type{<:ReferenceFE{D}}) where D
num_point_dims(reffe::ReferenceFE{D}) where DReturns D.
Gridap.ReferenceFEs.num_faces — Methodnum_faces(reffe::ReferenceFE)
num_faces(reffe::ReferenceFE,d::Integer)Gridap.ReferenceFEs.num_vertices — Methodnum_vertices(reffe::ReferenceFE)Gridap.ReferenceFEs.num_edges — Methodnum_edges(reffe::ReferenceFE)Gridap.ReferenceFEs.num_facets — Methodnum_facets(reffe::ReferenceFE)Gridap.ReferenceFEs.get_face_own_dofs — Methodget_face_own_dofs(reffe::ReferenceFE,d::Integer)Gridap.ReferenceFEs.get_face_own_dofs_permutations — Methodget_face_own_dofs_permutations(reffe::ReferenceFE,d::Integer)Gridap.ReferenceFEs.get_own_dofs_permutations — Methodget_own_dofs_permutations(reffe::ReferenceFE)Gridap.ReferenceFEs.get_face_dofs — Methodget_face_dofs(reffe::ReferenceFE,d::Integer)Gridap.ReferenceFEs.get_shapefuns — Methodget_shapefuns(reffe::ReferenceFE) -> FieldReturns the basis of shape functions (i.e. the canonical basis) associated with the reference FE. The result is encoded as a Field object.
Gridap.ReferenceFEs.compute_shapefuns — Methodcompute_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))
change_basis(prebasis,change)Generic reference elements
Gridap.ReferenceFEs.GenericRefFE — Typestruct GenericRefFE{D} <: ReferenceFE{D}
ndofs::Int
polytope::Polytope{D}
prebasis::Field
dofs::Dof
face_own_dofs::Vector{Vector{Int}}
face_own_dofs_permutations::Vector{Vector{Vector{Int}}}
face_dofs::Vector{Vector{Int}}
shapefuns::Field
endThis 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.
Gridap.ReferenceFEs.GenericRefFE — MethodGenericRefFE(
ndofs::Int,
polytope::Polytope{D},
prebasis::Field,
dofs::Dof,
face_own_dofs::Vector{Vector{Int}},
face_own_dofs_permutations::Vector{Vector{Vector{Int}}},
face_dofs::Vector{Vector{Int}},
shapefuns::Field=compute_shapefuns(dofs,prebasis)) where DConstructs a GenericRefFE object with the provided data.
Node-based reference Finite Elements
Interface
Gridap.ReferenceFEs.NodalReferenceFE — Typeabstract type NodalReferenceFE{D} <: ReferenceFE{D}Abstract type representing a node-based reference FE. We understand a node-based reference FE as one that uses the concept of node to locate dofs on the underlying polytope. Here, nodal-based does not necessary mean an interpolatory reference FE. We only assume that each dof is assigned to a node, whereas several dofs can share a same node in general.
The interface for this type is defined with the methods of ReferenceFE plus the following ones
Gridap.ReferenceFEs.get_node_coordinates — Methodget_node_coordinates(reffe::NodalReferenceFE)Gridap.ReferenceFEs.get_node_and_comp_to_dof — Methodget_node_and_comp_to_dof(reffe::NodalReferenceFE)Gridap.ReferenceFEs.get_face_own_nodes — Methodget_face_own_nodes(reffe::NodalReferenceFE)Gridap.ReferenceFEs.get_face_own_nodes_permutations — Methodget_face_own_nodes_permutations(reffe::NodalReferenceFE)Gridap.ReferenceFEs.get_face_nodes — Methodget_face_nodes(reffe::NodalReferenceFE)Gridap.ReferenceFEs.test_nodal_reference_fe — Functiontest_nodal_reference_fe(reffe::NodalReferenceFE)Gridap.ReferenceFEs.num_nodes — Methodnum_nodes(reffe::NodalReferenceFE)Gridap.ReferenceFEs.get_dof_to_node — Methodget_dof_to_node(reffe::NodalReferenceFE)Gridap.ReferenceFEs.get_own_nodes_permutations — Methodget_own_nodes_permutations(reffe::NodalReferenceFE)Gridap.ReferenceFEs.get_vertex_node — Methodget_vertex_node(reffe::NodalReferenceFE) -> Vector{Int}Gridap.ReferenceFEs.get_face_own_nodes — Methodget_face_own_nodes(reffe::NodalReferenceFE,d::Integer)Gridap.ReferenceFEs.get_face_own_nodes_permutations — Methodget_face_own_nodes_permutations(reffe::NodalReferenceFE,d::Integer)Gridap.ReferenceFEs.get_face_nodes — Methodget_face_nodes(reffe::NodalReferenceFE,d::Integer)GenericNodalRefFE
Gridap.ReferenceFEs.GenericNodalRefFE — Typestruct GenericNodalRefFE{D,T,V} <: NodalReferenceFE{D} reffe::GenericRefFE{D} nodecoordinates::Vector{Point{D,T}} nodeandcomptodof::Vector{V} faceownnodes::Vector{Vector{Int}} faceownnodespermutations::Vector{Vector{Vector{Int}}} face_nodes::Vector{Vector{Int}} end
Lagrangian reference elements
Gridap.ReferenceFEs.LagrangianRefFE — Typestruct LagrangianRefFE{D} <: NodalReferenceFE{D}
# private fields
endType representing a Lagrangian finite element.
For this type
get_dof_basis(reffe)returns aLagrangianDofBasisget_prebasis(reffe)returns aMonomialBasis
Gridap.ReferenceFEs.LagrangianRefFE — MethodLagrangianRefFE(
polytope::Polytope{D},
prebasis::MonomialBasis,
dofs::LagrangianDofBasis,
face_own_nodes::Vector{Vector{Int}},
own_nodes_permutations::Vector{Vector{Int}},
reffaces) where DLow level (inner) constructor of LagrangianRefFE.
Gridap.ReferenceFEs.LagrangianRefFE — MethodLagrangianRefFE(::Type{T},p::Polytope,orders) where T
LagrangianRefFE(::Type{T},p::Polytope,order::Int) where TBuilds 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.
compute_monomial_basis(::Type{T},p::Polytope,orders) where Tcompute_own_nodes(p::Polytope,orders)compute_face_orders(p::Polytope,face::Polytope,iface::Int,orders)
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.
Gridap.ReferenceFEs.compute_monomial_basis — Methodcompute_monomial_basis(::Type{T},p::Polytope,orders) where T -> MonomialBasisReturns the monomial basis of value type T and order per direction described by orders on top of the polytope p.
Gridap.ReferenceFEs.compute_own_nodes — Methodcompute_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.
Gridap.ReferenceFEs.compute_face_orders — Methodcompute_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.
Gridap.ReferenceFEs.compute_nodes — Methodcompute_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.
Gridap.ReferenceFEs.compute_own_nodes_permutations — Methodcompute_own_nodes_permutations(
p::Polytope, own_nodes_coordinates) -> Vector{Vector{Int}}Returns a vector of vectors with the permutations of the nodes owned by the interior of the polytope.
Gridap.ReferenceFEs.compute_lagrangian_reffaces — Methodcompute_lagrangian_reffaces(::Type{T},p::Polytope,orders) where TReturns 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.
Gridap.ReferenceFEs.get_dof_to_comp — Methodget_dof_to_comp(reffe::LagrangianRefFE)Gridap.ReferenceFEs.ReferenceFE — MethodReferenceFE{N}(reffe::LagrangianRefFE,iface::Integer) where NBase.:== — Method(==)(a::LagrangianRefFE{D}, b::LagrangianRefFE{D}) where DGridap.ReferenceFEs.get_reffaces — Methodget_reffaces(
::Type{ReferenceFE{d}},
reffe::LagrangianRefFE) where d -> Vector{LagrangianRefFE{d}}Gridap.ReferenceFEs.get_face_type — Methodget_face_type(reffe::LagrangianRefFE, d::Integer) -> Vector{Int}Gridap.ReferenceFEs.is_first_order — Methodis_first_order(reffe::NodalReferenceFE) -> BoolGridap.ReferenceFEs.is_affine — Methodis_affine(reffe::NodalReferenceFE) -> BoolQuery if the reffe leads to an afine map (true only for first order spaces on top of simplices)
Gridap.ReferenceFEs.is_P — Methodis_P(reffe::LagrangianRefFE)Gridap.ReferenceFEs.is_Q — Methodis_Q(reffe::LagrangianRefFE)
Gridap.ReferenceFEs.is_S — Methodis_S(reffe::LagrangianRefFE)
Serendipity reference elements
Gridap.ReferenceFEs.SerendipityRefFE — FunctionSerendipityRefFE(::Type{T},p::Polytope,order::Int) where T
SerendipityRefFE(::Type{T},p::Polytope,orders::Tuple) where TReturns 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
PDiscRefFE
Gridap.ReferenceFEs.PDiscRefFE — Typestruct PDiscRefFE{D} <: NodalReferenceFE{D}
# Private fields
endRaviartThomasRefFE
Gridap.ReferenceFEs.RaviartThomasRefFE — FunctionRaviartThomasRefFE(::Type{et},p::Polytope,order::Integer) where etPre-defined ReferenceFE instances
Gridap.ReferenceFEs.VERTEX1 — Constantconst VERTEX1 = LagrangianRefFE(Float64,VERTEX,1)Gridap.ReferenceFEs.SEG2 — Constantconst SEG2 = LagrangianRefFE(Float64,SEGMENT,1)Gridap.ReferenceFEs.TRI3 — Constantconst TRI3 = LagrangianRefFE(Float64,TRI,1)Gridap.ReferenceFEs.QUAD4 — Constantconst QUAD4 = LagrangianRefFE(Float64,QUAD,1)Gridap.ReferenceFEs.TET4 — Constantconst TET4 = LagrangianRefFE(Float64,TET,1)Gridap.ReferenceFEs.HEX8 — Constantconst HEX8 = LagrangianRefFE(Float64,HEX,1)