create.mesh.2D: Create a 2D triangular mesh

View source: R/mesh.R

create.mesh.2DR Documentation

Create a 2D triangular mesh

Description

This function is a wrapper of the Triangle library (http://www.cs.cmu.edu/~quake/triangle.html). It can be used to create a triangulation of the domain of interest starting from a list of points, to be used as triangles' vertices, and a list of segments, that define the domain boundary. The resulting mesh is a Constrained Delaunay triangulation. This is constructed in a way to preserve segments provided in the input segments without splitting them. This imput can be used to define the boundaries of the domain. If this imput is NULL, it generates a triangulation over the convex hull of the points. It is also possible to create a mesh.2D from the nodes locations and the connectivity matrix.

Usage

create.mesh.2D(nodes, nodesattributes = NA, segments = NA, holes = NA,
                     triangles = NA, order = 1, verbosity = 0)

Arguments

nodes

A #nodes-by-2 matrix containing the x and y coordinates of the mesh nodes.

nodesattributes

A matrix with #nodes rows containing nodes' attributes. These are passed unchanged to the output. If a node is added during the triangulation process or mesh refinement, its attributes are computed by linear interpolation using the attributes of neighboring nodes. This functionality is for instance used to compute the value of a Dirichlet boundary condition at boundary nodes added during the triangulation process.

segments

A #segments-by-2 matrix. Each row contains the row's indices in nodes of the vertices where the segment starts from and ends to. Segments are edges that are not splitted during the triangulation process. These are for instance used to define the boundaries of the domain. If this is input is NULL, it generates a triangulation over the convex hull of the points specified in nodes.

holes

A #holes-by-2 matrix containing the x and y coordinates of a point internal to each hole of the mesh. These points are used to carve holes in the triangulation, when the domain has holes.

triangles

A #triangles-by-3 (when order = 1) or #triangles-by-6 (when order = 2) matrix. This option is used when a triangulation is already available. It specifies the triangles giving the row's indices in nodes of the triangles' vertices and (when nodes = 2) also if the triangles' edges midpoints. The triangles' vertices and midpoints are ordered as described at
https://www.cs.cmu.edu/~quake/triangle.highorder.html. In this case the function create.mesh.2D is used to produce a complete mesh.2D object.

order

Either '1' or '2'. It specifies wether each mesh triangle should be represented by 3 nodes (the triangle' vertices) or by 6 nodes (the triangle's vertices and midpoints). These are respectively used for linear (order = 1) and quadratic (order = 2) Finite Elements. Default is order = 1.

verbosity

This can be '0', '1' or '2'. It indicates the level of verbosity in the triangulation process. When verbosity = 0 no message is returned during the triangulation. When verbosity = 2 the triangulation process is described step by step by displayed messages. Default is verbosity = 0.

Value

An object of the class mesh.2D with the following output:

  • nodesA #nodes-by-2 matrix containing the x and y coordinates of the mesh nodes.

  • nodesmarkersA vector of length #nodes, with entries either '1' or '0'. An entry '1' indicates that the corresponding node is a boundary node; an entry '0' indicates that the corresponding node is not a boundary node.

  • nodesattributesA matrix with #nodes rows containing nodes' attributes. These are passed unchanged from the input.

  • trianglesA #triangles-by-3 (when order = 1) or #triangles-by-6 (when order = 2) matrix. This option is used when a triangulation is already available. It specifies the triangles giving the indices in nodes of the triangles' vertices and (when nodes = 2) also if the triangles' edges midpoints. The triangles' vertices and midpoints are ordered as described at
    https://www.cs.cmu.edu/~quake/triangle.highorder.html.

  • segmentsmarkerA vector of length #segments with entries either '1' or '0'. An entry '1' indicates that the corresponding element in segments is a boundary segment; an entry '0' indicates that the corresponding segment is not a boundary segment.

  • edgesA #edges-by-2 matrix containing all the edges of the triangles in the output triangulation. Each row contains the row's indices in nodes, indicating the nodes where the edge starts from and ends to.

  • edgesmarkersA vector of lenght #edges with entries either '1' or '0'. An entry '1' indicates that the corresponding element in edge is a boundary edge; an entry '0' indicates that the corresponding edge is not a boundary edge.

  • neighborsA #triangles-by-3 matrix. Each row contains the indices of the three neighbouring triangles. An entry '-1' indicates that one edge of the triangle is a boundary edge.

  • holesA #holes-by-2 matrix containing the x and y coordinates of a point internal to each hole of the mesh. These points are used to carve holes in the triangulation, when the domain has holes.

  • orderEither '1' or '2'. It specifies wether each mesh triangle should be represented by 3 nodes (the triangle' vertices) or by 6 nodes (the triangle's vertices and midpoints). These are respectively used for linear (order = 1) and quadratic (order = 2) Finite Elements.

See Also

refine.mesh.2D, create.FEM.basis

Examples

library(fdaPDE)

## Upload the quasicirle2D data
data(quasicircle2D)
boundary_nodes = quasicircle2D$boundary_nodes
boundary_segments = quasicircle2D$boundary_segments
locations = quasicircle2D$locations
data = quasicircle2D$data

## Create mesh from boundary
## if the domain is convex it is sufficient to call:
mesh = create.mesh.2D(nodes = rbind(boundary_nodes, locations))
plot(mesh)

## if the domain is not convex, pass in addition the segments the compose the boundary:
mesh = create.mesh.2D(nodes = rbind(boundary_nodes, locations), segments = boundary_segments)

## Create mesh from data locations (without knowing the boundary)
mesh = create.mesh.2D(nodes = locations)
plot(mesh)
## In this case the domain is the convex hull of the data locations.
## Do this only if you do not have any information about the shape of the domain of interest.

fdaPDE documentation built on March 7, 2023, 5:28 p.m.