Grid & AbstractGrid
Grid
Ferrite.generate_grid — Functiongenerate_grid(celltype::Cell, nel::NTuple, [left::Vec, right::Vec)Return a Grid for a rectangle in 1, 2 or 3 dimensions. celltype defined the type of cells, e.g. Triangle or Hexahedron. nel is a tuple of the number of elements in each direction. left and right are optional endpoints of the domain. Defaults to -1 and 1 in all directions.
Ferrite.Node — TypeNode{dim, T}A Node is a point in space.
Fields
x::Vec{dim,T}: stores the coordinates
Ferrite.Cell — TypeCell{dim,N,M} <: AbstractCell{dim,N,M}A Cell is a subdomain defined by a collection of Nodes. The parameter dim refers here to the geometrical/ambient dimension, i.e. the dimension of the nodes in the grid and not the topological dimension of the cell (i.e. the dimension of the reference element obtained by default_interpolation). A Cell has N nodes and M faces. Note that a Cell is not defined geometrically by node coordinates, but rather topologically by node indices into the node vector of some grid.
Fields
nodes::Ntuple{N,Int}: N-tuple that stores the node ids. The ordering defines a cell's and its subentities' orientations.
Ferrite.CellIndex — TypeA CellIndex wraps an Int and corresponds to a cell with that number in the mesh
Ferrite.VertexIndex — TypeA VertexIndex wraps an (Int, Int) and defines a local vertex by pointing to a (cell, vert).
Ferrite.EdgeIndex — TypeA EdgeIndex wraps an (Int, Int) and defines a local edge by pointing to a (cell, edge).
Ferrite.FaceIndex — TypeA FaceIndex wraps an (Int, Int) and defines a local face by pointing to a (cell, face).
Ferrite.Grid — TypeGrid{dim, C<:AbstractCell, T<:Real} <: AbstractGrid}A Grid is a collection of Cells and Nodes which covers the computational domain, together with Sets of cells, nodes and faces. There are multiple helper structures to apply boundary conditions or define subdomains. They are gathered in the cellsets, nodesets, facesets, edgesets and vertexsets.
Fields
cells::Vector{C}: stores all cells of the gridnodes::Vector{Node{dim,T}}: stores thedimdimensional nodes of the gridcellsets::Dict{String,Set{Int}}: maps aStringkey to aSetof cell idsnodesets::Dict{String,Set{Int}}: maps aStringkey to aSetof global node idsfacesets::Dict{String,Set{FaceIndex}}: maps aStringto aSetofSet{FaceIndex} (global_cell_id, local_face_id)edgesets::Dict{String,Set{EdgeIndex}}: maps aStringto aSetofSet{EdgeIndex} (global_cell_id, local_edge_idvertexsets::Dict{String,Set{VertexIndex}}: maps aStringkey to aSetof local vertex idsboundary_matrix::SparseMatrixCSC{Bool,Int}: optional, only needed byonboundaryto check if a cell is on the boundary, see, e.g. Helmholtz example
Utility Functions
Ferrite.getcells — Functiongetcells(grid::AbstractGrid)
getcells(grid::AbstractGrid, v::Union{Int,Vector{Int}}
getcells(grid::AbstractGrid, setname::String)Returns either all cells::Collection{C<:AbstractCell} of a <:AbstractGrid or a subset based on an Int, Vector{Int} or String. Whereas the last option tries to call a cellset of the grid. Collection can be any indexable type, for Grid it is Vector{C<:AbstractCell}.
Ferrite.getncells — FunctionReturns the number of cells in the <:AbstractGrid.
Ferrite.getnodes — Functiongetnodes(grid::AbstractGrid)
getnodes(grid::AbstractGrid, v::Union{Int,Vector{Int}}
getnodes(grid::AbstractGrid, setname::String)Returns either all nodes::Collection{N} of a <:AbstractGrid or a subset based on an Int, Vector{Int} or String. The last option tries to call a nodeset of the <:AbstractGrid. Collection{N} refers to some indexable collection where each element corresponds to a Node.
Ferrite.getnnodes — FunctionReturns the number of nodes in the grid.
Ferrite.nnodes_per_cell — FunctionReturns the number of nodes of the i-th cell.
Ferrite.getcellset — Functiongetcellset(grid::AbstractGrid, setname::String)Returns all cells as cellid in a Set of a given setname.
Ferrite.getcellsets — Functiongetcellsets(grid::AbstractGrid)Returns all cellsets of the grid.
Ferrite.getnodeset — Functiongetnodeset(grid::AbstractGrid, setname::String)Returns all nodes as nodeid in a Set of a given setname.
Ferrite.getnodesets — Functiongetnodesets(grid::AbstractGrid)Returns all nodesets of the grid.
Ferrite.getfaceset — Functiongetfaceset(grid::AbstractGrid, setname::String)Returns all faces as FaceIndex in a Set of a given setname.
Ferrite.getfacesets — Functiongetfacesets(grid::AbstractGrid)Returns all facesets of the grid.
Ferrite.getedgeset — Functiongetedgeset(grid::AbstractGrid, setname::String)Returns all edges as EdgeIndex in a Set of a given setname.
Ferrite.getedgesets — Functiongetedgesets(grid::AbstractGrid)Returns all edge sets of the grid.
Ferrite.getvertexset — Functiongetedgeset(grid::AbstractGrid, setname::String)Returns all vertices as VertexIndex in a Set of a given setname.
Ferrite.getvertexsets — Functiongetvertexsets(grid::AbstractGrid)Returns all vertex sets of the grid.
Ferrite.compute_vertex_values — Functionfunction compute_vertex_values(grid::AbstractGrid, f::Function)
function compute_vertex_values(grid::AbstractGrid, v::Vector{Int}, f::Function)
function compute_vertex_values(grid::AbstractGrid, set::String, f::Function)Given a grid and some function f, compute_vertex_values computes all nodal values, i.e. values at the nodes, of the function f. The function implements two dispatches, where only a subset of the grid's node is used.
compute_vertex_values(grid, x -> sin(x[1]) + cos([2]))
compute_vertex_values(grid, [9, 6, 3], x -> sin(x[1]) + cos([2])) #compute function values at nodes with id 9,6,3
compute_vertex_values(grid, "right", x -> sin(x[1]) + cos([2])) #compute function values at nodes belonging to nodeset rightFerrite.transform! — Functiontransform!(grid::Abstractgrid, f::Function)Transform all nodes of the grid based on some transformation function f.
Ferrite.getcoordinates — Functiongetcoordinates(grid::AbstractGrid, cell)Return a vector with the coordinates of the vertices of cell number cell.
Ferrite.getcoordinates! — Functiongetcoordinates!(x::Vector{Vec{dim,T}}, grid::AbstractGrid, cell::Int)
getcoordinates!(x::Vector{Vec{dim,T}}, grid::AbstractGrid, cell::AbstractCell)Fills the vector x with the coordinates of a cell defined by either its cellid or the cell object itself.
Topology
Ferrite.ExclusiveTopology — TypeExclusiveTopology(cells::Vector{C}) where C <: AbstractCellExclusiveTopology saves topological (connectivity) data of the grid. The constructor works with an AbstractCell vector for all cells that dispatch vertices, faces and in 3D edges as well as the utility functions face_npoints and edge_npoints. The struct saves the highest dimensional neighborhood, i.e. if something is connected by a face and an edge only the face neighborhood is saved. The lower dimensional neighborhood is recomputed, if needed.
Fields
vertex_to_cell::Dict{Int,Vector{Int}}: global vertex id to all cells containing the vertexcell_neighbor::Vector{EntityNeighborhood{CellIndex}}: cellid to all connected cellsface_neighbor::SparseMatrixCSC{EntityNeighborhood,Int}:face_neighbor[cellid,local_face_id]-> neighboring facevertex_neighbor::SparseMatrixCSC{EntityNeighborhood,Int}:vertex_neighbor[cellid,local_vertex_id]-> neighboring vertexedge_neighbor::SparseMatrixCSC{EntityNeighborhood,Int}:edge_neighbor[cellid_local_vertex_id]-> neighboring edgevertex_vertex_neighbor::Dict{Int,EntityNeighborhood{VertexIndex}}: global vertex id -> all connected vertices by edge or faceface_skeleton::Vector{FaceIndex}: list of unique faces in the grid
Ferrite.getneighborhood — Functiongetneighborhood(top::ExclusiveTopology, grid::AbstractGrid, cellidx::CellIndex, include_self=false)
getneighborhood(top::ExclusiveTopology, grid::AbstractGrid, faceidx::FaceIndex, include_self=false)
getneighborhood(top::ExclusiveTopology, grid::AbstractGrid, vertexidx::VertexIndex, include_self=false)
getneighborhood(top::ExclusiveTopology, grid::AbstractGrid, edgeidx::EdgeIndex, include_self=false)Returns all directly connected entities of the same type, i.e. calling the function with a VertexIndex will return a list of directly connected vertices (connected via face/edge). If include_self is true, the given *Index is included in the returned list.
This feature is highly experimental and very likely subjected to interface changes in the future.
Ferrite.faceskeleton — Functionfaceskeleton(grid) -> Vector{FaceIndex}Returns an iterateable face skeleton. The skeleton consists of FaceIndex that can be used to reinit FaceValues.
Grid Sets Utility
Ferrite.addcellset! — Functionaddcellset!(grid::AbstractGrid, name::String, cellid::Union{Set{Int}, Vector{Int}})
addcellset!(grid::AbstractGrid, name::String, f::function; all::Bool=true)Adds a cellset to the grid with key name. Cellsets are typically used to define subdomains of the problem, e.g. two materials in the computational domain. The MixedDofHandler can construct different fields which live not on the whole domain, but rather on a cellset. all=true implies that f(x) must return true for all nodal coordinates x in the cell if the cell should be added to the set, otherwise it suffices that f(x) returns true for one node.
addcellset!(grid, "left", Set((1,3))) #add cells with id 1 and 3 to cellset left
addcellset!(grid, "right", x -> norm(x[1]) < 2.0 ) #add cell to cellset right, if x[1] of each cell's node is smaller than 2.0Ferrite.addfaceset! — Functionaddfaceset!(grid::AbstractGrid, name::String, faceid::Union{Set{FaceIndex},Vector{FaceIndex}})
addfaceset!(grid::AbstractGrid, name::String, f::Function; all::Bool=true)Adds a faceset to the grid with key name. A faceset maps a String key to a Set of tuples corresponding to (global_cell_id, local_face_id). Facesets are used to initialize Dirichlet structs, that are needed to specify the boundary for the ConstraintHandler. all=true implies that f(x) must return true for all nodal coordinates x on the face if the face should be added to the set, otherwise it suffices that f(x) returns true for one node.
addfaceset!(grid, "right", Set(((2,2),(4,2))) #see grid manual example for reference
addfaceset!(grid, "clamped", x -> norm(x[1]) ≈ 0.0) #see incompressible elasticity example for referenceFerrite.addnodeset! — Functionaddnodeset!(grid::AbstractGrid, name::String, nodeid::Union{Vector{Int},Set{Int}})
addnodeset!(grid::AbstractGrid, name::String, f::Function)Adds a nodeset::Dict{String, Set{Int}} to the grid with key name. Has the same interface as addcellset. However, instead of mapping a cell id to the String key, a set of node ids is returned.
Multithreaded Assembly
Ferrite.create_coloring — Functioncreate_coloring(g::Grid, cellset=1:getncells(g); alg::ColoringAlgorithm)Create a coloring of the cells in grid g such that no neighboring cells have the same color. If only a subset of cells should be colored, the cells to color can be specified by cellset.
Returns a vector of vectors with cell indexes, e.g.:
ret = [
[1, 3, 5, 10, ...], # cells for color 1
[2, 4, 6, 12, ...], # cells for color 2
]Two different algorithms are available, specified with the alg keyword argument:
alg = ColoringAlgorithm.WorkStream(default): Three step algorithm from WorkStream , albeit with a greedy coloring in the second step. Generally results in more colors thanColoringAlgorithm.Greedy, however the cells are more equally distributed among the colors.alg = ColoringAlgorithm.Greedy: greedy algorithm that works well for structured quadrilateral grids such as e.g. quadrilateral grids fromgenerate_grid.
The resulting colors can be visualized using vtk_cell_data_colors.
In a previous version of Ferrite this function returned a dictionary mapping cell ID to color numbers as the first argument. If you need this mapping you can create it using the following construct:
colors = create_coloring(...)
cell_colormap = Dict{Int,Int}(
cellid => color for (color, cellids) in enumerate(final_colors) for cellid in cellids
)