CompScienceMeshes.AbstractMeshType
`U::Int`: indicating dimension of the embedding space
`D1::Int`: one plus the manifold dimension
`T<:AbstractFloat`: the type of a vertex coordinate
CompScienceMeshes.MeshMethod
Converts an SComplex to Mesh.
Convenience function only. Use with caution.
(Only makes sense where the SComplex can originally be replaced with a Mesh.
In that case, the user is expected to use Mesh)
CompScienceMeshes.ReferenceSimplexType
ReferenceSimplex{Dimension, CoordType, NumVertices}

This domain is defined to bootstrap the quadrature generation strategy. The generic definition of numquads on a chart pulls back to the domain. For a limit set of reference domains, explicit quadrature rules are defined. The weights and points are then pushed forward to the configuaration space element over which integration is desired.

For more details see the implementation in quadpoints.jl

CompScienceMeshes.SegmentedAxisType

This type conforms to the mesh interface but is specialised for the case of a segment of the real axis subdived in equally sized intervals. Typical use is a discretisation of the time axis.

Base.:-Method
-mesh -> flipped_mesh

Create a mesh with opposite orientation.

Base.convertMethod
convert(::Type{NewT}, mesh::Mesh{U,D1,T}) where {NewT<:Real,U,D1,T}

Converts a Mesh with coordtype T to a Mesh with coordtype NewT.

Base.getindexMethod
getindex(simplex, I)
simplex[I]

Get the vertices at index I (scalar or array) defining the simplex

Base.lengthMethod
length(simplex)

Returns the number of vertices (equals dimension + 1)

Base.unionMethod
union(mesh1, mesh2, ...)

Create the topological union of two meshes. This requires them to be defined on the same vertex set. No geometric considerations are taken into account.

CompScienceMeshes.GSubdMeshMethod
GSubdMesh(mesh::Mesh)

Given a linear polygon mesh, construct a data structure for subdivision surfaces
CompScienceMeshes.Loop_subdivisionMethod
Loop_subdivsion(mesh) -> refinement

Construct a refinement of mesh by Loop subdivision scheme. Every face is subdived into four small faces and use weights to smooth the surface.

Only defined for 2D meshes.

CompScienceMeshes.barycentric_refinementMethod
barycentric refinement(mesh) -> refined_mesh

Create the mesh obtained by inserting an extra vertex in the barycenters of all cells and recusively creating fine cells by connecting the barycenter of a k-cell to the already constructed refined (k-1)-cells on its boundary.

CompScienceMeshes.bisecting_refinementMethod
bisecting_refinement(mesh) -> refinement

Construct a refinement of mesh by inserting a new vertex on every existing edge. Every face is subdived in four small faces.

Only defined for 2D meshes.

CompScienceMeshes.boundingboxMethod

Returns the boundingbox of a patch in terms of its center and halfsize.

function boundingbox{U,D,C,N,T}(p::Simplex{U,D,C,N,T}) -> center, halfsize
CompScienceMeshes.cartesianMethod
cartesian(neighborhood) -> point

Return the cartesian coordinates of the point where neighborhood is located.

CompScienceMeshes.cellpairsMethod
pairs = cellpairs(mesh, edges, dropjunctionpair=false)

Given a mesh and set of oriented edges from that mesh (as generated by skeleton), cellpairs will generate a 2 x K matrix, where K is the number of pairs and each column contains a pair of indices in the cell array of mesh that have one of the supplied edges in common.

Returns an array of pairs of indices, each pair corresponding to a pair of adjacent faces.

(If the mesh is oriented, the first row of facepairs will contain indices to the cell for which the corresponding edge has a positive relative orientation.

If a edge lies on the boundary of the mesh, and only has one neighboring cell, the second row of facepairs will contain -k with k the local index of the corresponding edge in its neighboring triangle.

If an edge has more than two neighboring cells (i.e. the edge is on a junction), all possible pairs of cells that have the junction edge in common are supplied. if dropjunctionpair == false then one of the possible pairs of cells is not recorded. This is done to avoid the creation of linearly dependent basis functions in the construction of boundary element methods for Maxwell's equations.)

CompScienceMeshes.cellsMethod
cells(mesh)

Return an iterable collection containing the cells making up the mesh.

CompScienceMeshes.centerMethod
center(simplex) -> neighborhood

Create the neighborhood at the center of the simplex, i.e. the point corresponding to parameter (1/(D+1), 1/(D+1), ...) where D is the simplex dimension.

CompScienceMeshes.chartMethod
chart(mesh, cell) -> cell_chart

Return a chart describing the supplied cell of mesh.

CompScienceMeshes.connectivityFunction
connectivity(faces, cells, op=sign)

Create a sparse matrix D of size numcells(cells) by numcells(faces) that contiains the connectivity info of the mesh. In particular D[m,k] is op(r) where r is the local index of face k in cell m. The sign of r is positive or negative depending on the relative orientation of face k in cell m.

For op=sign, the matrix returned is the classic connectivity matrix, i.e. the graph version of the exterior derivative.

CompScienceMeshes.dimensionMethod
dim = dimension(mesh)

Returns the dimension of the mesh. Note that this is the dimension of the cells, not of the surrounding space.

CompScienceMeshes.embeddingMethod
embedding(small::AbstractMesh, big::AbstractMesh) -> S::AbstractSparseArray
predicate(a,b) => extra constraints that must both be fullfiled by each
cell `a` and `b` that are found to exist in `small` and `big` respectively.
`a` and `b` are the indices of the corresponding edges. It is the user's
responsibility to ensure that the predicate allows for all cells of `small`
to be mapped to a cell in `big`. Specifying a predicate is useful for cases
where cells in `small` matches/overlaps with multiple cells in `big`
#TODO: change pred to accept simplices
CompScienceMeshes.find_neighborMethod
find_neighbor(faces,edges,F,EdgesIndices,orientation,Sedges)

Given a face find out the neighbor elements (share edge) and vertices.
CompScienceMeshes.flipMethod
flip(cell)

Change the orientation of a cell by interchanging the first to indices.

CompScienceMeshes.flip_normalMethod
flip_normal(simplex, sign)

Flips the normal of the simplex if sign is -1. Only on triangles embedded in 3D space

CompScienceMeshes.interiorFunction

Complement to boundary. This function selects those edges that have at least two faces adjacent. The case with more than two neighboring faces occurs on non-manifold structures (e.g. containing junctions)

CompScienceMeshes.interior_tpredicateMethod

Creates a predicate that can be used to check wheter an edge is interior to a surface (true) or on its boundary (false). This predicate is based on combinatorics. In particular it expects as argument a tuple of indices pointing into the vertex buffer of mesh.

CompScienceMeshes.interior_vpredicateMethod

Creates a predicate that can be used to check wheter a vertex is interior to a surface (true) or on its boundary (false). In particular it expects as argument an index pointing into the vertex buffer of mesh.

CompScienceMeshes.intersectionMethod
intersect(chartA, chartB) -> [chart1, chart2, ...]

Compute the intersection of two charts of equal dimension.

Compute an array of charts such that the disjoint union of those charts produces the intersection of the two charts provided as inputs. In particular the sum of integrals over the returned charts will equal the integral over the intersection of the two given charts.

CompScienceMeshes.intersection2Method
intersection2(triangleA, triangleB)

returns intersection in which the first operand inherits orientation from first argument and second argument inherits orientation of second argument.

CompScienceMeshes.isinclosureMethod
isinclosure(simplex, point) -> Bool

Determine whether point is in the closure of simplex. False positives are possible for points just outside of the simplex.

CompScienceMeshes.isinsideMethod
isinside(chart, point) -> Bool

Returns true is the given point is in the image of the given chart, false otherwise.

CompScienceMeshes.leftofMethod
inside(p,a,b)

Determines is p is on the interior side of the segment b of the boundary, assuming that the boundary is oriented counter-clockwise.

CompScienceMeshes.meshFunction
mesh(type, mdim, udim=mdim+1)

Returns an empty mesh with coordtype equal to type, of dimension mdim and embedded in a universe of dimension udim

CompScienceMeshes.meshcuboidMethod
meshcuboid(width, height, length, delta)

Creates a mesh for a cuboid of width (along the x-axis) width, height (along the y-axis) height and length (along the z-axis) length by parsing a .geo script incorporating these parameters into the GMSH mesher.

The target edge size is delta. physical => in {"TopPlate", "BottomPlate", "SidePlates", "OpenBox"} extracts and returns only the specified part of the cuboid

CompScienceMeshes.meshrectangleMethod
meshrectangle(width, height, delta, udim)

Create a mesh for a rectangle of width (along the x-axis) width and height (along the y-axis) height.

The target edge size is delta and the dimension of the embedding universe is udim (>= 2).

The mesh is oriented such that the normal is pointing down. This is subject to change.

CompScienceMeshes.meshsphereMethod
meshsphere(radius, delta)
meshsphere(;radius, h)

Create a mesh of a sphere of radius radius by parsing a .geo script incorporating these parameters into the GMSH mesher.

The target edge size is delta.

CompScienceMeshes.mirrorMethod
mirror(vertex, normal, anchor)

Mirror vertex across a plane defined by its normal and a containing point.

CompScienceMeshes.numverticesMethod
numvertices(mesh)

Returns the number of vertices in the mesh.

Note: this is the number of vertices in the vertex buffer and might include floatin vertices or vertices not appearing in any cell. In other words the following is not necessarily true:

    numvertices(mesh) == numcells(skeleton(mesh,0))
CompScienceMeshes.octreeMethod

Store the k-cells of a mesh in an octree.

function octree{U,D,T}(mesh::Mesh{U,D,T}, kcells::Array{Int,2})
CompScienceMeshes.overlap_gpredicateMethod
overlap_gpredicate(γ::Mesh) -> (patch -> Bool)

Create a predicate that for a given patch determinees if it overlaps with the provided target mesh γ.

CompScienceMeshes.permute_verticesMethod
permute_simplex(simplex,permutation)

Permutation is a Vector v which sets the v[i]-th vertex at the i-th place.

Return Simplex with permuted vertices list, tangents are recalculated, normal is kept the same

CompScienceMeshes.pointMethod
point(type, xs...)

Create point of default type and supplied precision for the coordinates

CompScienceMeshes.quadpointsMethod
quadpoints(refspace, charts, rules)

Computed a matrix of vectors containing (weight, point, value) triples that can be used in numerical integration over the elements described by the charts. Internally, this method used quadpoints(chart, rule) to retrieve the points and weights for a certain quadrature rule over chart.

CompScienceMeshes.quadpointsMethod
pw = quadpoints(chart, rule)

Returns a collection of (point, weight) tuples corresponding to the numerical quadrature rule defined on the domain of chart. The weights returned take into account the Jacobian determinant resulting from mapping from the reference domain to the configuration space.

Functions can be integrated like:

PW = quadpoints(chart, rule)
I = sum(pw[2]*f(pw[1]) for pw in PW)
CompScienceMeshes.read_TRI_meshFunction
read_TRI_mesh(mesh_filename) -> mesh::Mesh

Imports a surface mesh (stored in an ASCII file named mesh_filename) into a Mesh object (i.e. node list and element list).

NOTE: The contents of mesh_filename must include the file extension, and the file must be stored in the current directory.

CompScienceMeshes.read_gmsh3d_meshMethod
read_gmsh3d_mesh(iostream) -> mesh::Mesh

Reads the mesh nodes and elements stored in the input .msh file (io, output by GMSH) into arrays of node vectors and vertex vectors respectively.

Returns an object mesh of type Mesh, comprising both vector arrays.

CompScienceMeshes.read_gmsh_meshMethod
read_gmsh_mesh(iostream) -> mesh::Mesh

Reads the mesh nodes and elements stored in the input .msh file (io, output by GMSH) into arrays of node vectors and vertex vectors respectively.

Returns an object mesh of type Mesh, comprising both vector arrays.

CompScienceMeshes.readmeshMethod
readmesh(filename)

Reads a mesh in in format from filename. The format follows:

1
V C
x1_1    x1_2    ... x1_U
x2_1    x2_2    ... x2_U
...
xV_1    xV_2    ... xV_U
i1_1    i1_2    ... i1_D1
i2_1    i2_2    ... i2_D1
...
iC_1    iC_2    ... iC_D1

where U is the universedimension of the mesh, D1 the dimension of the mesh plus one, V the number of vertices, and C the number of cells in the mesh.

CompScienceMeshes.restrictionMethod
restriction(submesh, supermesh)

Computes the restriction matrix relative to a submesh submesh of supermesh.

The restriction matrix has size (m,n), where

m == numcells(submesh)
n == numcells(supermesh)

It has entries 1 at location [i,j] iff cell i of submesh equals cell j of supermesh. The remaining entries are zero.

This matrix is referred to as the restriction matrix because if it acts on an array of samples taken at the cells of supermesh is selects out the samples in the cells that are retained in submesh, taking into account any renumbering. Likewise, its transpose is sometimes referred to as the extension-by-zero operator because it maps arrays of samples taken in the cells of submesh into an array of samples taken in the cells of supermesh by inserting zeros at cells that were not retained in submesh.

CompScienceMeshes.simplexMethod
simplex(vertices)
simplex(v1, v2, ...)
simplex(vertices, Val{D})

Build a D-dimensional simplex. The vertices can be passed in an array (static or dynamic), or supplied separately. If the length of the array is not part of its type, the speed of the construction can be improved by supplying an extra Val{D} argument. In case it is not clear from the context whether the vertex array is dynamically or statically sized, use the third form as it will not incur notable performance hits.

Note that D is the dimension of the simplex, i.e. the number of vertices supplied minus one.

CompScienceMeshes.skeletonMethod
skeleton(mesh, dim)

Returns a mesh comprising the dim-dimensional sub cells of mesh. For example to retrieve the edges of a given surface mesh,

edges = skelton(mesh, 1)
CompScienceMeshes.submeshMethod

Returns a mesh on the same vertexbuffer as the input mesh. The submesh will be a mesh of dimension k containing all the k-cells that are in mesh and that fulfill the predicate pred.

pred is a function with signature pred(cell) -> Bool returning true if the simplex is to be added to the submesh under construction.

CompScienceMeshes.submeshMethod
submesh(selection, mesh)

Create a submesh from mesh comprising those elements that overlap with elements from selection. It is assumed that selection and mesh have the same dimension.

CompScienceMeshes.sutherlandhodgmanMethod
sutherlandhodgman(subject, clipper)

Compute the intersection of two coplanar triangles, potentially embedded in a higher dimensional space.

CompScienceMeshes.sutherlandhodgman_keep_clippingsMethod
sutherlandhodgman_keep_clippings(subject, clipper)

Compute both the intersection and the setminus of two convex polygons subject and clipper. The polygons should be supplied as vectors of points. The intersection is returned as a vector of points and the clippings as a vector of vectors of points.

CompScienceMeshes.trgaussMethod
trgauss(n) -> (u,w)

Returns the n-th triangle quadrature rule. Returns a Matrix u of size (Q,2) with Q the number of quadrature points and a Vector w of size (Q,) containing the quadrature weights.

CompScienceMeshes.universedimensionMethod
universedimension(mesh)

Returns the dimension of the surrounding space. Equals the number of coordinates required to describe a vertex.

CompScienceMeshes.vertextocellmapMethod
vertextocellmap(mesh) -> vertextocells, numneighbors

Computed an V×M array vertextocells where V is the number of vertices and M is the maximum number of cells adjacent to any given vertex such that vertextocells[v,i] is the index in the cells of mesh of the ith cell adjacent to teh v-th vertex. numneighbors[v] contains the number of cells adjacent to the v-th vertex.

This method allows e.g. for the efficient computation of the connectivity matrix of the mesh.

CompScienceMeshes.weldFunction
weld(mesh1, mesh2, ...;boundary=false) -> welded_mesh

Build a mesh by welding or pasting together the inputs. Vertices from different meshes that coincide up to the tolerance will be merged into one. The order cells appear in the output mesh is equal to the order in the inputs. boundary = true; will merge only vertices of different meshes if one of them exists on the 'open' boundary of the mesh.