ExtendableGrids.ElementInfoType
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.AT_NODESType
abstract type AT_NODES <: AssemblyType

causes interpolation at vertices of the grid (only for H1-conforming interpolations)

ExtendableGrids.BFaceCellsType
abstract type BFaceCells <: AbstractGridAdjacency

Adjacency describing cells per boundary or interior face

ExtendableGrids.BFaceEdgesType
abstract type BFaceEdges <: AbstractGridAdjacency

Adjacency describing edges per boundary or interior face

ExtendableGrids.BFaceNodesType
abstract type BFaceNodes <: AbstractGridAdjacency

Adjacency describing nodes per grid boundary face

ExtendableGrids.BFaceNormalsType
abstract type BFaceNormals <: AbstractGridComponent

Adjacency describing outer normals to boundary faces

ExtendableGrids.BFaceParentsType
abstract type BFaceParents <: AbstractGridIntegerArray1D

Grid component key type for storing parent bfaces

ExtendableGrids.BinnedPointListType
mutable struct BinnedPointList{T}

Binned point list structure allowing for fast check for already existing points.

This provides better performance for indendifying already inserted points than the naive linear search.

OTOH the implementation is still quite naive - it dynamically maintains a cuboid binning region with a fixed number of bins.

Probably tree based adaptive methods (a la octree) will be more efficient, however they will be harder to implement.

In an ideal world, we would maintain a dynamic Delaunay triangulation, which at once could be the starting point of mesh generation which will follow here anyway.

  • dim::Int32: Space dimension
  • tol::Any: Point distance tolerance. Points closer than tol (in Euclidean distance) will be identified, i.e. are collapsed to the first inserted.
  • binning_region_min::Vector: " The union of all bins is the binning region - a cuboid given by two of its corners. It is calculated dynamically depending on the inserted points.
  • binning_region_max::Vector

  • binning_region_increase_factor::Any: Increase factor of binning region (with respect to the cuboid defined by the coordinates of the binned points)

  • points::ElasticArrays.ElasticArray{T, 2, M, V} where {T, M, V<:DenseVector{T}}: The actual point list
  • bins::Array{Vector{Int32}}: The bins are vectors of indices of points in the point list We store them in a dim-dimensional array of length "numberofdirectional_bins^dim"
  • number_of_directional_bins::Int32: Number of bins in each space dimension
  • unbinned::Vector{Int32}: Some points will fall outside of the binning region. We collect them in vector of ubinned point indices
  • num_allowed_unbinned_points::Int32: Number of unbinned points tolerated without rebinning
  • max_unbinned_ratio::Any: Maximum ratio of unbinned points in point list
  • current_bin::Vector{Int32}: Storage of current point bin
ExtendableGrids.BinnedPointListMethod
 BinnedPointList(::Type{T}, dim;
                 tol = 1.0e-12,
                 number_of_directional_bins = 10,
                 binning_region_increase_factor = 0.01,
                 num_allowed_unbinned_points = 5,
                 max_unbinned_ratio = 0.05) where {T}

Create and initialize binned point list

ExtendableGrids.Cartesian3DType
abstract type Cartesian3D <: AbstractCoordinateSystem

2D cartesion coordinate system (unknowns x,y,z)

ExtendableGrids.CellParentsType
abstract type CellParents <: AbstractGridIntegerArray1D

Grid component key type for storing parent cells

ExtendableGrids.FaceParentsType
abstract type FaceParents <: AbstractGridIntegerArray1D

Grid component key type for storing parent faces (only for SubGrid relation when FaceNodes is instantiated)

ExtendableGrids.NodeParentsType
abstract type NodeParents <: AbstractGridIntegerArray1D

Grid component key type for storing node parents (=ids of nodes in ParentGrid) in an array

ExtendableGrids.NodePermutationType
abstract type NodePermutation <: AbstractGridIntegerArray1D

Key type describing the permutation of the nodes of a partitioned grid with respect to the unpartitioned origin.

If pgrid is the partitioned grid and grid is the unpartitioned origin, then

pgrid[Coordinates][:,pgrid[NodePermutation]]==grid[Coordinates]

ExtendableGrids.ON_BEDGESType
abstract type ON_BEDGES <: AssemblyType

causes assembly/interpolation on boundary edges of the grid (only in 3D)

ExtendableGrids.ON_BFACESType
abstract type ON_BFACES <: AssemblyType

causes assembly/interpolation on boundary faces of the grid

ExtendableGrids.ON_CELLSType
abstract type ON_CELLS <: AssemblyType

causes assembly/interpolation on cells of the grid

ExtendableGrids.ON_EDGESType
abstract type ON_EDGES <: AssemblyType

causes assembly/interpolation on edges of the grid (only in 3D)

ExtendableGrids.ON_FACESType
abstract type ON_FACES <: AssemblyType

causes assembly/interpolation on faces of the grid

ExtendableGrids.ON_IFACESType
abstract type ON_IFACES <: ON_FACES

causes assembly/interpolation on interior faces of the grid

ExtendableGrids.PColorPartitionsType
abstract type PColorPartitions <: AbstractGridIntegerArray1D

Key type describing colors of partitions. These correspond to a coloring of the neigborhood graphs of partitions such that operations (e.g. FEM assembly) on partitions of a given color can be performed in parallel.

grid[PColorPartitions] returns an integer vector describing the partition colors ("pcolors") of a grid. Let p=grid[PColorPartitions]. Then all partitions with numbers i ∈ p[c]:p[c+1]-1 have "color" c. See also pcolors.

ExtendableGrids.ParentGridType
abstract type ParentGrid <: AbstractGridComponent

Grid component key type for storing parent grid

ExtendableGrids.PartitionBFacesType
abstract type PartitionBFaces <: AbstractGridIntegerArray1D

Key type describing the bondary faces of a given partition.

grid[PartitionBFaces] returns an integer vector describing the boundary faces of a partition given by its number. Let pc=grid[PartitionCells]. Then all cells with index i ∈ pc[p]:pc[p+1]-1 belong to partition p.

ExtendableGrids.PartitionCellsType
abstract type PartitionCells <: AbstractGridIntegerArray1D

Key type describing the cells of a given partition.

grid[PartitionCells] returns an integer vector describing the cells of a partition given by its number. Let pc=grid[PartitionCells]. Then all cells with index i ∈ pc[p]:pc[p+1]-1 belong to partition p.

ExtendableGrids.PartitionEdgesType
abstract type PartitionEdges <: AbstractGridIntegerArray1D

Key type describing the edges of a given partition.

grid[PartitionEdges] returns an integer vector describing the edges of a partition given by its number. Let pe=grid[PartitionEdges]. Then all edges with index i ∈ pe[p]:pe[p+1]-1 belong to partition p.

ExtendableGrids.PartitionNodesType
abstract type PartitionNodes <: AbstractGridIntegerArray1D

Key type describing the nodes of a given partition.

grid[PartitionNodes] returns an integer vector describing the nodes of a partition given by its number. Let pn=grid[PartitionNodes]. Then all nodes with index i ∈ pn[p]:pn[p+1]-1 belong to partition p.

ExtendableGrids.PlainMetisPartitioningType
struct PlainMetisPartitioning <: AbstractPartitioningAlgorithm

Subdivide grid into npart partitions using Metis.partition and color the resulting partition neigborhood graph. This requires to import Metis.jl in order to trigger the corresponding extension.

This algorithm allows to control the overall number of partitions. The number of partitions per color comes from the subsequent partition graph coloring and in the moment cannot be controled.

Parameters:

  • npart::Int64: Number of partitions (default: 20)
ExtendableGrids.Polar1DType
abstract type Polar1D <: AbstractCoordinateSystem

1D polar coordinate system (unknown r)

ExtendableGrids.Polar2DType
abstract type Polar2D <: AbstractCoordinateSystem

2D polar coordinate system (unknowns r,ϕ)

ExtendableGrids.RecursiveMetisPartitioningType
struct RecursiveMetisPartitioning <: AbstractPartitioningAlgorithm

Subdivide grid into npart partitions using Metis.partition and calculate cell separators from this partitioning. The initial partitions get color 1, and the separator gets color 2. This is continued recursively with partitioning of the separator.

This algorithm allows to control the number of partitions in color 1 which coorespond to the bulk of the work.

Parameters:

  • npart::Int64: Number of color 1 partitions (default: 4)

  • maxdepth::Int64: Recursion depth (default: 1)

  • separatorwidth::Int64: Separator width (default: 2)

ExtendableGrids.RefinedGridType
abstract type RefinedGrid <: ParentGridRelation

Grid component key type for indicating that grid is a refinement of the parentgrid

ExtendableGrids.Spherical3DType
abstract type Spherical3D <: AbstractCoordinateSystem

3D spheriacal coordinate system (unknowns r,ϕ,θ)

ExtendableGrids.SubGridType
abstract type SubGrid{support} <: ParentGridRelation

Grid component key type for indicating that grid is a subgrid of the parentgrid

ExtendableGrids.SubgridVectorViewType
struct SubgridVectorView{Tv, Ti} <: AbstractArray{Tv, 1}

Vector view on subgrid

  • sysarray::AbstractVector

  • node_in_parent::Vector

ExtendableGrids.TokenStreamType
mutable struct TokenStream

Tokenstream allows to read tokenized data from file without keeping the file ocntent in memory.

  • input::IOStream: Input stream
  • tokens::Vector{SubString{String}}: Array of current tokens kept in memory.
  • itoken::Int64: Position of actual token in tokens array
  • lineno::Int64: Line number in IOStream
  • comment::Char: Comment character
  • dlm::Function: Function telling if given character is a delimiter.
ExtendableGrids.TokenStreamMethod
TokenStream(input::IOStream; comment, dlm) -> TokenStream

Create Tokenstream with IOStream argument.

ExtendableGrids.TokenStreamMethod
TokenStream(filename::String; comment, dlm) -> TokenStream

Create Tokenstream with file name argument.

ExtendableGrids.TrivialPartitioningType
struct TrivialPartitioning <: AbstractPartitioningAlgorithm

Trivial partitioning: all grid cells belong to single partition number 1.

ExtendableGrids.UnexpectedTokenErrorType
struct UnexpectedTokenError <: Exception

Error thrown when the token expected in expect! is not there.

  • found::String

  • expected::String

  • lineno::Int64

ExtendableGrids.VariableTargetAdjacencyType
struct VariableTargetAdjacency{T}

Adjacency struct. Essentially, this is the sparsity pattern of a matrix whose nonzero elements all have the same value in the CSC format.

ExtendableGrids.VariableTargetAdjacencyMethod
VariableTargetAdjacency(
    m::SparseArrays.SparseMatrixCSC{Tv<:Integer, Ti<:Integer}
) -> VariableTargetAdjacency{Ti} where Ti<:Integer

Create variable target adjacency from adjacency matrix

AbstractTrees.childrenMethod
children(T::Type) -> Union{Vector{Type}, Vector{Any}}

Define children for types.

Base.:==Method
==(a, b)
Comparison of two adjacencies
Base.:==Method
==(a, b)
Comparison of two adjacencies
Base.append!Method
append!(adj::SerialVariableTargetAdjacency, len) -> Vector

Append a column to adjacency.

Base.append!Method
append!(adj::VariableTargetAdjacency, column) -> Vector

Append a column to adjacency.

Base.delete!Method
delete!(
    grid::ExtendableGrid,
    T::Type{<:AbstractGridComponent}
) -> Dict{Type{<:AbstractGridComponent}, Any}

Remove grid component

Base.eofMethod
eof(tks::TokenStream) -> Bool

Check if all tokens have been consumed.

Base.get!Method
get!(
    grid::ExtendableGrid,
    T::Type{<:AbstractGridComponent}
) -> Any

To be called by getindex. This triggers lazy creation of non-existing gridcomponents

Base.getindexMethod
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.

Base.getindexMethod
getindex(
    aview::ExtendableGrids.SubgridVectorView,
    inode::Integer
) -> Any

Accessor method for subgrid vector view.

Base.getindexMethod
getindex(
    adj::SerialVariableTargetAdjacency,
    i,
    isource
) -> Any

Access adjacency as if it is a 2D Array

Base.getindexMethod
getindex(adj::VariableTargetAdjacency, i, isource) -> Any

Access adjacency as if it is a 2D Array

Base.getindexMethod
getindex(v::VectorOfConstants, i) -> Any

Access

Base.haskeyMethod
haskey(g::ExtendableGrid, k) -> Bool

Check if key is in grid

Base.insert!Method
 Base.insert!(binnedpointlist,p)

If another point with distance less the tol from p is in pointlist, return its index. Otherwise, insert point into pointlist. p may be a vector or a tuple.

Base.insert!Method
 Base.insert!(binnedpointlist,x,y,z)

Insert 3D point via coordinates.

Base.insert!Method
 Base.insert!(binnedpointlist,x)

Insert 1D point via coordinate.

Base.iterateMethod
iterate(
    v::VectorOfConstants,
    state
) -> Union{Nothing, Tuple{Any, Any}}

Iterator

Base.iterateMethod
iterate(v::VectorOfConstants) -> Tuple{Any, Int64}

Iterator

Base.keysMethod
keys(
    g::ExtendableGrid
) -> Base.KeySet{Type{<:AbstractGridComponent}, Dict{Type{<:AbstractGridComponent}, Any}}

Keys in grid

Base.lengthMethod
length(v::VectorOfConstants) -> Any

Length

Base.mapMethod
map(f::Function, grid::ExtendableGrid) -> Any

Map a function onto node coordinates of grid

Base.setindex!Method
setindex!(
    grid::ExtendableGrid,
    v,
    T::Type{<:AbstractGridComponent}
) -> Any

Set new grid component

Base.setindex!Method
setindex!(
    aview::ExtendableGrids.SubgridVectorView,
    v,
    inode::Integer
) -> ExtendableGrids.SubgridVectorView

Accessor method for subgrid vector view.

Base.showMethod
show(io::IO, adj::SerialVariableTargetAdjacency)

Show adjacency (in trasposed form; preliminary)

Base.showMethod
show(io::IO, adj::VariableTargetAdjacency)

Show adjacency (in trasposed form; preliminary)

Base.sizeMethod
size(a::ExtendableGrids.SubgridVectorView) -> Tuple{Int64}

Return size of vector view.

Base.sizeMethod
size(v::VectorOfConstants) -> Tuple{Any}

Size

Base.uniqueMethod
unique(v::VectorOfConstants) -> Vector

Shortcut for unique

Base.viewMethod
view(
    a::AbstractVector,
    subgrid::ExtendableGrid
) -> ExtendableGrids.SubgridVectorView

Create a view of the vector on a subgrid.

ExtendableGrids.RGB_refineMethod
RGB_refine(
    source_grid::ExtendableGrid{T, K},
    facemarkers::Vector{Bool};
    store_parents
) -> ExtendableGrid

generates a new ExtendableGrid by red-green-blue mesh refinement of triangular meshes, see e.g.

Carstensen, C. –An Adaptive Mesh-Refining Algorithm Allowing for an H^1 Stable L^2 Projection onto Courant Finite Element Spaces– Constr Approx 20, 549–564 (2004). https://doi.org/10.1007/s00365-003-0550-5

The bool array facemarkers determines which faces should be bisected. Note, that a closuring is performed such that the first face in every triangle with a marked face is also refined.

ExtendableGrids._findpointMethod
_findpoint(binnedpointlist, index, p)

Find point in index list (by linear search) Return its index, or zero if not found

ExtendableGrids._rebin_all_points!Method
_rebin_all_points!(bpl)

Re-calculate binning if there are too many unbinned points This amounts to two steps:

  • Enlarge binning area in order to include all points
  • Re-calculate all point bins
ExtendableGrids.asparseMethod
asparse(a::Matrix) -> SparseArrays.SparseMatrixCSC{Int64}

Create sparse incidence matrix from adjacency

ExtendableGrids.asparseMethod
asparse(
    a::VariableTargetAdjacency
) -> SparseArrays.SparseMatrixCSC{Int64}

Create sparse incidence matrix from adjacency

ExtendableGrids.assemble_bfacesMethod
function assemble_bfaces(simplices, dim, nn, Ti)

Assemble the BoundaryFaces corresponding to the simplices passed. In this function, the faces are encoded for performance reasons. If a large grid with many nodes is used, Ti has to be chosen accordingly (e.g. Int128), or encode=false has to be passed to seal!. simplices is a $(dim+1) x 'number cells'$ matrix and nn is the total number of nodes. We can not guarantee, that the orientation of the BoundaryFaces is correct.

ExtendableGrids.assemble_bfaces_directMethod
function assemble_bfaces_direct(simplices, dim, Ti)

Assemble the BoundaryFaces corresponding to the simplices passed. In this function, the faces are not encoded. This may make sense for grids with many nodes. For smaller grids it can lead to performance losses. simplices is a $(dim+1) x 'number cells'$ matrix and nn is the total number of nodes. We can not guarantee, that the orientation of the BoundaryFaces is correct.

ExtendableGrids.barycentric_refineMethod
barycentric_refine(
    source_grid::ExtendableGrid{T, K};
    store_parents
) -> ExtendableGrid

generates a new ExtendableGrid by barycentric refinement of each cell in the source grid

barycentric refinement is available for these ElementGeometries

  • Quadrilateral2D (first split into Triangle2D)
  • Triangle2D
ExtendableGrids.bedgemask!Method
bedgemask!(
    grid::ExtendableGrid,
    xa,
    xb,
    ireg::Int64;
    tol
) -> ExtendableGrid

Edit region numbers of grid boundary edges via line mask. This only works for 3D grids.

ExtendableGrids.bfacemask!Method
bfacemask!(grid::ExtendableGrid,
                maskmin,
                maskmax,
                ireg;
                allow_new=true,
                tol=1.0e-10)

Edit region numbers of grid boundary facets via rectangular mask. If allow_new is true (default), new facets are added.

ireg may be an integer or a function ireg(current_region).

A zero region number removes boundary faces.

Examples: Rectangle-with-multiple-regions

ExtendableGrids.check_partitioningMethod
check_partitioning(grid; 
                   verbose=true, 
                   cellpartonly=false)

Check correctness of cell partitioning, necessesary for parallel assembly:

  • Check if every node belongs to one of the cell partitions
  • Check if no node belongs to two cell partitions of the same color at once

If cellpartonly==false check correctness of node partitioning necessary for parallel sparse matrix multiplication and ILU preconditioning

  • Check if no node belongs to two node partitions of the same color at once
  • Check if no node is a neighbor of nodes from two node partitions of the same color
ExtendableGrids.decodeMethod
function decode(y::Integer, nn::Integer, dim::Integer)

Decode y to the vector x. x has the length dim. The en/-decoding is similar to using the base-nn number system. For details of the encoding, see the documentation of the function encode.

ExtendableGrids.dopartitionMethod
dopartition(grid, alg)

(Internal utility function) Core function for partitioning grid cells which dispatches over partitioning algorithms. Partitioning extensions should add methods to this function.

ExtendableGrids.encodeMethod
function encode(x::Vector, nn::Integer)

Encode th vector x into an Int y. The en/-decoding is similar to using the base-nn number system. Example: $[x₁, x₂, x₃] → (x₁-1) + (x₂-1)*nn + (x₃-1)*nn²$``

ExtendableGrids.expecttokenMethod
expecttoken(tks::TokenStream, expected::String) -> Bool

Expect keyword token.

If token is missing, an UnexpectedTokenError is thrown If the token has been found, reading will continue at the position after the token found.

ExtendableGrids.faces_of_ndim_simplexMethod
function faces_of_ndim_simplex(x::Vector, dim::Integer, nn::Integer)

Return all faces of a n-dim simplex. The orientation is not guaranteed to be right. x contains the nodes of the simplex. nn is the total number of nodes. The faces (=the nodes contained in the face), are encoded to Integers (of nn's type).

ExtendableGrids.faces_of_ndim_simplex_directMethod
function faces_of_ndim_simplex(x::Vector, dim::Integer, nn::Integer)

Return all faces of a n-dim simplex. The orientation is not guaranteed to be right. x contains the nodes of the simplex. nn is the total number of nodes. The faces (=the nodes contained in the face), are not encoded to Integers.

ExtendableGrids.findpointMethod
findpoint(binnedpointlist, p)

Find point in binned point list. Return its index in the point list if found, otherwise return 0.

ExtendableGrids.gFindBruteForce!Method
icellfound=gFindBruteForce!(xref,cellfinder,p; icellstart=1,eps=1.0e-14)

Find cell containing point p starting with cell number icellstart.

Returns cell number if found, zero otherwise. Upon return, xref contains the barycentric coordinates of the point in the sequence dim+1, 1...dim

Warning

Currently implemented for simplex grids only.

ExtendableGrids.gFindLocal!Method
icellfound=GFindLocal!(xref,cellfinder,p; icellstart=1,eps=1.0e-14, trybrute=true)

Find cell containing point p starting with cell number icellstart.

Returns cell number if found, zero otherwise. If trybrute==true try gFindBruteForce! before giving up. Upon return, xref contains the barycentric coordinates of the point in the sequence dim+1, 1...dim

Warning

Currently implemented for simplex grids only.

ExtendableGrids.geomspaceMethod
geomspace(a, b, ha, hb; tol, maxiterations) -> Any

(Try to) create a subdivision of interval (a,b) stored in the returned array X such that

  • X[1]==a, X[end]==b
  • (X[2]-X[1])<=ha+tol*(b-a)
  • (X[end]-X[end-1])<=hb+tol*(b-a)
  • There is a number q such that X[i+1]-X[i] == q*(X[i]-X[i-1])
  • X is the array with the minimal possible number of points with the above property

Caveat: the algorithm behind this is tested for many cases but unproven.

Returns an Array containing the points of the subdivision.

ExtendableGrids.gettokenMethod
gettoken(
    tks::TokenStream
) -> Union{Nothing, SubString{String}}

Get next token from tokenstream.

ExtendableGrids.glueMethod
c=glue(a,b)

Glue together two vectors a and b resulting in a vector c. They last element of a shall be equal (up to tol) to the first element of b. The result fulfills length(c)=length(a)+length(b)-1

ExtendableGrids.glueMethod
glue(g1,g2;
     g1regions=1:num_bfaceregions(g1),
     g2regions=1:num_bfaceregions(g2),
     interface=0,
     warnonly = false,
     tol=1.0e-10,
     naive=false)

Merge two grids along their common boundary facets.

  • g1: First grid to be merged
  • g2: Second grid to be merged
  • g1regions: boundary regions to be used from grid1. Default: all.
  • g2regions: boundary regions to be used from grid2. Default: all.
  • interface: if nonzero, create interface region in new grid, otherwise, ignore
  • strict: Assume all bfaces form specfied regions shall be matched, throw error on failure
  • tol: Distance below which two points are seen as identical. Default: 1.0e-10
  • naive: use naive quadratic complexity matching (for checking backward compatibility). Default: false

Deprecated:

  • breg: old notation for interface
ExtendableGrids.grid_triangleMethod
grid_triangle(coords::AbstractArray{T,2}) where {T}

Generates a single triangle with the given coordinates, that should be a 2 x 3 array with the coordinates of the three vertices, e.g. coords = [0.0 0.0; 1.0 0.0; 0.0 1.0]'.

ExtendableGrids.grid_unitcubeMethod
grid_unitcube(EG::Type{<:Hexahedron3D}; scale = [1,1,1], shift = [0,0,0])

Unit cube as one cell with six boundary regions (bottom, front, right, back, left, top)

ExtendableGrids.grid_unitcubeMethod
grid_unitcube(::Type{Tetrahedron3D}; scale = [1,1,1], shift = [0,0,0])

Unit cube as six tets with six boundary regions (bottom, front, right, back, left, top)

ExtendableGrids.grid_unitsquareMethod
grid_unitsquare(EG::Type{<:Quadrilateral2D}; scale = [1,1], shift = [0,0])

Unit square as one cell with four boundary regions (bottom, right, top, left)

ExtendableGrids.grid_unitsquareMethod
grid_unitsquare(::Type{<:Triangle2D}; scale = [1,1], shift = [0,0])

Unit square as two triangles with four boundary regions (bottom, right, top, left)

ExtendableGrids.induce_edge_partitioning!Method
induce_edge_partitioning!(grid; trivial)

(internal) Induce edge partitioning from cell partitioning of grid. The algorithm assumes that nodes get the partition number from the partition numbers of the cells having this node in common. If these are differnt, the highest number is taken.

This method triggers creation of rather complex edge information and should be called only if this information is really necessary.

ExtendableGrids.induce_node_partitioning!Method
induce_node_partitioning!(
    grid,
    cn,
    nc;
    trivial,
    keep_nodepermutation
)

(internal) Induce node partitioning from cell partitioning of grid. The algorithm assumes that nodes get the partition number from the partition numbers of the cells having this node in common. If these are differnt, the highest number is taken.

Node partitioning should support parallel matrix-vector products with SparseMatrixCSC. The current algorithm assumes that nodes get the partition number from the partition numbers of the cells having this node in common. If these are differnt, the highest number is taken.

Simply inducing node partition numbers from cell partition numbers does not always fulfill the condition that there is no node which is neigbour of nodes from two different partition with the same color.

This situation is detected and corrected by joining respective critical partitions, sacrificing a bit of parallel efficiency for correctness.

ExtendableGrids.instantiateMethod
instantiate(grid::ExtendableGrid, ::Type{PColorPartitions})

If not given otherwise, instantiate partition data with trivial partitioning.

ExtendableGrids.instantiateMethod
instantiate(grid::ExtendableGrid, ::Type{PartitionBFaces})

If not given otherwise, instantiate partition data with trivial partitioning.

ExtendableGrids.instantiateMethod
instantiate(grid::ExtendableGrid, ::Type{PartitionCells})

If not given otherwise, instantiate partition data with trivial partitioning.

ExtendableGrids.instantiateMethod
instantiate(grid::ExtendableGrid, ::Type{PartitionNodes})

If not given otherwise, instantiate partition data with trivial partitioning.

ExtendableGrids.instantiateMethod
instantiate(grid::ExtendableGrid, ::Type{PartitionEdges})

If not given otherwise, instantiate partition data with trivial partitioning.

ExtendableGrids.interpolateMethod
u_to=interpolate(grid_to, u_from, grid_from;eps=1.0e-14,trybrute=true)

Piecewise linear interpolation of function u_from on grid grid_from to grid_to. Works for matrices with second dimension corresponding to grid nodes and for vectors.

Warning

May be slow on non-convex domains. If trybrute==false it may even fail.

Warning

Currently implemented for simplex grids only.

ExtendableGrids.isconsistentMethod
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.linspaceMethod
linspace(a, b, n) -> Any

Resurrect linspace despite https://github.com/JuliaLang/julia/pull/25896#issuecomment-363769368

ExtendableGrids.makevarMethod
makevar(a::Array{T, 2}) -> VariableTargetAdjacency

Turn fixed target adjacency into variable target adjacency

ExtendableGrids.naiveinsert!Method
naiveinsert(binnedpointlist, p)

Insert via linear search, without any binning. Just for being able to check of all of the above was worth the effort...

ExtendableGrids.num_partitions_per_colorMethod
num_partitions_per_color(grid)

Return a vector containing the number of partitions for each of the colors of the grid partitioning. These define the maximum number of parallel threads for each color.

ExtendableGrids.num_targetsMethod
num_targets(
    adj::SerialVariableTargetAdjacency,
    isource
) -> Any

Number of targets for given source

ExtendableGrids.partgraphMethod
partgraph(cellpartitions, ncellpartitions, cellcelladj)

(internal) Create neigbourhood graph for given partitioning.

ExtendableGrids.partitionMethod
     partition(grid::ExtendableGrid,
               alg::AbstractPartitioningAlgorithm;
               nodes = false,
               keep_nodepermutation = false,
               edges = false )

Partition cells of grid according to alg, such that the neigborhood graph of partitions is colored in such a way, that all partitions with a given color can be worked on in parallel. Cells are renumbered such that cell numbers for a given partition are numbered contiguously.

Return the resulting grid.

Useful for parallel FEM assembly and cellwise FVM assembly.

Keyword arguments:

  • nodes: if true, induce node partitioning from cell partitioning. Used for node/edgewise FVM assembly. In addition the resulting partitioning supports parallel matrix-vector products with SparseMatrixCSC. Nodes are renumbered compared to the original grid.
  • keep_nodepermutation: if true, keep the node permutation with respect to the original grid in grid[NodePermutation].
  • edges: if true, induce partitioning of edges from cell partitioning. Used for node/edgewise FVM assembly. This step creates a number of relatively expensive additional adjacencies.

Access:

A parallel loop over grid cells thus looks like

for color in pcolors(grid)
    @threads for part in pcolor_partitions(grid, color)
                for cell in partition_cells(grid, part)
                 ...
                end
             end
end

Without a call to partition, all these functions return trivial data such that the above sample code stays valid.

Note

partition must be called before obtaining any other adjacencies of a grid.

Currently, partitioning does not cover the boundary, boundary cells belong to one big trivial partition.

ExtendableGrids.prepare_edges!Method
prepare_edges!(grid)

Prepare edge adjacencies (celledges, edgecells, edgenodes)

Currently depends on ExtendableSparse, we may want to remove this adjacency.

ExtendableGrids.rect!Method
rect!(grid,maskmin,maskmax; 
      region=1, 
      bregion=1, 
      bregions=nothing, 
      tol=1.0e-10)

Place a rectangle into a rectangular grid. It places a cellmask according to maskmin and maskmax, and introduces boundary faces via `bfacesmask! at all sides of the mask area. It is checked that the coordinate values in the mask match (with tolerance) corresponding directional coordinates of the grid.

If bregions is given it is assumed to be a vector corresponding to the number of sides, im the sequence w,e in 1D. s,e,n,w in 2D and s,e,n,w,b,t in 3D.

bregion or elements of bregions can be numbers or functions ireg(current_region).

Examples: Subgrid-from-rectangle, Rect2d-with-bregion-function, Cross3d

ExtendableGrids.reference_domainFunction
    reference_domain(EG::Type{<:AbstractElementGeometry}, T::Type{<:Real} = Float64; scale = [1,1,1], shift = [0,0,0]) -> ExtendableGrid{T,Int32}

Generates an ExtendableGrid{T,Int32} for the reference domain of the specified Element Geometry. With scale and shift the coordinates can be manipulated.

ExtendableGrids.reorder_cellsMethod
reorder_cells(
    grid,
    cellpartitions,
    ncellpartitions,
    colpart
)

(Internal utility function) Create cell permutation such that all cells belonging to one partition are numbered contiguously, return grid with reordered cells.

ExtendableGrids.ringsectorMethod
ringsector(rad,ang; eltype=Triangle2D)

Sector of ring or full ring (if ang[begin]-ang[end]≈2π)

ExtendableGrids.seal!Method
function seal!(grid::ExtendableGrid; bfaceregions=[], encode=true, Ti=Int64)

Take an (simplex-) ExtendableGrid and compute and add the BoundaryFaces. A so called incomplete ExtendableGrid can e.g. be read from an msh file using the Gmsh.jl-extension of the ExtendableGrids package and the function $simplexgrid_from_gmsh(filename::String; incomplete=true)$. If a non empty vector is passed as bfaceregions, this vector is used for the 'BFaceRegions'. If bfaceregions is empty, all BoundaryFaces get the region number 1.

For performance reasons, the faces (=the nodes contained in the face) can be encoded (see the function $encode(x::Vector, nn::Integer)$) to Integers encoding_type. To do this, encode=true is used. But for each encoding_type there is a limit on the number of nodes:

- For Int64  and a 2d grid: 3*10^9 nodes
- For Int64  and a 3d grid: 2*10^6 nodes
- For Int128 and a 2d grid: 1.3*10^19 nodes
- For Int128 and a 3d grid: 5.5*10^12 nodes

If encode=false is passed, there is no limit (besides the MaxValue of the Integer type used).

ExtendableGrids.seemingly_equalMethod
seemingly_equal(grid1, grid2; sort=false, confidence=:full

Recursively check seeming equality of two grids. Seemingly means that long arrays are only compared via random samples.

Keyword args:

  • sort: if true, sort grid points
  • confidence: 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.simplexgridMethod
simplexgrid(X,Y,Z; bregions=[1,2,3,4,5,6],cellregion=1)

Constructor for 3D grid from coordinate arrays. Boundary region numbers:

locationnumber
south1
east2
north3
west4
bottom5
top6

The keyword arguments allow to overwrite the default region numbers.

ExtendableGrids.simplexgridMethod
simplexgrid(X,Y; bregions=[1,2,3,4],cellregion=1)

Constructor for 2D grid from coordinate arrays.

Boundary region numbers count counterclockwise:

locationnumber
south1
east2
north3
west4

The keyword arguments allow to overwrite the default region numbers.

ExtendableGrids.simplexgridMethod
simplexgrid(X; bregions=[1,2],cellregion=1)

Constructor for 1D grid.

Construct 1D grid from an array of node cordinates. It creates two boundary regions with index 1 at the left end and index 2 at the right end by default.

The keyword arguments allow to overwrite the default region numbers.

Primal grid holding unknowns: marked by o, dual grid marking control volumes: marked by |.

 o-----o-----o-----o-----o-----o-----o-----o-----o
 |--|-----|-----|-----|-----|-----|-----|-----|--|
ExtendableGrids.simplexgridMethod
simplexgrid(grid2d::ExtendableGrid, coordZ; bot_offset=0,cell_offset=0,top_offset=0, bface_offset=0)

Create tensor product of 2D grid and 1D coordinate array.

Cellregions and outer facet regions are taken over from 2D grid and added to cell_offset and bface_offset, respectively. Top an bottom facet regions are detected from the cell regions and added to bot_offset resp. top_offset.

ExtendableGrids.simplexgridMethod
simplexgrid(
    file::String;
    format,
    kwargs...
) -> ExtendableGrid{Float64, Int32}

Read grid from file. Supported formats:

  • "*.sg": pdelib sg files. Format versions:
    • format=v"2.0": long version with some unnecessary data
    • format=v"2.1": shortened version only with cells, cellnodes, cellregions, bfacenodes, bfaceregions
    • format=v"2.2": like 2.1, but additional info on cell and node partitioning. Edge partitioning is not stored in the file and may be re-established by induce_edge_partitioning!.
  • "*.geo": gmsh geometry description (requires using Gmsh)
  • "*.msh": gmsh mesh (requires using Gmsh)
ExtendableGrids.simplexgridMethod
function simplexgrid(coord::Array{Tc,2},
                     cellnodes::Array{Ti,2},
                     cellregions,
                     bfacenodes,
                     bfaceregions
                     ) where {Tc,Ti}

Create d-dimensional simplex grid from five arrays.

  •    coord: d ``\times`` n_points matrix of coordinates
  • cellnodes: d+1 $\times$ n_tri matrix of triangle - point incidence
  • cellregions: n_tri vector of cell region markers
  • bfacenodes: d $\times$ n_bf matrix of boundary facet - point incidences
  • bfaceregions: n_bf vector of boundary facet region markers

Coordinate type Tc index type Ti are detected from the first two parameters. cellregions, bfaceregions, bfacenodes are converted to have the same element type as cellnodes.

ExtendableGrids.simplexgridMethod
function simplexgrid(coord::Array{Tc,2},
                     cellnodes::Array{Ti,2},
                     cellregions,
                     bfacenodes,
                     bfaceregions,
                     bedgenodes,
                     bedgeregions
                     ) where {Tc,Ti}

Create simplex grid from coordinates, cell-nodes-adjancency, cell-region-numbers, boundary-face-nodes adjacency, boundary-face-region-numbers, boundary-edge-nodes, and boundary-edge-region-numbers arrays.

The index type Ti is detected from cellnodes, all other arrays besides coord are converted to this index type.

ExtendableGrids.split_grid_intoMethod
split_grid_into(
    source_grid::ExtendableGrid{T, K},
    targetgeometry::Type{<:AbstractElementGeometry};
    store_parents
) -> ExtendableGrid

generates a new ExtendableGrid by splitting each cell into subcells of the specified targetgeometry

split rules exist for

  • Quadrilateral2D into Triangle2D
  • Hexahedron3D into Tetrahedron3D
ExtendableGrids.subgridMethod
subgrid(parent,                                                             
        subregions::AbstractArray;                                          
        transform::T=function(a,b) @views a.=b[1:length(a)] end,                                      
        boundary=false,                                                     
        coordinatesystem=codim1_coordinatesystem(parent[CoordinateSystem]), 
        project=true) where T

Create subgrid from list of regions.

  • parent: parent grid
  • subregions: Array of subregions which define the subgrid
  • 'support': support of subgrid, default is ONCELLS but can be also ONFACES or ON_BFACES to create codimension 1 subgrid from face/bfaces region
  • boundary: if true, create codimension 1 subgrid from boundary regions (same as support = ON_BFACES)
  • transform (kw parameter): transformation function between grid and subgrid coordinates acting on one point.
  • coordinatesystem: if boundary==true, specify coordinate system for the boundary. Default: if parent coordinatesystem is cartesian, just the cooresponding codim1 coordinatesystem, otherwise: nothing, requiring user specification for use of e.g. CellFinder with the subgrid.
  • project: project coordinates onto subgrid dimension

A subgrid is of type ExtendableGrid and stores two additional components: ParentGrid and NodeParents

ExtendableGrids.tricircumcenter!Method
tricircumcenter!(circumcenter, a, b, c)

Find the circumcenter of a triangle.

Derived from C source of Jonathan R Shewchuk <jrs@cs.cmu.edu>

Modified to return absolute coordinates.

ExtendableGrids.tryfixMethod
tryfix(
    a::Union{Array{T, 2}, VariableTargetAdjacency{T}}
) -> Any

Try to turn variable target adjacency into fixed target adjacency

ExtendableGrids.trytokenMethod
trytoken(tks::TokenStream, expected::String) -> Bool

Try for keyword token.

It token is missing, the token read is put back into stream, a value of false is returned and the next try/gettoken command continues at the same position,

Otherwise, true is returned, and reading continues after the token found.

ExtendableGrids.uniform_refineMethod
uniform_refine(
    source_grid::ExtendableGrid{T, K};
    store_parents
) -> ExtendableGrid

generates a new ExtendableGrid by uniform refinement of each cell in the given grid

uniform refinement rules are available for these AbstractElementGeometries:

  • Line1D (bisection into two subsegments)
  • Triangle2D (red refinement into four subtriangles)
  • Quadrilateral2D (into four subquadrilaterals)
  • Tetrahedron (into eight subtetrahedrons)
  • Hexahedron (into eight subhexahedrons)

if multiple geometries are in the mesh uniform refinement will only work if all refinement rules refine faces and edges (in 3D) equally (so no hanging nodes are created)

ExtendableGrids.update!Method
update!(
    grid::ExtendableGrid,
    T::Type{<:AbstractGridComponent}
) -> Any

Reinstantiate grid component (only if it exists)

ExtendableGrids.veryformMethod
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.

ExtendableGrids.veryformMethod
veryform(grid::ExtendableGrid{Tc,Ti},v,T::Type{<:AbstractGridAdjacency}) where{Tc,Ti}

Check proper type of adjacencies upon insertion

ExtendableGrids.writeVTKMethod
writeVTK(
    filename::String,
    grid::ExtendableGrid{Tc, Ti};
    kwargs...
) -> Vector{String}

exports grid and optional provided data as a vtk file

  • filename: filename of the exported file
  • grid: grid

Each '(key, value)' pair adds another data entry to the vtk file via WriteVTK functionality.