tin: TIN object

Description Usage Arguments Details Value Note Examples

View source: R/tin.R

Description

Initialise a TIN mesh object for use within TELEMAC.

Usage

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
tin(x, ...)

## S3 method for class 'character'
tin(x, ...)

## S3 method for class 'matrix'
tin(x, ..., ikle, ipobo)

## S3 method for class 'list'
tin(x, ..., s, s_brk, a, q = 30)

## S3 method for class 't2d_tin'
print(x, ...)

Arguments

x

Either: a character string providing the name of a SELAFIN (*.slf) or a Gmsh mesh (*.msh, format version 2) file of which the mesh will be extracted; a matrix of mesh points with 2 columns containing the x and y coordinates (arguments ikle and ipobo are required); a list with boundary and breakline definitions (see Details).

...

Arguments passed to or from other methods. If x is a list, further arguments passed to triangulate.

ikle

If x is a matrix of points: A matrix with 3 columns of indices referring to rows in x; each row represents a mesh element (triangle). In TELEMAC termed IKLE.

ipobo

If x is a matrix of points: A vector of indices referring to rows in x, each marking a boundary point. In TELEMAC termed IPOBO.

s

numeric, if x is a list: OPTIONAL value giving the resolution of vertices along the boundary line for triangulation. If not given, the points are used as they are supplied, otherwise line_spacing is called to ensure equal spacing of vertices with the given segment lengths.

s_brk

As s but for breaklines.

a

numeric, maximum triangle area; passed to triangulate. Default: squared spacing of points (either given as s or inferred from x$boundary).

q

numeric, minimum triangle angle; passed to triangulate. Default: 30 degrees.

Details

If x is a list this function creates a Triangulated Irregular Network (TIN) using function triangulate. The following list elements are required to perform the triangulation:

boundary

A matrix, data.frame, SpatialLines* or sf object with two columns, each row defining a point along the outer catchment boundary. Points are connected one-by-one to a line starting with the first point, i.e. make sure points are in the right order! The first and last point will be connected to close the boundary.

breaklines

OPTIONAL, a matrix, data.frame, SpatialLines* or sf object with three columns x and y, the x and y coordinates of vertices along the breaklines, and line, an identifier to identify individual breaklines.

Value

An object of class t2d_tin, which is a list with the following elements:

points

A matrix with the x and y coordinates (as columns) of mesh points.

triangles

A matrix with 3 columns of indices referring to rows in points; each row represents a mesh element (triangle).

edges

A matrix with 2 columns of indices referring to rows in points, the node points; each row represents an edge / segment of a triangle.

boundaries

A vector of indices referring to rows in points, each marking a point of the outer catchment boundary.

breaklines

A matrix with 2 columns of indices referring to rows in points, the vertices of the breaklines (used for mesh refinement during triangulation).

Note

Duplicated mesh points are silently removed.

Make sure breaklines do not intersect as this is not supported by the Triangle algorithm. A possible workaround to split intersecting breaklines in R using sf is shown in the examples.

If you want to construct a t2d_tin object and get the error Boundary points do not comply with requirements: [...] the reason might be that breaklines are too close to the boundary causing that points of the breaklines are used as boundary points which eventually results in a discontinuous outer boundary. Try to increase the distance of breaklines to the catchment boundary.

Examples

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
### BASIC FUNCTIONALITY ###
library(sf)

# load boundary as sf linestring
bnd <- st_read(system.file("dem/boundary_lagos.gpkg", package = "telemac"))

# create t2d_tin object
tin_obj <- tin(list(boundary = bnd), s = 90, a = 100^2, q = 30)

# inspection
tin_obj
str(tin_obj)
plot(tin_obj, pch = ".")

### DEALING WITH INTERSECTING BREAKLINES ###
library(sf)
library(tidyverse)

# example boundary and with intersecting breaklines
test_bnd <- st_linestring(
  matrix(c(seq(0,100,5),  rep(100,21), seq(100,0,-5),     rep(0,21),
           rep(0,21), seq(0,100,5),   rep(100,21), seq(100,0,-5)),
         ncol = 2)
) %>% st_sfc()
test_brk <- list(
  st_linestring(matrix(c(seq(0,100,5), rep(50,21)), ncol = 2)),
  st_linestring(matrix(c(rep(50,21), seq(0,100,5)), ncol = 2)),
  st_linestring(matrix(c(seq(30,60,5), rep(60,11),
                         rep(20,7), seq(20,70,5)), ncol = 2))) %>% st_sfc()

# get intersection points and define buffer of 2 around these points
pt_inters <- c(test_bnd, test_brk) %>%
  st_intersection() %>%
  st_collection_extract(type = "POINT") %>%
  st_buffer(2)

plot(test_bnd)
plot(test_brk, add = TRUE)
plot(pt_inters, add = TRUE)

# split breaklines
test_brk_unique <- st_difference(st_union(test_brk), st_union(pt_inters))
plot(test_bnd)
plot(test_brk_unique, add = TRUE)

# create mesh
tin_obj <- tin(list(boundary = test_bnd, breaklines = test_brk_unique),
               s = 2, s_brk = 2, a = 4, q = 30)
plot(tin_obj, pch = ".")

tpilz/telemac documentation built on Feb. 10, 2022, 2:12 p.m.