Extendable grid
An ExtendableGrid in form of a dictionary with types as keys and type stable value access. This means that grid components are accessed as dict entries, e.g. grid[Coordinates] . The rationale of this approach is explained here.
Notations
A grid is assumed to be a subset of components of a polyhedral complex in d-dimensional space. We distinguish the following element classes characterized by their dimension:
| Element class | Meaning |
|---|---|
| Node | 0-dimensional node |
| Edge | 1-dimensional line connecting two neighboring nodes |
| Face | codimension 1 object separating a cell from outer space or neighboring cell |
| Cell | codimension 0 object |
| BFace | Face situated at inner or domain boundary |
| Region | number to be used to characterize subdomains, contacts etc. |
Grid components
Grid components are accessed like Dict entries, the keys must be subtypes of AbstractGridComponent.
Basic set of grid components
Upon construction, an ExtendableGrid needs to be provided with the basic set of grid components denoted by the following component type keys:
| Component type key | Meaning |
|---|---|
Coordinates | Coordinates of the vertices of the grid cells |
CellNodes | Adjacency describing the nodes of grid cell |
CellGeometries | Abstract array of subtypes of AbstractElementGeometry describing the geometry of each cell |
CellRegions | Abstract array of integers describing region numbers |
BFaceNodes | Adjacency structure describing the nodes corresponding to each grid boundary face |
BFaceGeometries | Abstract array of subtypes of AbstractElementGeometry describing the geometry of each boundary face |
BFaceRegions | Abstract array of integers describig region numbers |
CoordinateSystem | Abstract type describing the coordinate system to be used |
Hierarchy of component type keys
The list of components can be printed using the gridcomponents method.
AbstractGridComponent
├─ AbstractElementGeometries
│ ├─ BEdgeGeometries
│ ├─ BFaceGeometries
│ ├─ CellGeometries
│ ├─ EdgeGeometries
│ ├─ FaceGeometries
│ ├─ UniqueBEdgeGeometries
│ ├─ UniqueBFaceGeometries
│ ├─ UniqueCellGeometries
│ ├─ UniqueEdgeGeometries
│ └─ UniqueFaceGeometries
├─ AbstractElementRegions
│ ├─ BEdgeRegions
│ ├─ BFaceRegions
│ ├─ CellRegions
│ ├─ EdgeRegions
│ └─ FaceRegions
├─ AbstractGridAdjacency
│ ├─ BEdgeAssemblyGroups
│ ├─ BEdgeNodes
│ ├─ BFaceAssemblyGroups
│ ├─ BFaceCells
│ ├─ BFaceEdges
│ ├─ BFaceNodes
│ ├─ CellAssemblyGroups
│ ├─ CellEdgeSigns
│ ├─ CellEdges
│ ├─ CellFaceOrientations
│ ├─ CellFaceSigns
│ ├─ CellFaces
│ ├─ CellNodes
│ ├─ EdgeAssemblyGroups
│ ├─ EdgeCells
│ ├─ EdgeNodes
│ ├─ FaceAssemblyGroups
│ ├─ FaceCells
│ ├─ FaceEdgeSigns
│ ├─ FaceEdges
│ ├─ FaceNodes
│ ├─ NodeCells
│ ├─ NodeEdges
│ └─ NodeFaces
├─ AbstractGridFloatArray1D
│ ├─ BEdgeVolumes
│ ├─ BFaceVolumes
│ ├─ CellVolumes
│ ├─ EdgeVolumes
│ ├─ FaceVolumes
│ ├─ XCoordinates
│ ├─ YCoordinates
│ └─ ZCoordinates
├─ AbstractGridFloatArray2D
│ ├─ Coordinates
│ ├─ EdgeTangents
│ ├─ FaceNormals
│ └─ VoronoiFaceCenters
├─ AbstractGridIntegerArray1D
│ ├─ BEdgeEdges
│ ├─ BEdgeParents
│ ├─ BFaceCellPos
│ ├─ BFaceFaces
│ ├─ BFaceParents
│ ├─ CellParents
│ ├─ EdgeParents
│ ├─ FaceParents
│ ├─ NodeInParent
│ ├─ NodeParents
│ ├─ NodePatchGroups
│ ├─ NodePermutation
│ ├─ PColorPartitions
│ ├─ PartitionBFaces
│ ├─ PartitionCells
│ ├─ PartitionEdges
│ └─ PartitionNodes
├─ AbstractGridIntegerArray2D
├─ BFaceNormals
├─ CoordinateSystem
├─ AbstractGridFloatConstant
├─ AbstractGridIntegerConstant
│ ├─ NumBEdgeRegions
│ ├─ NumBFaceRegions
│ └─ NumCellRegions
├─ ParentGrid
└─ ParentGridRelation
├─ RefinedGrid
└─ SubGridAdditional components
Additional components can be added by defining a subtype of AbstractGridComponent or a fitting subtype thereof, and assigning the value to the corresponding Dict entry:
g=simplexgrid([1,2,3,4.0])
abstract type MyComponent <: AbstractGridComponent end
g[MyComponent]=13
show(g)ExtendableGrid{Float64, Int32}(dim=1, nnodes=4, ncells=3, nbfaces=2)Alternatively, component creation can be performed lazily. For this purpose one needs to define an instantiate method:
abstract type NodeCells <: AbstractGridAdjacency end
ExtendableGrids.instantiate(grid, ::Type{NodeCells})=atranspose(grid[CellNodes])
g=simplexgrid([1,2,3,4.0])
show(g[NodeCells])1
1 2
2 3
3Grid API
ExtendableGrids.ElementInfo — Type
const ElementInfo{T}=Union{Vector{T},VectorOfConstants{T}}Union type for element information arrays. If all elements have the same information, it can be stored in an economical form as a VectorOfConstants.
ExtendableGrids.AbstractElementGeometries — Type
abstract type AbstractElementGeometries <: AbstractGridComponentArray of element geometry information.
sourceExtendableGrids.AbstractElementRegions — Type
abstract type AbstractElementRegions <: AbstractGridComponentArray of element region number information.
sourceExtendableGrids.AbstractGridAdjacency — Type
abstract type AbstractGridAdjacency <: AbstractGridComponentAny kind of adjacency between grid components
sourceExtendableGrids.AbstractGridComponent — Type
abstract type AbstractGridComponent <: AbstractExtendableGridApexTypeApex type for grid components.
sourceExtendableGrids.AbstractGridFloatArray1D — Type
abstract type AbstractGridFloatArray1D <: AbstractGridComponent1D Array of floating point data
sourceExtendableGrids.AbstractGridFloatArray2D — Type
abstract type AbstractGridFloatArray2D <: AbstractGridComponent2D Array of floating point data
sourceExtendableGrids.BEdgeRegions — Type
abstract type BEdgeRegions <: AbstractElementRegionsBoundary edge region number per boundary edge
sourceExtendableGrids.BFaceGeometries — Type
Description of boundary face geometries
abstract type BFaceGeometries <: AbstractElementGeometriessourceExtendableGrids.BFaceNodes — Type
abstract type BFaceNodes <: AbstractGridAdjacencyAdjacency describing nodes per grid boundary face
sourceExtendableGrids.BFaceRegions — Type
ExtendableGrids.CellNodes — Type
ExtendableGrids.CellRegions — Type
ExtendableGrids.Coordinates — Type
ExtendableGrids.NumBEdgeRegions — Type
abstract type NumBEdgeRegions <: ExtendableGrids.AbstractGridIntegerConstantNumber of boundary edge regions
sourceExtendableGrids.NumBFaceRegions — Type
abstract type NumBFaceRegions <: ExtendableGrids.AbstractGridIntegerConstantNumber of boundary face regions
sourceExtendableGrids.NumCellRegions — Type
abstract type NumCellRegions <: ExtendableGrids.AbstractGridIntegerConstantNumber of cell regions
sourceBase.delete! — Method
delete!(
grid::ExtendableGrid,
T::Type{<:AbstractGridComponent}
) -> Dict{Type{<:AbstractGridComponent}, Any}
Remove grid component
sourceBase.getindex — Method
Base.getindex(grid::ExtendableGrid,T::Type{<:AbstractGridComponent})Generic method for obtaining grid component.
This method is mutating in the sense that non-existing grid components are created on demand.
Due to the fact that components are stored as Any the return value triggers type instability. To prevent this, specialized methods must be (and are) defined.
sourceBase.haskey — Method
Base.map — Method
map(f,grid)Map function f returning a number onto node coordinates of grid. Returns a vector of length corresponding to the number of nodes of the grid. The function can take either a vector or a numbers as arguments. E.g. for a two-dimensional grid g, both
map(X->X[1]+X[2], g)and
map((x,y)->x+y, g)are possible.
sourceBase.setindex! — Method
setindex!(
grid::ExtendableGrid,
v,
T::Type{<:AbstractGridComponent}
) -> Any
Set new grid component
sourceExtendableGrids.coord_type — Method
ExtendableGrids.dim_grid — Method
ExtendableGrids.dim_space — Method
ExtendableGrids.gridcomponents — Method
gridcomponents()
Print the hierarchy of grid component key types (subtypes of AbstractGridComponent. This includes additionally user defined subptypes.
ExtendableGrids.index_type — Method
ExtendableGrids.instantiate — Function
"Hook" for methods instantiating lazy components.
sourceExtendableGrids.instantiate — Method
instantiate(grid, _::Type{NumBEdgeRegions}) -> Any
Instantiate number of boundary edge regions
sourceExtendableGrids.instantiate — Method
ExtendableGrids.instantiate — Method
ExtendableGrids.isconsistent — Method
isconsistent(grid; warnonly=false)Check consistency of grid: a grid is consistent if
- Grid has no dangling nodes
- ... more to be added
If grid is consistent, return true, otherwise throw an error, or, if warnoly==true, return false.
ExtendableGrids.num_bedgeregions — Method
ExtendableGrids.num_bedges — Method
ExtendableGrids.num_bfaceregions — Method
ExtendableGrids.num_bfaces — Method
ExtendableGrids.num_cellregions — Method
ExtendableGrids.num_cells — Method
ExtendableGrids.num_nodes — Method
ExtendableGrids.seemingly_equal — Method
ExtendableGrids.seemingly_equal — Method
seemingly_equal(grid1, grid2; sort=false, confidence=:fullRecursively check seeming equality of two grids. Seemingly means that long arrays are only compared via random samples.
Keyword args:
sort: if true, sort grid pointsconfidence: Confidence level:- :low : Point numbers etc are the same
- :full : all arrays are equal (besides the coordinate array, the arrays only have to be equal up to permutations)
ExtendableGrids.trim! — Method
trim!(grid; keep)
Remove precomputed grid components for lightweight storage of the grid.
Use the keep list to exclude grid components from trimming.
ExtendableGrids.trim — Method
ExtendableGrids.update! — Method
update!(
grid::ExtendableGrid,
T::Type{<:AbstractGridComponent}
) -> Any
Reinstantiate grid component (only if it exists)
sourceExtendableGrids.veryform — Method
veryform(
grid::ExtendableGrid,
v,
_::Type{<:AbstractGridComponent}
) -> Any
Default veryform method.
"veryform" means "verify and/or transform" and is called to check and possibly transform components to be added to the grid via setindex!.
The default method just passes data through.
sourceExtendableGrids.veryform — Method
veryform(grid::ExtendableGrid{Tc,Ti},v,T::Type{<:AbstractGridAdjacency}) where{Tc,Ti}Check proper type of adjacencies upon insertion
source