fm_bary {fmesher} | R Documentation |
Compute barycentric coordinates
Description
Identify knot intervals or triangles and compute barycentric coordinates
Usage
fm_bary(...)
## S3 method for class 'fm_bary'
fm_bary(bary, ..., extra_class = NULL)
## S3 method for class 'list'
fm_bary(bary, ..., extra_class = NULL)
## S3 method for class 'tbl_df'
fm_bary(bary, ..., extra_class = NULL)
## S3 method for class 'fm_mesh_1d'
fm_bary(mesh, loc, method = c("linear", "nearest"), restricted = FALSE, ...)
## S3 method for class 'fm_mesh_2d'
fm_bary(mesh, loc, crs = NULL, ..., max_batch_size = NULL)
## S3 method for class 'fm_mesh_3d'
fm_bary(mesh, loc, ..., max_batch_size = NULL)
## S3 method for class 'fm_lattice_2d'
fm_bary(mesh, loc, crs = NULL, ...)
## S3 method for class 'fm_lattice_Nd'
fm_bary(mesh, loc, ...)
Arguments
... |
Arguments forwarded to sub-methods. |
bary |
An |
extra_class |
character; If non- |
mesh |
|
loc |
Points for which to identify the containing interval/triangle, and
corresponding barycentric coordinates. May be a vector (for 1d) or a matrix
of raw coordinates, |
method |
character; method for defining the barycentric coordinates, "linear" (default) or "nearest" |
restricted |
logical, used for |
crs |
Optional crs information for |
max_batch_size |
integer; maximum number of points to process in a
single batch. This speeds up calculations by avoiding repeated large
internal memory allocations and data copies. The default, |
Value
A fm_bary
object, a tibble
with columns index
; either
vector of triangle indices (triangle meshes),
vector of knot indices (1D meshes, either for edges or individual knots), or
vector of lower left box indices (2D lattices),
and where
, a matrix of barycentric coordinates.
Methods (by class)
-
fm_bary(fm_bary)
: Returns thebary
input unchanged -
fm_bary(list)
: Converts alist
bary
tofm_bary
. In the list elements are unnamed, the namesindex
andwhere
are assumed. -
fm_bary(tbl_df)
: Converts atibble::tibble()
bary
tofm_bary
-
fm_bary(fm_mesh_1d)
: Return anfm_bary
object with elementsindex
(edge index vector pointing to the first knot of each edge) andwhere
(barycentric coordinates, 2-column matrices). Usefm_bary_simplex()
to obtain the corresponding endpoint knot indices.For
method = "nearest"
,index
contains the index of the nearest mesh knot, andwhere
is a single-column all-ones matrix. -
fm_bary(fm_mesh_2d)
: Anfm_bary
object with columnsindex
(vector of triangle indices) andwhere
(3-column matrix of barycentric coordinates). Points that were not found giveNA
entries inindex
andwhere
. -
fm_bary(fm_mesh_3d)
: Anfm_bary
object with columnsindex
(vector of triangle indices) andwhere
(4-column matrix of barycentric coordinates). Points that were not found giveNA
entries inindex
andwhere
. -
fm_bary(fm_lattice_2d)
: Anfm_bary
object with columnsindex
(vector of lattice cell indices) andwhere
(4-column matrix of barycentric coordinates). Points that are outside the lattice are givenNA
entries inindex
andwhere
. -
fm_bary(fm_lattice_Nd)
: Anfm_bary
object with columnsindex
(vector of lattice cell indices) andwhere
2^d
-column matrix of barycentric coordinates). Points that are outside the lattice are givenNA
entries inindex
andwhere
.
See Also
fm_bary_simplex()
, fm_bary_loc()
Examples
bary <- fm_bary(fm_mesh_1d(1:4), seq(0, 5, by = 0.5))
bary
str(fm_bary(fmexample$mesh, fmexample$loc_sf))
m <- fm_mesh_3d(
rbind(
c(1, 0, 0),
c(0, 1, 0),
c(0, 0, 1),
c(0, 0, 0)
),
matrix(c(1, 2, 3, 4), 1, 4)
)
b <- fm_bary(m, matrix(c(1, 1, 1) / 4, 1, 3))
str(fm_bary(fmexample$mesh, fmexample$loc_sf))