Grid constructors
Tensor product simplex grids
ExtendableGrids.simplexgrid
— Functionfunction simplexgrid(coord::Array{Tc,2},
cellnodes::Array{Ti,2},
cellregions=ones(Ti, size(coord,2)),
bfacenodes=zeros(Ti,size(coord,1),0),
bfaceregions=zeros(Ti,length(bfacenodes))
) 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: (optional) n_tri vector of cell region markers
- bfacenodes: (optional) d $\times$ n_bf matrix of boundary facet - point incidences
- bfaceregions: (optional) 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
.
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.
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
|--|-----|-----|-----|-----|-----|-----|-----|--|
simplexgrid(X,Y; bregions=[1,2,3,4],cellregion=1)
Constructor for 2D grid from coordinate arrays.
Boundary region numbers count counterclockwise:
location | number |
---|---|
south | 1 |
east | 2 |
north | 3 |
west | 4 |
The keyword arguments allow to overwrite the default region numbers.
simplexgrid(X,Y,Z; bregions=[1,2,3,4,5,6],cellregion=1)
Constructor for 3D grid from coordinate arrays. Boundary region numbers:
location | number |
---|---|
south | 1 |
east | 2 |
north | 3 |
west | 4 |
bottom | 5 |
top | 6 |
The keyword arguments allow to overwrite the default region numbers.
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
.
simplexgrid(
file::String;
format,
kwargs...
) -> Union{Nothing, ExtendableGrid}
Read grid from file. Supported formats:
- "*.sg": pdelib sg files. Format versions:
format=v"2.0"
: long version with some unnecessary dataformat=v"2.1"
: shortened version only with cells, cellnodes, cellregions, bfacenodes, bfaceregionsformat=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 byinduce_edge_partitioning!
.
- "*.geo": gmsh geometry description (requires
using Gmsh
) - "*.msh": gmsh mesh (requires
using Gmsh
)
ExtendableGrids.glue
— Functionc=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
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
Various special grids
ExtendableGrids.grid_lshape
— Methodgrid_lshape(::Type{<:Triangle2D}; scale = [1,1], shift = [0,0])
Lshape domain
ExtendableGrids.grid_triangle
— Methodgrid_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_unitcube
— Methodgrid_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_unitcube
— Methodgrid_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_unitsquare
— Methodgrid_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_unitsquare
— Methodgrid_unitsquare(::Type{<:Triangle2D}; scale = [1,1], shift = [0,0])
Unit square as two triangles with four boundary regions (bottom, right, top, left)
ExtendableGrids.grid_unitsquare_mixedgeometries
— Methodgrid_unitsquare_mixedgeometries()
Unit suqare as mixed triangles and squares with four boundary regions (bottom, right, top, left)
ExtendableGrids.reference_domain
— Function 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.ringsector
— Methodringsector(rad,ang; eltype=Triangle2D)
Sector of ring or full ring (if ang[begin]-ang[end]≈2π
)