inla.mesh.2d.Rd since
23.08.18.
Use fmesher::fm_mesh_2d_inla() instead.
Create a triangle mesh based on initial point locations, specified or automatic boundaries, and mesh quality parameters.
inla.mesh.2d(
loc = NULL,
loc.domain = NULL,
offset = NULL,
n = NULL,
boundary = NULL,
interior = NULL,
max.edge = NULL,
min.angle = NULL,
cutoff = 1e-12,
max.n.strict = NULL,
max.n = NULL,
plot.delay = NULL,
crs = NULL
)Matrix of point locations to be used as initial triangulation
nodes. Can alternatively be a SpatialPoints or
SpatialPointsDataFrame object.
Matrix of point locations used to determine the domain
extent. Can alternatively be a SpatialPoints or
SpatialPointsDataFrame object.
The automatic extension distance. One or two values, for an inner and an optional outer extension. If negative, interpreted as a factor relative to the approximate data diameter (default=-0.10???)
The number of initial nodes in the automatic extensions (default=16)
A list of one or two inla.mesh.segment() objects
describing domain boundaries.
An inla.mesh.segment() object describing desired
interior edges.
The largest allowed triangle edge length. One or two values.
The smallest allowed triangle angle. One or two values. (Default=21)
The minimum allowed distance between points. Point at most as far apart as this are replaced by a single vertex prior to the mesh refinement step.
The maximum number of vertices allowed, overriding
min.angle and max.edge (default=-1, meaning no limit). One or
two values, where the second value gives the number of additional vertices
allowed for the extension.
The maximum number of vertices allowed, overriding
max.edge only (default=-1, meaning no limit). One or two values,
where the second value gives the number of additional vertices allowed for
the extension.
On Linux (and Mac if appropriate X11 libraries are
installed), specifying a nonnegative numeric value activates a rudimentary
plotting system in the underlying fmesher program, showing the
triangulation algorithm at work, with waiting time factor plot.delay
between each step.
On all systems, specifying any negative value activates displaying the result after each step of the multi-step domain extension algorithm.
An optional CRS or inla.CRS object
An inla.mesh object.
loc <- matrix(runif(10 * 2), 10, 2)
if (require("splancs")) {
boundary <- list(
inla.nonconvex.hull(loc, 0.1, 0.15),
inla.nonconvex.hull(loc, 0.2, 0.2)
)
offset <- NULL
} else {
boundary <- NULL
offset <- c(0.1, 0.2)
}
#> Loading required package: splancs
#>
#> Spatial Point Pattern Analysis Code in S-Plus
#>
#> Version 2 - Spatial and Space-Time analysis
#> Warning: `inla.nonconvex.hull()` was deprecated in INLA 23.08.18.
#> ℹ Please use `fmesher::fm_nonconvex_hull()` instead.
#> ℹ For more information, see
#> https://inlabru-org.github.io/fmesher/articles/inla_conversion.html
#> ℹ To silence these deprecation messages in old legacy code, set
#> `inla.setOption(fmesher.evolution.warn = FALSE)`.
#> ℹ To ensure visibility of these messages in package tests, also set
#> `inla.setOption(fmesher.evolution.verbosity = 'warn')`.
mesh <- inla.mesh.2d(loc, boundary = boundary, offset = offset, max.edge = c(0.05, 0.1))
#> Warning: `inla.mesh.2d()` was deprecated in INLA 23.08.18.
#> ℹ Please use `fmesher::fm_mesh_2d_inla()` instead.
#> ℹ For more information, see
#> https://inlabru-org.github.io/fmesher/articles/inla_conversion.html
#> ℹ To silence these deprecation messages in old legacy code, set
#> `inla.setOption(fmesher.evolution.warn = FALSE)`.
#> ℹ To ensure visibility of these messages in package tests, also set
#> `inla.setOption(fmesher.evolution.verbosity = 'warn')`.
plot(mesh)