Gridap.Geometry
DiscreteModels
In Gridap, a DiscreteModel is the main object representing a discretized domain, where the finite element problem is defined adn solved. It has three main components:
- A
GridTopology, which defines the topology of the mesh. That is the polytopes that make up the mesh as well as their connectivity. This object also provides the connectivity between all d-dimensional entities of the mesh, such as vertices, edges, faces, and cells. - A
Grid, which defines the geometric properties of the mesh. That is, the coordinates of the vertices and the geometric map between the reference space of each polytope and the physical space. - A
FaceLabeling, which classifies all the d-dimensional entities of the mesh into (non-overlapping) physical entities. These entities are then associated one or multiple tags which can be used to impose boundary conditions or define physical regions on the mesh.
To get a deeper understanding of these three components, we encourage the reader to check the low-level tutorial on geometry. The API for each of these components is documented in the following sections.
GridTopology API
Gridap.Geometry.GridTopology — Typeabstract type GridTopology{Dc,Dp}Abstract type representing the topological information associated with a grid.
The GridTopology interface is defined by overloading the methods:
get_faces(g::GridTopology,dimfrom::Integer,dimto::Integer)get_polytopes(g::GridTopology)get_cell_type(g::GridTopology)get_vertex_coordinates(g::GridTopology)
The GridTopology interface has the following traits
and tested with this function:
Gridap.Geometry.GridTopology — MethodGridTopology(grid::Grid)
GridTopology(grid::Grid, cell_to_vertices::Table, vertex_to_node::Vector)Gridap.Geometry.OrientationStyle — MethodOrientationStyle(::Type{<:GridTopology})
OrientationStyle(::GridTopology)Oriented() if has oriented faces, NonOriented() otherwise (default).
Gridap.Geometry.RegularityStyle — MethodRegularityStyle(::Type{<:GridTopology})
RegularityStyle(::GridTopology)Regular() if no hanging-nodes default), Irregular() otherwise.
Gridap.Geometry.compute_cell_faces — Methodcompute_cell_faces(g::GridTopology)Gridap.Geometry.compute_face_vertices — Methodcompute_face_vertices(g::GridTopology)Gridap.Geometry.compute_reffaces — Methodcompute_reffaces(g::GridTopology)Gridap.Geometry.compute_reffaces — Methodcompute_reffaces(::Type{Polytope{d}}, g::GridTopology) where dGridap.Geometry.get_cell_faces — Methodget_cell_faces(g::GridTopology)Defaults to
compute_cell_faces(g)Gridap.Geometry.get_cell_permutations — Methodget_cell_permutations(top::GridTopology)
get_cell_permutations(top::GridTopology,d::Integer)Returns an cell-wise array of permutations. For each cell, the entry contains the permutations of the local d-faces of the cell w.r.t the global d-faces of the mesh.
If no target dimension is specified, all dimensions are returned.
Gridap.Geometry.get_cell_polytopes — Methodget_cell_polytopes(topo::GridTopology)Gridap.Geometry.get_cell_type — Methodget_cell_type(g::GridTopology)Gridap.Geometry.get_cell_vertices — Methodget_cell_vertices(g::GridTopology)Gridap.Geometry.get_isboundary_face — Methodget_isboundary_face(g::GridTopology)
get_isboundary_face(g::GridTopology,d::Integer)Returns a vector of booleans indicating if the face is a boundary face. Boundary faces are defined in the following way:
- If
d = D-1, i.e facets, a boundary facet is a facet that is adjacent to only one cell. - Otherwise, a face is a boundary face if it belongs to a boundary facet.
If no target dimension is specified, all dimensions are returned.
Gridap.Geometry.get_polytopes — Methodget_polytopes(g::GridTopology)Gridap.Geometry.get_reffaces_offsets — Methodget_reffaces_offsets(topo::GridTopology)Gridap.Geometry.is_oriented — Methodis_oriented(::Type{<:GridTopology}) -> Bool
is_oriented(a::GridTopology) -> BoolGridap.Geometry.is_regular — Methodis_regular(::Type{<:GridTopology}) -> Bool
is_regular(a::GridTopology) -> BoolGridap.Geometry.num_cells — Methodnum_cells(g::GridTopology)Gridap.Geometry.restrict — Methodrestrict(topo::GridTopology, cell_to_parent_cell::AbstractVector{<:Integer})
restrict(topo::GridTopology, parent_cell_to_mask::AbstractVector{Bool})Gridap.Geometry.test_grid_topology — Methodtest_grid_topology(top::GridTopology)Gridap.ReferenceFEs.get_dimrange — Methodget_dimrange(g::GridTopology,d::Integer)Gridap.ReferenceFEs.get_dimranges — Methodget_dimranges(g::GridTopology)Gridap.ReferenceFEs.get_face_coordinates — Methodget_face_coordinates(g::GridTopology)
get_face_coordinates(g::GridTopology,d::Integer)Gridap.ReferenceFEs.get_face_type — Methodget_face_type(g::GridTopology,d::Integer)By default, it calls to compute_reffaces.
Gridap.ReferenceFEs.get_face_type — Methodget_face_type(topo::GridTopology)Gridap.ReferenceFEs.get_face_vertices — Methodget_face_vertices(g::GridTopology,d::Integer)Gridap.ReferenceFEs.get_face_vertices — Methodget_face_vertices(g::GridTopology)Defaults to
compute_face_vertices(g)Gridap.ReferenceFEs.get_facedims — Methodget_facedims(g::GridTopology)Gridap.ReferenceFEs.get_faces — Methodget_faces(g::GridTopology,dimfrom::Integer,dimto::Integer)Gridap.ReferenceFEs.get_offset — Methodget_offset(g::GridTopology,d::Integer)Gridap.ReferenceFEs.get_offsets — Methodget_offsets(g::GridTopology)Gridap.ReferenceFEs.get_reffaces — Methodget_reffaces(topo::GridTopology)Gridap.ReferenceFEs.get_reffaces — Methodget_reffaces(::Type{Polytope{d}}, g::GridTopology) where dBy default, it calls to compute_reffaces.
Gridap.ReferenceFEs.get_vertex_coordinates — Methodget_vertex_coordinates(g::GridTopology)Gridap.ReferenceFEs.is_n_cube — Methodis_n_cube(p::GridTopology) -> BoolGridap.ReferenceFEs.is_simplex — Methodis_simplex(p::GridTopology) -> BoolGridap.ReferenceFEs.num_cell_dims — Methodnum_cell_dims(::GridTopology) -> Int
num_cell_dims(::Type{<:GridTopology}) -> IntGridap.ReferenceFEs.num_dims — Methodnum_dims(::GridTopology) -> Int
num_dims(::Type{<:GridTopology}) -> IntEquivalent to num_cell_dims.
Gridap.ReferenceFEs.num_edges — Methodnum_edges(g::GridTopology)Gridap.ReferenceFEs.num_faces — Methodnum_faces(g::GridTopology,d::Integer)
num_faces(g::GridTopology)Gridap.ReferenceFEs.num_facets — Methodnum_facets(g::GridTopology)Gridap.ReferenceFEs.num_point_dims — Methodnum_point_dims(::GridTopology) -> Int
num_point_dims(::Type{<:GridTopology}) -> IntGridap.ReferenceFEs.num_vertices — Methodnum_vertices(g::GridTopology)Grid API
Gridap.Geometry.Grid — Typeabstract type Grid{Dc,Dp}Abstract type that represents mesh of a domain of parametric dimension Dc and physical dimension Dp.
The interface of Grid is defined by overloading:
get_node_coordinates(trian::Grid)get_cell_node_ids(trian::Grid)get_reffes(trian::Grid)get_cell_type(trian::Grid)
The Grid interface has the following traits
The interface of Grid is tested with
Gridap.Geometry.Grid — MethodGrid(reffe::LagrangianRefFE)Gridap.Geometry.Grid — MethodGrid(::Type{ReferenceFE{d}},p::Polytope) where dGridap.Geometry.OrientationStyle — MethodOrientationStyle(::Type{<:Grid})
OrientationStyle(::Grid)Oriented() if has oriented faces, NonOriented() otherwise (default).
Gridap.Geometry.RegularityStyle — MethodRegularityStyle(::Type{<:Grid})
RegularityStyle(::Grid)Regular() if no hanging-nodes default), Irregular() otherwise.
Gridap.Geometry.compute_linear_grid — Methodcompute_linear_grid(reffe::LagrangianRefFE)Gridap.Geometry.compute_reference_grid — Methodcompute_reference_grid(p::LagrangianRefFE, nelems::Integer)Gridap.Geometry.compute_reference_grid — Methodcompute_reference_grid(p::Polytope,nelems)Gridap.Geometry.compute_reference_grid — Methodcompute_reference_grid(p::LagrangianRefFE, partition::NTuple{D,Integer})Gridap.Geometry.get_cell_map — Methodget_cell_map(trian::Grid) -> Vector{<:Field}Gridap.Geometry.get_cell_node_ids — Methodget_cell_node_ids(trian::Grid)Gridap.Geometry.get_cell_polytopes — Methodget_cell_polytopes(trian::Grid) -> Vector{<:Polytope}Gridap.Geometry.get_cell_ref_coordinates — MethodGridap.Geometry.get_cell_reffe — Methodget_cell_reffe(trian::Grid) -> Vector{<:LagrangianRefFE}It is not desirable to iterate over the resulting array for large number of cells if the underlying reference FEs are of different Julia type.
Gridap.Geometry.get_cell_shapefuns — Methodget_cell_shapefuns(trian::Grid) -> Vector{<:Field}Gridap.Geometry.get_cell_type — Methodget_cell_type(trian::Grid) -> AbstractVector{<:Integer}Gridap.Geometry.get_reffes — Methodget_reffes(trian::Grid) -> Vector{LagrangianRefFE}Gridap.Geometry.is_oriented — Methodis_oriented(::Type{<:Grid}) -> Bool
is_oriented(a::Grid) -> BoolGridap.Geometry.is_regular — Methodis_regular(::Type{<:Grid}) -> Bool
is_regular(a::Grid) -> BoolGridap.Geometry.num_cells — Methodnum_cells(trian::Grid) -> IntGridap.Geometry.test_grid — Methodtest_grid(trian::Grid)Gridap.ReferenceFEs.get_facet_normal — Methodget_facet_normal(trian::Grid)Gridap.ReferenceFEs.get_node_coordinates — Methodget_node_coordinates(trian::Grid) -> AbstractArray{<:Point{Dp}}Gridap.ReferenceFEs.is_first_order — Methodis_first_order(trian::Grid) -> BoolGridap.ReferenceFEs.num_cell_dims — Methodnum_cell_dims(::Grid) -> Int
num_cell_dims(::Type{<:Grid}) -> IntGridap.ReferenceFEs.num_dims — Methodnum_dims(::Grid) -> Int
num_dims(::Type{<:Grid}) -> IntEquivalent to num_cell_dims.
Gridap.ReferenceFEs.num_nodes — Methodnum_nodes(trian::Grid) -> IntGridap.ReferenceFEs.num_point_dims — Methodnum_point_dims(::Grid) -> Int
num_point_dims(::Type{<:Grid}) -> IntGridap.ReferenceFEs.simplexify — Methodsimplexify(grid::Grid;kwargs...)FaceLabeling API
Gridap.Geometry.FaceLabeling — Typestruct FaceLabeling <: GridapType
d_to_dface_to_entity::Vector{Vector{Int32}}
tag_to_entities::Vector{Vector{Int32}}
tag_to_name::Vector{String}
endGridap.Geometry.FaceLabeling — MethodGridap.Geometry.FaceLabeling — MethodFaceLabeling(d_to_num_dfaces::Vector{Int})
FaceLabeling(topo::GridTopology)
FaceLabeling(topo::GridTopology,cell_to_tag::Vector{<:Integer},tag_to_name::Vector{String})Base.merge! — MethodBase.merge!(a::FaceLabeling,b::FaceLabeling...)Gridap.Geometry.add_tag! — Methodadd_tag!(lab::FaceLabeling,name::String,entities::Vector{<:Integer})Gridap.Geometry.add_tag_from_tag_filter! — Methodadd_tag_from_tag_filter!(lab::FaceLabeling, name::String, filter::Function)Adds a new tag name, by including all entities selected by a filter function. The filter function must have signature
filter(entity_tags::Vector{Int}) -> Boolwhere entity_tags are the tags of a particular geometrical entity.
Gridap.Geometry.add_tag_from_tags! — Methodadd_tag_from_tags!(lab::FaceLabeling, name::String, tags::Vector{Int})
add_tag_from_tags!(lab::FaceLabeling, name::String, tags::Vector{String})
add_tag_from_tags!(lab::FaceLabeling, name::String, tag::Int)
add_tag_from_tags!(lab::FaceLabeling, name::String, tag::String)Adds a new tag name, given by the union of the tags tags.
Gridap.Geometry.add_tag_from_tags_complementary! — Methodadd_tag_from_tags_complementary!(lab::FaceLabeling, name::String, tags::Vector{Int})
add_tag_from_tags_complementary!(lab::FaceLabeling, name::String, tags::Vector{String})
add_tag_from_tags_complementary!(lab::FaceLabeling, name::String, tag::Int)
add_tag_from_tags_complementary!(lab::FaceLabeling, name::String, tag::String)Adds a new tag name, given by the complementary of the tags tags.
Gridap.Geometry.add_tag_from_tags_intersection! — Methodadd_tag_from_tags_intersection!(labels::FaceLabeling, name::String, tags::Vector{Int})
add_tag_from_tags_intersection!(labels::FaceLabeling, name::String, tags::Vector{String})Adds a new tag name, given by the intersection of the tags tags.
Gridap.Geometry.add_tag_from_tags_setdiff! — Methodadd_tag_from_tags_setdiff!(lab::FaceLabeling, name::String, tags_include::Vector{Int}, tags_exclude::Vector{Int})
add_tag_from_tags_setdiff!(lab::FaceLabeling, name::String, tags_include::Vector{String}, tags_exclude::Vector{String})Adds a new tag name, given by the set difference between tags_include and tags_exclude.
Gridap.Geometry.face_labeling_from_cell_tags — Methodface_labeling_from_cell_tags(topo::GridTopology,cell_to_tag::Vector{<:Integer},tag_to_name::Vector{String})Creates a FaceLabeling object from a GridTopology and a vector containing a tag for each cell.
You can use this function to add custom cell-defined tags to a pre-existing FaceLabeling, by calling merge! on the two FaceLabeling objects, i.e
labels = FaceLabeling(topo)
additional_labels = face_labeling_from_cell_tags(topo,cell_to_tag,tag_to_name)
merge!(labels,additional_labels)Gridap.Geometry.face_labeling_from_vertex_filter — Methodface_labeling_from_vertex_filter(topo::GridTopology, name::String, filter::Function)Creates a FaceLabeling object from a GridTopology and a coordinate-based filter function. The filter function takes in vertex coordinates and returns a boolean values. A geometrical entity is tagged if all its vertices pass the filter.
You can use this function to add custom vertex-defined tags to a pre-existing FaceLabeling, by calling merge! on the two FaceLabeling objects, i.e
labels = FaceLabeling(topo)
additional_labels = face_labeling_from_vertex_filter(topo,name,filter)
merge!(labels,additional_labels)Gridap.Geometry.get_cell_entity — MethodGridap.Geometry.get_face_entity — Methodget_face_entity(lab::FaceLabeling,d::Integer)Gridap.Geometry.get_face_entity — Methodget_face_entity(lab::FaceLabeling)Gridap.Geometry.get_face_mask — Methodget_face_mask(labeling::FaceLabeling,tags::Vector{Int},d::Integer)
get_face_mask(labeling::FaceLabeling,tags::Vector{String},d::Integer)
get_face_mask(labeling::FaceLabeling,tag::Int,d::Integer)
get_face_mask(labeling::FaceLabeling,tag::String,d::Integer)Gridap.Geometry.get_face_tag — Methodget_face_tag(labeling::FaceLabeling,tags::Vector{Int},d::Integer)
get_face_tag(labeling::FaceLabeling,tags::Vector{String},d::Integer)
get_face_tag(labeling::FaceLabeling,tag::Int,d::Integer)
get_face_tag(labeling::FaceLabeling,tag::String,d::Integer)
get_face_tag(labeling::FaceLabeling,d::Integer)The first of the given tags appearing in the face is taken. If there is no tag on a face, this face will have a value equal to UNSET. If not tag or tags are provided, all the tags in the model are considered
Gridap.Geometry.get_face_tag_index — Methodget_face_tag_index(labeling::FaceLabeling,tags::Vector{Int},d::Integer)
get_face_tag_index(labeling::FaceLabeling,tags::Vector{String},d::Integer)
get_face_tag_index(labeling::FaceLabeling,tag::Int,d::Integer)
get_face_tag_index(labeling::FaceLabeling,tag::String,d::Integer)Like get_face_tag by provides the index into the array tags instead of the tag stored in tags.
Gridap.Geometry.get_tag_entities — Methodget_tag_entities(lab::FaceLabeling,tag::Integer)
get_tag_entities(lab::FaceLabeling,tag::String)Gridap.Geometry.get_tag_entities — Methodget_tag_entities(lab::FaceLabeling)Gridap.Geometry.get_tag_from_name — Methodget_tag_from_name(lab::FaceLabeling,name::String)Gridap.Geometry.get_tag_from_name — Methodget_tag_from_name(lab::FaceLabeling)Gridap.Geometry.get_tag_name — Methodget_tag_name(lab::FaceLabeling,tag::Integer)Gridap.Geometry.get_tag_name — Methodget_tag_name(lab::FaceLabeling)Gridap.Geometry.get_tags_from_names — Methodget_tags_from_names(lab::FaceLabeling,names::Vector{String})Gridap.Geometry.num_cells — Methodnum_cells(lab::FaceLabeling)Gridap.Geometry.num_entities — Methodnum_entities(lab::FaceLabeling)Gridap.Geometry.num_tags — Methodnum_tags(lab::FaceLabeling)Gridap.Geometry.restrict — Methodrestrict(labels::FaceLabeling,d_to_dface_to_parent_dface)Gridap.ReferenceFEs.num_cell_dims — Methodnum_cell_dims(lab::FaceLabeling)Gridap.ReferenceFEs.num_dims — Methodnum_dims(lab::FaceLabeling)Gridap.ReferenceFEs.num_edges — Methodnum_edges(lab::FaceLabeling)Gridap.ReferenceFEs.num_faces — Methodnum_faces(lab::FaceLabeling,d::Integer)Gridap.ReferenceFEs.num_faces — Methodnum_faces(lab::FaceLabeling)Gridap.ReferenceFEs.num_facets — Methodnum_facets(lab::FaceLabeling)Gridap.ReferenceFEs.num_vertices — Methodnum_vertices(lab::FaceLabeling)DiscreteModel API
Gridap.Geometry.DiscreteModel — Typeabstract type DiscreteModel{Dc,Dp} <: GridAbstract type holding information about a physical grid, the underlying grid topology, and a labeling of the grid faces. This is the information that typically provides a mesh generator, and it is what one needs to perform a simulation.
The DiscreteModel interface is defined by overloading the methods:
get_grid(model::DiscreteModel)get_grid_topology(model::DiscreteModel)get_face_labeling(g::DiscreteModel)
The interface is tested with this function:
Gridap.Geometry.DiscreteModel — MethodGridap.Geometry.DiscreteModel — MethodGridap.Geometry.Grid — MethodGrid(::Type{ReferenceFE{d}},model::DiscreteModel) where dGridap.Geometry.GridTopology — MethodGridap.Geometry.DiscreteModelFromFile — MethodDiscreteModelFromFile(filename::AbstractString)Gridap.Geometry.compute_face_nodes — Methodcompute_face_nodes(model::DiscreteModel,d::Integer)Gridap.Geometry.compute_face_nodes — Methodcompute_face_nodes(model::DiscreteModel)Gridap.Geometry.compute_face_own_nodes — Methodcompute_face_own_nodes(model::DiscreteModel,d::Integer)Gridap.Geometry.compute_face_own_nodes — Methodcompute_face_own_nodes(model::DiscreteModel)Gridap.Geometry.compute_node_face_owner — Methodcompute_node_face_owner(g::DiscreteModel)Gridap.Geometry.compute_reffaces — Methodcompute_reffaces(g::DiscreteModel)Gridap.Geometry.compute_reffaces — Methodcompute_reffaces(::Type{ReferenceFE{d}}, g::DiscreteModel) where dGridap.Geometry.compute_vertex_node — Methodcompute_vertex_node(g::DiscreteModel)Gridap.Geometry.get_face_labeling — Methodget_face_labeling(g::DiscreteModel)Gridap.Geometry.get_grid — Methodget_grid(model::DiscreteModel)Gridap.Geometry.get_grid_topology — Methodget_grid_topology(model::DiscreteModel)Gridap.Geometry.get_node_face_owner — Methodget_node_face_owner(g::DiscreteModel)Gridap.Geometry.get_polytopes — Methodget_polytopes(model::DiscreteModel)Gridap.Geometry.get_reffaces_offsets — Methodget_reffaces_offsets(model::DiscreteModel)Gridap.Geometry.num_cells — Methodnum_cells(g::DiscreteModel)Gridap.Geometry.test_discrete_model — Methodtest_discrete_model(model::DiscreteModel)Gridap.ReferenceFEs.get_face_nodes — Methodget_face_nodes(g::DiscreteModel,d::Integer)Gridap.ReferenceFEs.get_face_nodes — Methodget_face_nodes(g::DiscreteModel)Gridap.ReferenceFEs.get_face_own_nodes — Methodget_face_own_nodes(g::DiscreteModel,d::Integer)Gridap.ReferenceFEs.get_face_own_nodes — Methodget_face_own_nodes(g::DiscreteModel)Gridap.ReferenceFEs.get_face_type — Methodget_face_type(g::DiscreteModel,d::Integer)Index to the vector get_reffaces(ReferenceFE{d},g)
Gridap.ReferenceFEs.get_face_type — Methodget_face_type(model::DiscreteModel)Gridap.ReferenceFEs.get_reffaces — Methodget_reffaces(model::DiscreteModel)Gridap.ReferenceFEs.get_reffaces — Methodget_reffaces(::Type{ReferenceFE{d}},model::DiscreteModel) where dGridap.ReferenceFEs.get_vertex_node — Methodget_vertex_node(g::DiscreteModel)Gridap.ReferenceFEs.num_cell_dims — Methodnum_cell_dims(model::DiscreteModel)Gridap.ReferenceFEs.num_dims — Methodnum_dims(model::DiscreteModel)Gridap.ReferenceFEs.num_edges — Methodnum_edges(g::DiscreteModel)Gridap.ReferenceFEs.num_faces — Methodnum_faces(g::DiscreteModel,d::Integer)
num_faces(g::DiscreteModel)Gridap.ReferenceFEs.num_facets — Methodnum_facets(g::DiscreteModel)Gridap.ReferenceFEs.num_nodes — Methodnum_nodes(g::DiscreteModel)Gridap.ReferenceFEs.num_point_dims — Methodnum_point_dims(model::DiscreteModel)Gridap.ReferenceFEs.num_vertices — Methodnum_vertices(g::DiscreteModel)Gridap.ReferenceFEs.simplexify — Methodsimplexify(model::DiscreteModel)Types of models
UnstructuredDiscreteModels
These are the main types of models used for unstructured meshes. They have their own type of topology and grid.
Gridap.Geometry.UnstructuredDiscreteModel — Typestruct UnstructuredDiscreteModel{Dc,Dp,Tp,B} <: DiscreteModel{Dc,Dp}
grid::UnstructuredGrid{Dc,Dp,Tp,B}
grid_topology::UnstructuredGridTopology{Dc,Dp,Tp,B}
face_labeling::FaceLabeling
endGridap.Geometry.UnstructuredDiscreteModel — MethodGridap.Geometry.UnstructuredDiscreteModel — MethodUnstructuredDiscreteModel(grid::Grid)Gridap.Geometry.UnstructuredGrid — Typestruct UnstructuredGrid{Dc,Dp,Tp,Ti,O} <: Grid{Dc,Dp}
node_coordinates::Vector{Point{Dp,Tp}}
cell_node_ids::Table{Ti,Int32}
reffes::Vector{<:LagrangianRefFE{Dc}}
cell_types::Vector{Int8}
endGridap.Geometry.UnstructuredGrid — MethodUnstructuredGrid(x::AbstractArray{<:Point})Gridap.Geometry.UnstructuredGrid — MethodUnstructuredGrid(grid::Grid;kwargs...)Gridap.Geometry.UnstructuredGrid — MethodUnstructuredGrid(reffe::LagrangianRefFE)Build a grid with a single cell that is the given reference FE itself
Gridap.Geometry.UnstructuredGrid — Methodfunction UnstructuredGrid(
node_coordinates::Vector{Point{Dp,Tp}},
cell_node_ids::Table{Ti},
reffes::Vector{<:LagrangianRefFE{Dc}},
cell_types::Vector,
orientation_style::OrientationStyle=NonOriented()) where {Dc,Dp,Tp,Ti}
endLow-level inner constructor.
Gridap.Geometry.UnstructuredGrid — MethodUnstructuredGrid(::Type{ReferenceFE{d}},p::Polytope) where dGridap.Geometry.UnstructuredGridTopology — TypeUnstructuredGridTopology(
vertex_coordinates::Vector{<:Point},
d_to_dface_vertices::Vector{<:Table},
cell_type::Vector{<:Integer},
polytopes::Vector{<:Polytope},
orientation_style::OrientationStyle=NonOriented())Gridap.Geometry.UnstructuredGridTopology — TypeUnstructuredGridTopology(
vertex_coordinates::Vector{<:Point},
cell_vertices::Table,
cell_type::Vector{<:Integer},
polytopes::Vector{<:Polytope},
orientation_style::OrientationStyle=NonOriented())Gridap.Geometry.UnstructuredGridTopology — Typestruct UnstructuredGridTopology{Dc,Dp,T,O} <: GridTopology{Dc,Dp}
# private fields
endGridap.Geometry.UnstructuredGridTopology — MethodUnstructuredGridTopology(topo::GridTopology)Gridap.Geometry.UnstructuredGridTopology — MethodUnstructuredGridTopology(grid::UnstructuredGrid)
UnstructuredGridTopology(grid::UnstructuredGrid,cell_to_vertices::Table,vertex_to_node::AbstractVector)CartesianDiscreteModels
Cartesian models are specific structures to represent cartesian domains, possibly mapped by a prescribed function. They offer optimizations with respect to the unstructured models and are quite lighweight memory-wise. They implement their own type of grid and can be converted to unstructured models if needed.
Gridap.Geometry.CartesianDiscreteModel — Typestruct CartesianDiscreteModel{D,T,F} <: DiscreteModel{D,D}
# Private Fields
endGridap.Geometry.CartesianDiscreteModel — MethodCartesianDiscreteModel(args...)Same args needed to construct a CartesianDescriptor
Gridap.Geometry.CartesianDiscreteModel — MethodCartesianDiscreteModel(desc::CartesianDescriptor{D,T,F}, cmin::CartesianIndex, cmax::CartesianIndex)
Builds a CartesianDiscreteModel object which represents a subgrid of a (larger) grid represented by desc. This subgrid is described by its D-dimensional minimum (cmin) and maximum (cmax) CartesianIndex identifiers.
Gridap.Geometry.CartesianDiscreteModel — MethodCartesianDiscreteModel(desc::CartesianDescriptor)
Gridap.Geometry._find_ncube_face_neighbor_deltas — Methodfindncubefaceneighbor_deltas(p::ExtrusionPolytope{D}) -> Vector{CartesianIndex}
Given an n-cube type ExtrusionPolytope{D}, returns V=Vector{CartesianIndex} with as many entries as the number of faces in the boundary of the Polytope. For an entry facelid in this vector, V[facelid] returns what has to be added to the CartesianIndex of a cell in order to obtain the CartesianIndex of the cell neighbour of K across the face F with local ID face_lid.
Gridap.Geometry.get_cartesian_descriptor — Methodget_cartesian_descriptor(model::CartesianDiscreteModel)Gridap.Geometry.CartesianDescriptor — Typestruct CartesianDescriptor{D,T,F<:Function}
origin::Point{D,T}
sizes::NTuple{D,T}
partition::NTuple{D,Int}
map::F
endStruct that stores the data defining a Cartesian grid.
Gridap.Geometry.CartesianDescriptor — MethodCartesianDescriptor(
domain,
partition;
map::Function=identity,
isperiodic::NTuple{D,Bool}=tfill(false,Val{D}))domain and partition are 1D indexable collections of arbitrary type.
Gridap.Geometry.CartesianDescriptor — MethodCartesianDescriptor(
origin::Point{D},
sizes::NTuple{D},
partition;
map::Function=identity,
isperiodic::NTuple{D,Bool}=tfill(false,Val{D})) where Dpartition is a 1D indexable collection of arbitrary type.
Gridap.Geometry.CartesianDescriptor — MethodCartesianDescriptor(
pmin::Point{D},
pmax::Point{D},
partition;
map::Function=identity,
isperiodic::NTuple{D,Bool}=tfill(false,Val{D})) where Dpartition is a 1D indexable collection of arbitrary type.
Gridap.Geometry.CartesianGrid — Typestruct CartesianGrid{D,T,F} <: Grid{D,D}
# private fields
endGridap.Geometry.CartesianGrid — MethodCartesianGrid(args...;kwargs...)Same args needed to construct a CartesianDescriptor
Gridap.Geometry.CartesianGrid — MethodCartesianGrid(desc::CartesianDescriptor)Gridap.Geometry.get_cartesian_descriptor — Methodget_cartesian_descriptor(grid::CartesianGrid)Get the descriptor of the Cartesian grid
PolytopalDiscreteModels
Polytopal models are used to represent meshes made of arbitrarily shaped polytopes. They are directly build on the physical domain and therefore do require a geometric map. They can only be used with polytopal methods such as DG, HDG or HHO. More expensive than regular unstructured models, but also more flexible. They implement their own type of grid and topology.
Gridap.Geometry.PolytopalDiscreteModel — Typestruct PolytopalDiscreteModel{Dc,Dp,Tp,Tn} <: DiscreteModel{Dc,Dp}
grid::PolytopalGrid{Dc,Dp,Tp,Tn}
grid_topology::PolytopalGridTopology{Dc,Dp,Tp}
labels::FaceLabeling
endDiscrete model for polytopal grids.
Constructors:
PolytopalDiscreteModel(model::DiscreteModel)
PolytopalDiscreteModel(grid::PolytopalGrid,grid_topology::PolytopalGridTopology,labels::FaceLabeling)Gridap.Geometry.PolytopalGrid — Typestruct PolytopalGrid{Dc,Dp,Tp,Tn} <: Grid{Dc,Dp}
node_coordinates::Vector{Point{Dp,Tp}}
cell_node_ids::Table{Int32,Vector{Int32},Vector{Int32}}
polytopes::Vector{GeneralPolytope{Dc,Dp,Tp,Nothing}}
facet_normal::Tn
endGrid for polytopal meshes.
Constructors:
PolytopalGrid(grid::Grid)
PolytopalGrid(topo::PolytopalGridTopology)
PolytopalGrid(node_coordinates,cell_node_ids,polytopes[,facet_normal=nothing])Gridap.Geometry.PolytopalGridTopology — Typestruct PolytopalGridTopology{Dc,Dp,Tp} <: GridTopology{Dc,Dp}
vertex_coordinates::Vector{Point{Dp,Tp}}
n_m_to_nface_to_mfaces::Matrix{Table{Int32,Vector{Int32},Vector{Int32}}}
polytopes::Vector{GeneralPolytope{Dc,Dp,Tp,Nothing}}
endGrid topology for polytopal grids.
Constructors:
PolytopalGridTopology(topo::UnstructuredGridTopology)
PolytopalGridTopology(vertex_coordinates,cell_vertices,polytopes)Gridap.Geometry.voronoi — Methodvoronoi(model::DiscreteModel) -> PolytopalDiscreteModel
voronoi(topo::GridTopology) -> PolytopalGridTopologyGiven a mesh, computes it's associated Voronoi mesh. NOTE: Only working in 2D.
Other models
Gridap.Geometry.DiscreteModelPortion — TypeGridap.Geometry.DiscreteModelPortion — MethodGridap.Geometry.restrict — Methodrestrict(model::DiscreteModel, cell_to_parent_cell::AbstractVector{<:Integer})
restrict(model::DiscreteModel, parent_cell_to_mask::AbstractArray{Bool})Gridap.Geometry.MappedDiscreteModel — TypeMappedDiscreteModel
Represent a model with a MappedGrid grid. See also MappedGrid.
Gridap.Geometry.GridPortion — Typestruct GridPortion{Dc,Dp,G} <: Grid{Dc,Dp}
parent::G
cell_to_parent_cell::Vector{Int32}
node_to_parent_node::Vector{Int32}
endGridap.Geometry.GridPortion — MethodGridPortion(parent::Grid{Dc,Dp},cell_to_parent_cell::Vector{Int32}) where {Dc,Dp}Gridap.Geometry.restrict — Methodrestrict(grid::Grid, cell_to_parent_cell::AbstractVector{<:Integer})
restrict(grid::Grid, parent_cell_to_mask::AbstractArray{Bool})Triangulations
Given a model, which holds information for all dimensions and on the whole domain, we can take d-dimensional slices of it (or part of it) where we will define our finite element spaces and/or integrate our weakforms. These slices are called Triangulations.
Body-Fitted triangulations
The most basic type of triangulation is the BodyFittedTriangulation. It represents a d-dimensional set of faces attached to d-dimensional faces of the model. In particular, it is the type used to represent bulk triangulations (which contain a subset of the cells of the model).
Gridap.Geometry.Triangulation — Typeabstract type Triangulation{Dt,Dp}A discredited physical domain associated with a DiscreteModel{Dm,Dp}.
Dt and Dm can be different.
The (mandatory) Triangulation interface can be tested with
Gridap.Geometry.best_target — Methodbest_target(trian1::Triangulation,trian2::Triangulation)If possible, returns a Triangulation to which CellDatum objects can be transferred from trian1 and trian2. Can be trian1, trian2 or a new Triangulation.
Gridap.Geometry.is_change_possible — Methodis_change_possible(strian::Triangulation,ttrian::Triangulation)Returns true if CellDatum objects can be transferred from strian to ttrian.
Gridap.Geometry.pos_neg_data — MethodGridap.Geometry.test_triangulation — Methodtest_triangulation(trian::Triangulation)Boundary and Skeleton triangulations
To perform integration of bulk variables on a set of faces (for instance the integration of Neumann terms), we somehow need to link these faces to one of their neighboring cells (otherwise the bulk variables are not uniquely valued on the faces). This is done through the BoundaryTriangulation object. It is generally used for integrals on the physical boundary (where faces have a single neighboring cell), but offer quite a lot more flexibility. Conceptually, a BoundaryTriangulation is given by a BodyFittedTriangulation of faces and a FaceToCellGlue that links each face to the selected neighboring cell. To perform integration on interior faces, which have two neighboring cells, we can use the SkeletonTriangulation object. Each SkeletonTriangulation is a wrapper around two BoundaryTriangulations, one for each side of the face. It provides straighforward ways to performs jumps and means over the interior faces.
Gridap.Geometry.BoundaryTriangulation — TypeGridap.Geometry.BoundaryTriangulation — MethodBoundaryTriangulation(model::DiscreteModel,face_to_mask::Vector{Bool})
BoundaryTriangulation(model::DiscreteModel)Gridap.Geometry.BoundaryTriangulation — MethodBoundaryTriangulation(model::DiscreteModel,labeling::FaceLabeling;tags::Vector{Int})
BoundaryTriangulation(model::DiscreteModel,labeling::FaceLabeling;tags::Vector{String})
BoundaryTriangulation(model::DiscreteModel,labeling::FaceLabeling;tag::Int)
BoundaryTriangulation(model::DiscreteModel,labeling::FaceLabeling;tag::String)Gridap.Geometry.BoundaryTriangulation — MethodBoundaryTriangulation(model::DiscreteModel,tags::Vector{Int})
BoundaryTriangulation(model::DiscreteModel,tags::Vector{String})
BoundaryTriangulation(model::DiscreteModel,tag::Int)
BoundaryTriangulation(model::DiscreteModel,tag::String)Gridap.Geometry.SkeletonPair — Typestruct SkeletonPair{L,R} <: GridapType
plus::L
minus::R
endGridap.Geometry.SkeletonTriangulation — Typestruct SkeletonTriangulation{Dc,Dp,B} <: Triangulation{Dc,Dp}
plus::B
minus::B
endGridap.Geometry.SkeletonTriangulation — MethodSkeletonTriangulation(model::DiscreteModel,face_to_mask::Vector{Bool})
SkeletonTriangulation(model::DiscreteModel)Gridap.Geometry.InterfaceTriangulation — MethodPatch triangulations
We provide an API to integrate and solve problems on patches of cells. This API revolves around two objects: The PatchTopology defines the topology of the patches, and holds information on all the d-dimensional entities belonging to each patch. One can then take d-dimensional slices of these patches using PatchTriangulations.
Gridap.Geometry.PatchTopology — Typestruct PatchTopology{Dc,Dp} <: GridapType
topo :: GridTopology{Dc,Dp}
d_to_patch_to_dfaces :: Vector{Table{Int32,Vector{Int32},Vector{Int32}}}
endFields:
topo: Underlying grid topologyd_to_patch_to_dfaces: For each dimensiond, a table mapping each patch to it's d-faces.
Gridap.Geometry.PatchTriangulation — Typestruct PatchTriangulation{Dc,Dp} <: Triangulation{Dc,Dp}
trian :: Triangulation{Dc,Dp}
ptopo :: PatchTopology
glue :: PatchGlue{Dc}
endWrapper around a Triangulation, for patch-based assembly.
Fields:
trian: Underlying triangulation. In general, this can be a non-injective triangulation.ptopo: Patch topology, to which the triangulation faces can be mappedglue: Patch glue, mapping triangulation faces to the patches
Gridap.Geometry.PatchTriangulation — MethodPatchTriangulation(model::DiscreteModel,ptopo::PatchTopology;tags=nothing)
PatchTriangulation(::Type{ReferenceFE{D}},model::DiscreteModel,ptopo::PatchTopology;tags=nothing)Gridap.Geometry.PatchBoundaryTriangulation — MethodPatchBoundaryTriangulation(model::DiscreteModel,ptopo::PatchTopology{Dc};tags=nothing)Gridap.Geometry.PatchSkeletonTriangulation — MethodPatchSkeletonTriangulation(model::DiscreteModel,ptopo::PatchTopology{Dc};tags=nothing)Gridap.Geometry.get_patch_faces — Methodget_patch_faces(ptopo::PatchTopology,d::Integer) -> patch_facesGridap.Geometry.patch_extend — Methodpatch_extend(PD::PatchTopology{Dc},patch_to_data,Df=Dc) -> pface_to_dataGridap.Geometry.patch_reindex — Methodpatch_reindex(ptopo::PatchTopology{Dc},face_to_data,Df=Dc) -> pface_to_data