Title: | Manipulations of Triangular Meshes Based on the 'VCGLIB' API |
---|---|
Description: | Operations on triangular meshes based on 'VCGLIB'. This package integrates nicely with the R-package 'rgl' to render the meshes processed by 'Rvcg'. The Visualization and Computer Graphics Library (VCG for short) is an open source portable C++ templated library for manipulation, processing and displaying with OpenGL of triangle and tetrahedral meshes. The library, composed by more than 100k lines of code, is released under the GPL license, and it is the base of most of the software tools of the Visual Computing Lab of the Italian National Research Council Institute ISTI <https://vcg.isti.cnr.it/>, like 'metro' and 'MeshLab'. The 'VCGLIB' source is pulled from trunk <https://github.com/cnr-isti-vclab/vcglib> and patched to work with options determined by the configure script as well as to work with the header files included by 'RcppEigen'. |
Authors: | Stefan Schlager [aut, cre, cph], Girinon Francois [ctb], Tim Schaefer [ctb], Zhengjia Wang [ctb] |
Maintainer: | Stefan Schlager <[email protected]> |
License: | GPL (>= 2) | file LICENSE |
Version: | 0.24 |
Built: | 2024-10-30 08:30:03 UTC |
Source: | https://github.com/zarquon42b/rvcg |
Provides meshing functionality from vcglib (meshlab) for R. E.g. mesh smoothing, mesh decimation, closest point search.
Package: | Rvcg |
Type: | Package |
Version: | 0.24 |
Date: | 2024-09-19 |
License: | GPL |
LazyLoad: | yes |
Stefan Schlager
Maintainer: Stefan Schlager <[email protected]>
To be announced
Useful links:
Report bugs at https://github.com/zarquon42b/Rvcg/issues
check the orientation of a mesh assuming that expansion along normals increases centroid size
checkFaceOrientation(x, offset = NULL)
checkFaceOrientation(x, offset = NULL)
x |
mesh of class mesh3d |
offset |
numeric: amount to offset the mesh along the vertex normals. If NULL a reasonable value will be estimated. |
assuming that a correctly (i.e outward) oriented mesh increases its centroid size when 'growing' outwards, this function tests whether this is the case.
returns TRUE if mesh is oriented correctly and FALSE otherwise
data(dummyhead) ## now we invert faces inwards checkFaceOrientation(dummyhead.mesh) if (requireNamespace("Morpho", quietly = TRUE)) { dummyinward <- Morpho::invertFaces(dummyhead.mesh) checkFaceOrientation(dummyinward) }
data(dummyhead) ## now we invert faces inwards checkFaceOrientation(dummyhead.mesh) if (requireNamespace("Morpho", quietly = TRUE)) { dummyinward <- Morpho::invertFaces(dummyhead.mesh) checkFaceOrientation(dummyinward) }
A triangular mesh representing a dummyhead - called by data(dummyhead)
dummyhead.mesh
: triangular mesh representing a dummyhead.
dummyhead.lm
: landmarks on mesh 'dummyhead'
A triangular mesh representing a human face - called by data(humface)
humface
: triangular mesh representing a human face.
humfaceClean
: triangular mesh representing a human face but without errors or isolated pieces.
humface.lm
: landmarks on mesh 'humface'- called by data(humface)
print number of vertices and triangular faces of a mesh
meshInfo(x)
meshInfo(x)
x |
triangular mesh |
checks for existance and validity of vertices, faces and vertex normals of an object of class "mesh3d"
meshintegrity(mesh, facecheck = FALSE, normcheck = FALSE)
meshintegrity(mesh, facecheck = FALSE, normcheck = FALSE)
mesh |
object of class mesh3d |
facecheck |
logical: check the existence of valid triangular faces |
normcheck |
logical: check the existence of valid normals |
if mesh data are valid, the mesh is returned, otherwise it stops with an error message.
get number of vertices from a mesh
nfaces(x)
nfaces(x)
x |
triangular mesh |
integer: number of triangular faces
get number of vertices from a mesh
nverts(x)
nverts(x)
x |
triangular mesh |
integer: number of vertices
Find all intersections by tracing rays through mesh #
raysearchMulti(x, mesh, maxtol = 1e+15, threads = 1, offset = 0.001)
raysearchMulti(x, mesh, maxtol = 1e+15, threads = 1, offset = 0.001)
x |
a triangular mesh of class 'mesh3d' or a list containing vertices and vertex normals (fitting the naming conventions of 'mesh3d'). In the second case x must contain x$vb = 3 x n matrix containing 3D-coordinates and x$normals = 3 x n matrix containing normals associated with x$vb. |
mesh |
triangular mesh to be intersected. |
maxtol |
maximum distance to search along ray |
threads |
number of threads used during search. |
offset |
amount to offset the hit points along the ray to find the next intersection. This is needed to avoid finding the same intersection over and over again. |
This function iteratively uses vcgRaySearch
to find all intersections of rays and a given surface mesh.
list with following items:
intersects |
a list containing the result of |
hits |
Vector containging number of intersections for each ray |
## Not run: require(Morpho); require(rgl) data(humface) humface1 <- scalemesh(humface,size=1.1) mesh <- mergeMeshes(humface,humface1) #get normals of landmarks x <- vcgClost(humface.lm, humface) # offset landmarks along their normals for a negative amount of -5mm x$vb[1:3,] <- x$vb[1:3,]+x$normals[1:3,]*-5 myint <- raysearchMulti(x,mesh) wire3d(mesh,col="white") spheres3d(vert2points(x),radius=0.5,col=3) plotNormals(x,length=55,lwd=2) for (i in 1:length(myint$intersects)) { spheres3d(vert2points(myint$intersects[[i]])[which(as.logical(myint$intersects[[i]]$quality)),] ,col=i) } ## End(Not run)
## Not run: require(Morpho); require(rgl) data(humface) humface1 <- scalemesh(humface,size=1.1) mesh <- mergeMeshes(humface,humface1) #get normals of landmarks x <- vcgClost(humface.lm, humface) # offset landmarks along their normals for a negative amount of -5mm x$vb[1:3,] <- x$vb[1:3,]+x$normals[1:3,]*-5 myint <- raysearchMulti(x,mesh) wire3d(mesh,col="white") spheres3d(vert2points(x),radius=0.5,col=3) plotNormals(x,length=55,lwd=2) for (i in 1:length(myint$intersects)) { spheres3d(vert2points(myint$intersects[[i]])[which(as.logical(myint$intersects[[i]]$quality)),] ,col=i) } ## End(Not run)
create a search structure from a matrix of coordinates and one of directional vectors to be processed by vcgRaySearch
setRays(coords, dirs)
setRays(coords, dirs)
coords |
k x 3 matrix (or a vector of length 3) containing the starting points of the rays |
dirs |
k x 3 matrix (or a vector of length 3) containing the directons of the rays. The i-th row of |
an object of class "mesh3d" (without faces) and the vertices representing the starting points of the rays and the normals storing the directions.
compute surface area of a triangular mesh
vcgArea(mesh, perface = FALSE)
vcgArea(mesh, perface = FALSE)
mesh |
triangular mesh of class mesh3d |
perface |
logical: if TRUE, a list containing the overall area, as well as the individual per-face area are reported. |
surface area of mesh
data(humface) vcgArea(humface)
data(humface) vcgArea(humface)
Ball pivoting surface reconstruction
vcgBallPivoting( x, radius = 0, clustering = 0.2, angle = pi/2, deleteFaces = FALSE )
vcgBallPivoting( x, radius = 0, clustering = 0.2, angle = pi/2, deleteFaces = FALSE )
x |
k x 3 matrix or object of class mesh3d |
radius |
The radius of the ball pivoting (rolling) over the set of points. Gaps that are larger than the ball radius will not be filled; similarly the small pits that are smaller than the ball radius will be filled. 0 = autoguess. |
clustering |
Clustering radius (fraction of ball radius). To avoid the creation of too small triangles, if a vertex is found too close to a previous one, it is clustered/merged with it. |
angle |
Angle threshold (radians). If we encounter a crease angle that is too large we should stop the ball rolling. |
deleteFaces |
in case x is a mesh and |
triangular face of class mesh3d
if (requireNamespace("Morpho", quietly = TRUE)) { require(Morpho) data(nose) nosereko <- vcgBallPivoting(shortnose.lm) }
if (requireNamespace("Morpho", quietly = TRUE)) { require(Morpho) data(nose) nosereko <- vcgBallPivoting(shortnose.lm) }
get barycenters of all faces of a triangular mesh
vcgBary(mesh)
vcgBary(mesh)
mesh |
triangular mesh of class "mesh3d" |
n x 3 matrix containing 3D-coordinates of the barycenters (where n is the number of faces in mesh
.
data(humface) bary <- vcgBary(humface) ## Not run: require(rgl) points3d(bary,col=2) wire3d(humface) ## End(Not run)
data(humface) bary <- vcgBary(humface) ## Not run: require(rgl) points3d(bary,col=2) wire3d(humface) ## End(Not run)
Detect faces and vertices at the borders of a mesh and mark them.
vcgBorder(mesh)
vcgBorder(mesh)
mesh |
triangular mesh of class "mesh3d" |
bordervb |
logical: vector containing boolean value for each vertex, if it is a border vertex. |
borderit |
logical: vector containing boolean value for each face, if it is a border vertex. |
Stefan Schlager
data(humface) borders <- vcgBorder(humface) ## view border vertices ## Not run: require(rgl) points3d(t(humface$vb[1:3,])[which(borders$bordervb == 1),],col=2) wire3d(humface) require(rgl) ## End(Not run)
data(humface) borders <- vcgBorder(humface) ## view border vertices ## Not run: require(rgl) points3d(t(humface$vb[1:3,])[which(borders$bordervb == 1),],col=2) wire3d(humface) require(rgl) ## End(Not run)
Apply several cleaning algorithms to surface meshes
vcgClean(mesh, sel = 0, tol = 0, silent = FALSE, iterate = FALSE)
vcgClean(mesh, sel = 0, tol = 0, silent = FALSE, iterate = FALSE)
mesh |
triangular mesh of class 'mesh3d' |
sel |
integer vector selecting cleaning type (see "details"), |
tol |
numeric value determining Vertex Displacement Ratio used for splitting non-manifold vertices. |
silent |
logical, if TRUE no console output is issued. |
iterate |
logical: if TRUE, vcgClean is repeatedly run until nothing more is to be cleaned (see details). |
the vector sel determines which operations are performed in which order. E.g. removing degenerate faces may generate unreferenced vertices, thus the ordering of cleaning operations is important, multiple calls are possible (sel=c(1,3,1) will remove unreferenced vertices twice). available options are:
0 = only duplicated vertices and faces are removed
1 = remove unreferenced vertices
2 = Remove non-manifold Faces
3 = Remove degenerate faces
4 = Remove non-manifold vertices
5 = Split non-manifold vertices by threshold
6 = merge close vertices (radius=tol
)
7 = coherently orient faces
CAVEAT: sel=6 will not work keep vertex colors
cleaned mesh with an additional entry
remvert |
vector of length = number of vertices before cleaning. Entries = 1 indicate that this vertex was removed; 0 otherwise. |
data(humface) cleanface <- humface ##add duplicated faces cleanface$it <- cbind(cleanface$it, cleanface$it[,1:100]) ## add duplicated vertices cleanface$vb <- cbind(cleanface$vb,cleanface$vb[,1:100]) ## ad unreferenced vertices cleanface$vb <- cbind(cleanface$vb,rbind(matrix(rnorm(18),3,6),1)) cleanface <- vcgClean(cleanface, sel=1)
data(humface) cleanface <- humface ##add duplicated faces cleanface$it <- cbind(cleanface$it, cleanface$it[,1:100]) ## add duplicated vertices cleanface$vb <- cbind(cleanface$vb,cleanface$vb[,1:100]) ## ad unreferenced vertices cleanface$vb <- cbind(cleanface$vb,rbind(matrix(rnorm(18),3,6),1)) cleanface <- vcgClean(cleanface, sel=1)
For a set of 3D-coordinates/triangular mesh, the closest matches on a target surface are determined and normals at as well as distances to that point are calculated.
vcgClost( x, mesh, sign = TRUE, barycentric = FALSE, smoothNormals = FALSE, borderchk = FALSE, tol = 0, facenormals = FALSE, ... )
vcgClost( x, mesh, sign = TRUE, barycentric = FALSE, smoothNormals = FALSE, borderchk = FALSE, tol = 0, facenormals = FALSE, ... )
x |
k x 3 matrix containing 3D-coordinates or object of class "mesh3d". |
mesh |
triangular surface mesh stored as object of class "mesh3d". |
sign |
logical: if TRUE, signed distances are returned. |
barycentric |
logical: if TRUE, barycentric coordinates of the hit points are returned. |
smoothNormals |
logical: if TRUE, laplacian smoothed normals are used. |
borderchk |
logical: request checking if the hit face is at the border of the mesh. |
tol |
maximum distance to search. If distance is beyond that, the original point will be kept and the distance set to NaN. If tol = 0, tol is set to 2*diagonal of the bounding box of |
facenormals |
logical: if TRUE only the facenormal of the face the closest point has hit is returned, the weighted average of the surrounding vertex normals otherwise. |
... |
additional parameters, currently unused. |
returns an object of class "mesh3d" with:
vb |
4 x n matrix containing n vertices as homolougous coordinates. |
normals |
4 x n matrix containing vertex normals. |
quality |
numeric vector containing distances to target. |
it |
3 x m integer matrix containing vertex indices forming triangular faces.Only available, when x is a mesh. |
border |
integer vector of length n: if borderchk = TRUE, for each clostest point the value will be 1 if the hit face is at the border of the target mesh and 0 otherwise. |
barycoords |
3 x m Matrix containing barycentric coordinates of closest points; only available if barycentric=TRUE. |
faceptr |
vector of face indeces on which the closest points are located |
If large part of the reference mesh are far away from the target
surface, calculation can become very slow. In that case, the function
vcgClostKD
will be significantly faster.
Stefan Schlager
Baerentzen, Jakob Andreas. & Aanaes, H., 2002. Generating Signed Distance Fields From Triangle Meshes. Informatics and Mathematical Modelling.
data(humface) clost <- vcgClost(humface.lm, humface)
data(humface) clost <- vcgClost(humface.lm, humface)
For a set of 3D-coordinates/triangular mesh, the closest matches on a target surface are determined (by using KD-tree search) and normals at as well as distances to that point are calculated.
vcgClostKD( x, mesh, sign = TRUE, barycentric = FALSE, smoothNormals = FALSE, borderchk = FALSE, k = 50, nofPoints = 16, maxDepth = 64, angdev = NULL, weightnorm = FALSE, facenormals = FALSE, threads = 1, ... )
vcgClostKD( x, mesh, sign = TRUE, barycentric = FALSE, smoothNormals = FALSE, borderchk = FALSE, k = 50, nofPoints = 16, maxDepth = 64, angdev = NULL, weightnorm = FALSE, facenormals = FALSE, threads = 1, ... )
x |
k x 3 matrix containing 3D-coordinates or object of class "mesh3d". |
mesh |
triangular surface mesh stored as object of class "mesh3d". |
sign |
logical: if TRUE, signed distances are returned. |
barycentric |
logical: if TRUE, barycentric coordinates of the hit points are returned. |
smoothNormals |
logical: if TRUE, laplacian smoothed normals are used. |
borderchk |
logical: request checking if the hit face is at the border of the mesh. |
k |
integer: check the kdtree for the |
nofPoints |
integer: number of points per cell in the kd-tree (don't change unless you know what you are doing!) |
maxDepth |
integer: depth of the kd-tree (don't change unless you know what you are doing!) |
angdev |
maximum deviation between reference and target normals. If the none of the k closest triangles match this criterion, the closest point on the closest triangle is returned but the corresponding distance in $quality is set to 1e5. |
weightnorm |
logical if angdev is set, this requests the normal of the closest points to be estimated by weighting the surrounding vertex normals. Otherwise, simply the hit face's normal is used (faster but slightly less accurate) |
facenormals |
logical: if TRUE only the facenormal of the face the closest point has hit is returned, the weighted average of the surrounding vertex normals otherwise. |
threads |
integer: threads to use in closest point search. |
... |
additional parameters, currently unused. |
returns an object of class "mesh3d" with:
vb |
4 x n matrix containing n vertices as homolougous coordinates. |
normals |
4 x n matrix containing vertex normals. |
quality |
numeric vector containing distances to target. |
it |
3 x m integer matrix containing vertex indices forming triangular faces.Only available, when x is a mesh. |
border |
integer vector of length n: if borderchk = TRUE, for each clostest point the value will be 1 if the hit face is at the border of the target mesh and 0 otherwise. |
barycoords |
3 x m Matrix containing barycentric coordinates of closest points; only available if barycentric=TRUE. |
Other than vcgClost
this does not search a grid, but first uses a KD-tree search to find the k
closest barycenters for each point and then searches these faces for the closest match.
Stefan Schlager
Baerentzen, Jakob Andreas. & Aanaes, H., 2002. Generating Signed Distance Fields From Triangle Meshes. Informatics and Mathematical Modelling.
search a KD-tree from Barycenters for multiple closest point searches on a mesh
vcgClostOnKDtreeFromBarycenters( x, query, k = 50, sign = TRUE, barycentric = FALSE, borderchk = FALSE, angdev = NULL, weightnorm = FALSE, facenormals = FALSE, threads = 1 )
vcgClostOnKDtreeFromBarycenters( x, query, k = 50, sign = TRUE, barycentric = FALSE, borderchk = FALSE, angdev = NULL, weightnorm = FALSE, facenormals = FALSE, threads = 1 )
x |
object of class "vcgKDtreeWithBarycenters" |
query |
matrix or triangular mesh containing coordinates |
k |
integer: check the kdtree for the |
sign |
logical: if TRUE, signed distances are returned. |
barycentric |
logical: if TRUE, barycentric coordinates of the hit points are returned. |
borderchk |
logical: request checking if the hit face is at the border of the mesh. |
angdev |
maximum deviation between reference and target normals. If the none of the k closest triangles match this criterion, the closest point on the closest triangle is returned but the corresponding distance in $quality is set to 1e5. |
weightnorm |
logical if angdev is set, this requests the normal of the closest points to be estimated by weighting the surrounding vertex normals. Otherwise, simply the hit face's normal is used (faster but slightly less accurate) |
facenormals |
logical: if TRUE only the facenormal of the face the closest point has hit is returned, the weighted average of the surrounding vertex normals otherwise. |
threads |
integer: threads to use in closest point search. |
returns an object of class "mesh3d" with:
vb |
4 x n matrix containing n vertices as homolougous coordinates. |
normals |
4 x n matrix containing vertex normals. |
quality |
numeric vector containing distances to target. |
it |
3 x m integer matrix containing vertex indices forming triangular faces.Only available, when x is a mesh. |
border |
integer vector of length n: if borderchk = TRUE, for each clostest point the value will be 1 if the hit face is at the border of the target mesh and 0 otherwise. |
barycoords |
3 x m Matrix containing barycentric coordinates of closest points; only available if barycentric=TRUE. |
Stefan Schlager
vcgCreateKDtreeFromBarycenters, vcgSearchKDtree, vcgCreateKDtree
create a KD-tree
vcgCreateKDtree(mesh, nofPointsPerCell = 16, maxDepth = 64)
vcgCreateKDtree(mesh, nofPointsPerCell = 16, maxDepth = 64)
mesh |
matrix or triangular mesh containing coordinates |
nofPointsPerCell |
number of points per kd-cell |
maxDepth |
maximum tree depth |
returns an object of class vcgKDtree containing external pointers to the tree and the target points
data(humface) mytree <- vcgCreateKDtree(humface)
data(humface) mytree <- vcgCreateKDtree(humface)
create a KD-tree from Barycenters for multiple closest point searches on a mesh
vcgCreateKDtreeFromBarycenters(mesh, nofPointsPerCell = 16, maxDepth = 64)
vcgCreateKDtreeFromBarycenters(mesh, nofPointsPerCell = 16, maxDepth = 64)
mesh |
matrix or triangular mesh containing coordinates |
nofPointsPerCell |
number of points per kd-cell |
maxDepth |
maximum tree depth |
returns an object of class vcgKDtreeWithBarycenters containing external pointers to the tree, the barycenters and the target mesh
vcgClostOnKDtreeFromBarycenters, vcgSearchKDtree, vcgCreateKDtree
## Not run: data(humface);data(dummyhead) barytree <- vcgCreateKDtreeFromBarycenters(humface) closest <- vcgClostOnKDtreeFromBarycenters(barytree,dummyhead.mesh,k=50,threads=1) ## End(Not run)
## Not run: data(humface);data(dummyhead) barytree <- vcgCreateKDtreeFromBarycenters(humface) closest <- vcgClostOnKDtreeFromBarycenters(barytree,dummyhead.mesh,k=50,threads=1) ## End(Not run)
calculate curvature of faces/vertices of a triangular mesh using various methods.
vcgCurve(mesh)
vcgCurve(mesh)
mesh |
triangular mesh (object of class 'mesh3d') |
gaussvb |
per vertex gaussian curvature |
meanvb |
per vertex mean curvature |
RMSvb |
per vertex RMS curvature |
gaussitmax |
per face maximum gaussian curvature of adjacent vertices |
borderit |
per face information if it is on the mesh's border (0=FALSE, 1=TRUE) |
bordervb |
per vertex information if it is on the mesh's border (0=FALSE, 1=TRUE) |
meanitmax |
per face maximum mean curvature of adjacent vertices |
K1 |
Principal Curvature 1 |
K2 |
Principal Curvature 2 |
data(humface) curv <- vcgCurve(humface) ##visualise per vertex mean curvature ## Not run: require(Morpho) meshDist(humface,distvec=curv$meanvb,from=-0.2,to=0.2,tol=0.01) ## End(Not run)
data(humface) curv <- vcgCurve(humface) ##visualise per vertex mean curvature ## Not run: require(Morpho) meshDist(humface,distvec=curv$meanvb,from=-0.2,to=0.2,tol=0.01) ## End(Not run)
Compute pseudo-geodesic distances on a triangular mesh
vcgDijkstra(x, vertpointer, maxdist = NULL)
vcgDijkstra(x, vertpointer, maxdist = NULL)
x |
triangular mesh of class |
vertpointer |
integer: references indices of vertices on the mesh, typically only a single query vertex. |
maxdist |
positive scalar double, the maximal distance to travel along the mesh when computing distances. Leave at |
returns a vector of shortest distances for each of the vertices to one of the vertices referenced in vertpointer
. If maxdis
t is in use (not NULL
), the distance values for vertices outside the requested maxdist
are not computed and appear as 0
.
Make sure to have a clean manifold mesh. Note that this computes the length of the pseudo-geodesic path (following the edges) between the two vertices.
## Compute geodesic distance between all mesh vertices and the first vertex of a mesh data(humface) geo <- vcgDijkstra(humface,1) if (interactive()) { require(Morpho);require(rgl) meshDist(humface,distvec = geo) spheres3d(vert2points(humface)[1,],col=2) }
## Compute geodesic distance between all mesh vertices and the first vertex of a mesh data(humface) geo <- vcgDijkstra(humface,1) if (interactive()) { require(Morpho);require(rgl) meshDist(humface,distvec = geo) spheres3d(vert2points(humface)[1,],col=2) }
Compute normalized face normals for a mesh.
vcgFaceNormals(mesh)
vcgFaceNormals(mesh)
mesh |
triangular mesh of class 'mesh3d', from |
3xn numeric matrix of face normals for the mesh, where n
is the number of faces.
data(humface); hf_facenormals <- vcgFaceNormals(humface);
data(humface); hf_facenormals <- vcgFaceNormals(humface);
Compute geodesic path and path length between vertices on a mesh
vcgGeodesicPath(x, source, targets, maxdist = 1e+06)
vcgGeodesicPath(x, source, targets, maxdist = 1e+06)
x |
triangular mesh of class |
source |
scalar positive integer, the source vertex index. |
targets |
positive integer vector, the target vertex indices. |
maxdist |
numeric, the maximal distance to travel along the mesh edges during geodesic distance computation. |
named list with two entries as follows. 'paths'
: list of integer vectors, representing the paths. 'geodist'
: double vector, the geodesic distances from the source vertex to all vertices in the graph.
Currently no reachability checks are performed, so you have to be sure that the mesh is connected, or at least that the source and target vertices are reachable from one another.
data(humface) p = vcgGeodesicPath(humface,50,c(500,5000)) p$paths[[1]]; # The path 50..500 p$geodist[500]; # Its path length.
data(humface) p = vcgGeodesicPath(humface,50,c(500,5000)) p$paths[[1]]; # The path 50..500 p$geodist[500]; # Its path length.
Compute pseudo-geodesic distance between two points on a mesh
vcgGeodist(x, pt1, pt2)
vcgGeodist(x, pt1, pt2)
x |
triangular mesh of class |
pt1 |
3D coordinate on mesh or index of vertex |
pt2 |
3D coordinate on mesh or index of vertex |
returns the geodesic distance between pt1
and pt2
.
Make sure to have a clean manifold mesh. Note that this computes the length of the pseudo-geodesic path (following the edges) between the two vertices closest to these points.
data(humface) pt1 <- humface.lm[1,] pt2 <- humface.lm[5,] vcgGeodist(humface,pt1,pt2)
data(humface) pt1 <- humface.lm[1,] pt2 <- humface.lm[5,] vcgGeodist(humface,pt1,pt2)
Extract all edges from a mesh and retrieve adjacent faces and vertices
vcgGetEdge(mesh, unique = TRUE)
vcgGetEdge(mesh, unique = TRUE)
mesh |
triangular mesh of class 'mesh3d' |
unique |
logical: if TRUE each edge is only reported once, if FALSE, all occurences are reported. |
returns a dataframe containing:
vert1 |
integer indicating the position of the first vertex belonging to this edge |
vert2 |
integer indicating the position of the second vertex belonging to this edge |
facept |
integer pointing to the (or a, if unique = TRUE) face adjacent to the edge |
border |
integer indicating if the edge is at the border of the mesh. 0 = no border, 1 = border |
require(rgl) data(humface) edges <-vcgGetEdge(humface) ## Not run: ## show first edge lines3d(t(humface$vb[1:3,])[c(edges$vert1[1],edges$vert2[2]),],col=2,lwd=3) shade3d(humface, col=3) ## now find the edge - hint: it is at the neck. ## End(Not run)
require(rgl) data(humface) edges <-vcgGetEdge(humface) ## Not run: ## show first edge lines3d(t(humface$vb[1:3,])[c(edges$vert1[1],edges$vert2[2]),],col=2,lwd=3) shade3d(humface, col=3) ## now find the edge - hint: it is at the neck. ## End(Not run)
Import common mesh file formats and store the results in an object of class "mesh3d" - momentarily only triangular meshes are supported.
vcgImport( file, updateNormals = TRUE, readcolor = FALSE, clean = TRUE, silent = FALSE )
vcgImport( file, updateNormals = TRUE, readcolor = FALSE, clean = TRUE, silent = FALSE )
file |
character: file to be read. |
updateNormals |
logical: if TRUE and the imported file contais faces, vertex normals will be (re)calculated. Otherwise, normals will be a matrix containing zeros. |
readcolor |
if TRUE, vertex colors and texture (face and vertex) coordinates will be processed - if available, otherwise all vertices will be colored white. |
clean |
if TRUE, duplicated and unreferenced vertices as well as duplicate faces are removed (be careful when importing point clouds). |
silent |
logical, if TRUE no console output is issued. |
Object of class "mesh3d"
with:
vb |
4 x n matrix containing n vertices as homolougous coordinates |
it |
3 x m matrix containing vertex indices forming triangular faces |
normals |
4 x n matrix containing vertex normals (homologous coordinates) |
in case the imported files contains face or vertex quality, these will be stored as vectors named $quality (for vertex quality) and $facequality
if the imported file contains vertex colors and readcolor = TRUE, these will be saved in $material$color according to "mesh3d" specifications.
currently only meshes with either color or texture can be processed. If both are present, the function will mark the mesh as non-readable.
Stefan Schlager
data(humface) vcgPlyWrite(humface) readit <- vcgImport("humface.ply")
data(humface) vcgPlyWrite(humface) readit <- vcgImport("humface.ply")
Remove isolated pieces from a surface mesh, selected by a minimum amount of faces or of a diameter below a given threshold. Also the option only to keep the largest piece can be selected or to split a mesh into connected components.
vcgIsolated( mesh, facenum = NULL, diameter = NULL, split = FALSE, keep = 0, silent = FALSE )
vcgIsolated( mesh, facenum = NULL, diameter = NULL, split = FALSE, keep = 0, silent = FALSE )
mesh |
triangular mesh of class "mesh3d". |
facenum |
integer: all connected pieces with less components are removed. If not specified or 0 and diameter is NULL, then only the component with the most faces is kept. |
diameter |
numeric: all connected pieces smaller diameter are removed
removed. |
split |
logical: if TRUE, a list with all connected components (optionally matching requirements facenum/diameter) of the mesh will be returned. |
keep |
integer: if split=T, |
silent |
logical, if TRUE no console output is issued. |
returns the reduced mesh.
Stefan Schlager
## Not run: data(humface) cleanface <- vcgIsolated(humface) ## End(Not run)
## Not run: data(humface) cleanface <- vcgIsolated(humface) ## End(Not run)
Create Isosurface from 3D-array using Marching Cubes algorithm
vcgIsosurface( vol, threshold, from = NULL, to = NULL, spacing = NULL, origin = NULL, direction = NULL, IJK2RAS = diag(c(-1, -1, 1, 1)), as.int = FALSE )
vcgIsosurface( vol, threshold, from = NULL, to = NULL, spacing = NULL, origin = NULL, direction = NULL, IJK2RAS = diag(c(-1, -1, 1, 1)), as.int = FALSE )
vol |
an integer valued 3D-array |
threshold |
threshold for creating the surface |
from |
numeric: the lower threshold of a range (overrides |
to |
numeric: the upper threshold of a range (overrides |
spacing |
numeric 3D-vector: specifies the voxel dimensons in x,y,z direction. |
origin |
numeric 3D-vector: origin of the original data set, will transpose the mesh onto that origin. |
direction |
a 3x3 direction matrix |
IJK2RAS |
4x4 IJK2RAS transformation matrix |
as.int |
logical: if TRUE, the array will be stored as integer (might decrease RAM usage) |
returns a triangular mesh of class "mesh3d"
#this is the example from the package "misc3d" x <- seq(-2,2,len=50) g <- expand.grid(x = x, y = x, z = x) v <- array(g$x^4 + g$y^4 + g$z^4, rep(length(x),3)) storage.mode(v) <- "integer" ## Not run: mesh <- vcgIsosurface(v,threshold=10) require(rgl) wire3d(mesh) ##now smooth it a little bit wire3d(vcgSmooth(mesh,"HC",iteration=3),col=3) ## End(Not run)
#this is the example from the package "misc3d" x <- seq(-2,2,len=50) g <- expand.grid(x = x, y = x, z = x) v <- array(g$x^4 + g$y^4 + g$z^4, rep(length(x),3)) storage.mode(v) <- "integer" ## Not run: mesh <- vcgIsosurface(v,threshold=10) require(rgl) wire3d(mesh) ##now smooth it a little bit wire3d(vcgSmooth(mesh,"HC",iteration=3),col=3) ## End(Not run)
Isotropically remesh a triangular surface mesh
vcgIsotropicRemeshing( x, TargetLen = 1, FeatureAngleDeg = 10, MaxSurfDist = 1, iterations = 3, Adaptive = FALSE, split = TRUE, collapse = TRUE, swap = TRUE, smooth = TRUE, project = TRUE, surfDistCheck = TRUE )
vcgIsotropicRemeshing( x, TargetLen = 1, FeatureAngleDeg = 10, MaxSurfDist = 1, iterations = 3, Adaptive = FALSE, split = TRUE, collapse = TRUE, swap = TRUE, smooth = TRUE, project = TRUE, surfDistCheck = TRUE )
x |
mesh of class |
TargetLen |
numeric: edge length of the target surface |
FeatureAngleDeg |
define Crease angle (in degree). |
MaxSurfDist |
Max. surface distance |
iterations |
ToDo |
Adaptive |
enable adaptive remeshing |
split |
enable refine step |
collapse |
enable collapse step |
swap |
enable dge swap |
smooth |
enable smoothing |
project |
enable reprojection step |
surfDistCheck |
check distance to surface |
returns the remeshed surface mesh
## Not run: data(humface) resampledMesh <- vcgIsotropicRemeshing(humface,TargetLen=2.5) ## End(Not run)
## Not run: data(humface) resampledMesh <- vcgIsotropicRemeshing(humface,TargetLen=2.5) ## End(Not run)
perform kdtree search for 3D-coordinates.
vcgKDtree(target, query, k, nofPoints = 16, maxDepth = 64, threads = 1)
vcgKDtree(target, query, k, nofPoints = 16, maxDepth = 64, threads = 1)
target |
n x 3 matrix with 3D coordinates or mesh of class "mesh3d". These coordinates are to be searched. |
query |
m x 3 matrix with 3D coordinates or mesh of class "mesh3d". We seach the closest coordinates in |
k |
number of neighbours to find |
nofPoints |
integer: number of points per cell in the kd-tree (don't change unless you know what you are doing!) |
maxDepth |
integer: depth of the kd-tree (don't change unless you know what you are doing!) |
threads |
integer: threads to use in closest point search. |
a list with
index |
integer matrices with indeces of closest points |
distances |
corresponding distances |
fast Kmean clustering for 1D, 2D and 3D data
vcgKmeans(x, k = 10, iter.max = 10, getClosest = FALSE, threads = 0)
vcgKmeans(x, k = 10, iter.max = 10, getClosest = FALSE, threads = 0)
x |
matrix containing coordinates or mesh3d |
k |
number of clusters |
iter.max |
maximum number of iterations |
getClosest |
logical: if TRUE the indices of the points closest to the k-centers are sought. |
threads |
integer: number of threads to use |
returns a list containing
centers |
cluster center |
class |
vector with cluster association for each coordinate |
If getClosest=TRUE
selected |
vector with indices of points closest to the centers |
require(Rvcg);require(rgl) data(humface) set.seed(42) clust <- vcgKmeans(humface,k=1000,threads=1)
require(Rvcg);require(rgl) data(humface) set.seed(42) clust <- vcgKmeans(humface,k=1000,threads=1)
calculates the average edge length of a triangular mesh, iterating over all faces.
vcgMeshres(mesh)
vcgMeshres(mesh)
mesh |
triangular mesh stored as object of class "mesh3d" |
res |
average edge length (a.k.a. mesh resolution) |
edgelength |
vector containing lengths for each edge |
Stefan Schlager
data(humface) mres <- vcgMeshres(humface) #histogram of edgelength distribution hist(mres$edgelength) #visualise average edgelength points( mres$res, 1000, pch=20, col=2, cex=2)
data(humface) mres <- vcgMeshres(humface) #histogram of edgelength distribution hist(mres$edgelength) #visualise average edgelength points( mres$res, 1000, pch=20, col=2, cex=2)
Implementation of the command line tool "metro" to evaluate the difference between two triangular meshes.
vcgMetro( mesh1, mesh2, nSamples = 0, nSamplesArea = 0, vertSamp = TRUE, edgeSamp = TRUE, faceSamp = TRUE, unrefVert = FALSE, samplingType = c("SS", "MC", "SD"), searchStruct = c("SGRID", "AABB", "OCTREE", "HGRID"), from = 0, to = 0, colormeshes = FALSE, silent = FALSE )
vcgMetro( mesh1, mesh2, nSamples = 0, nSamplesArea = 0, vertSamp = TRUE, edgeSamp = TRUE, faceSamp = TRUE, unrefVert = FALSE, samplingType = c("SS", "MC", "SD"), searchStruct = c("SGRID", "AABB", "OCTREE", "HGRID"), from = 0, to = 0, colormeshes = FALSE, silent = FALSE )
mesh1 |
triangular mesh (object of class 'mesh3d'). |
mesh2 |
triangular mesh (object of class 'mesh3d'). |
nSamples |
set the required number of samples if 0, this will be set to approx. 10x the face number. |
nSamplesArea |
set the required number of samples per area unit, override nSamples. |
vertSamp |
logical: if FALSE, disable vertex sampling. |
edgeSamp |
logical: if FALSE, disable edge sampling. |
faceSamp |
logical: if FALSE, disable face sampling. |
unrefVert |
logical: if FALSE, ignore unreferred vertices. |
samplingType |
set the face sampling mode. options are: SS (similar triangles sampling), SD (subdivision sampling), MC (montecarlo sampling). |
searchStruct |
set search structures to use. options are: SGIRD (static Uniform Grid), OCTREE, AABB (AxisAligned Bounding Box Tree), HGRID (Hashed Uniform Grid). |
from |
numeric: minimum value for color mapping. |
to |
numeric: maximum value for color mapping. |
colormeshes |
if TRUE, meshes with vertices colored according to distance are returned |
silent |
logical: if TRUE, output to console is suppressed. |
ForwardSampling , BackwardSampling
|
lists containing information about forward (mesh1 to mesh2) and backward (mesh2 to mesh1) sampling with the following entries |
maxdist
maximal Hausdorff distance
meandist
mean Hausdorff distance
RMSdist
RMS of the Hausdorff distances
area
mesh area (of mesh1
in ForwardSampling
and mesh2
in BackwardSampling
)
RMSdist
RMS of the Hausdorff distances
nvbsamples
number of vertices sampled
nsamples
number of samples
distances1 , distances2
|
vectors containing vertex distances from mesh1 to mesh2 and mesh2 to mesh1. |
forward_hist , backward_hist
|
Matrices tracking the sampling results |
if colormeshes == TRUE
mesh1 , mesh2
|
meshes with color coded distances and an additional entry called quality containing the sampled per-vertex distances |
this is a straightforward implementation of the command line tool metro http://vcglib.net/metro.html
P. Cignoni, C. Rocchini and R. Scopigno. Metro: measuring error on simplified surfaces. Computer Graphics Forum, Blackwell Publishers, vol. 17(2), June 1998, pp 167-174
if (requireNamespace("Morpho", quietly = TRUE)) { require(Morpho) data(humface) data(dummyhead) ## align humface to dummyhead.mesh humfalign <- rotmesh.onto(humface,humface.lm,dummyhead.lm) samp <- vcgMetro(humfalign$mesh,dummyhead.mesh,faceSamp=FALSE,edgeSamp=FALSE) ## create heatmap using Morpho's meshDist function } ## Not run: ## create custom heatmaps based on distances mD <- meshDist(humfalign$mesh,distvec=samp$distances1) ## End(Not run)
if (requireNamespace("Morpho", quietly = TRUE)) { require(Morpho) data(humface) data(dummyhead) ## align humface to dummyhead.mesh humfalign <- rotmesh.onto(humface,humface.lm,dummyhead.lm) samp <- vcgMetro(humfalign$mesh,dummyhead.mesh,faceSamp=FALSE,edgeSamp=FALSE) ## create heatmap using Morpho's meshDist function } ## Not run: ## create custom heatmaps based on distances mD <- meshDist(humfalign$mesh,distvec=samp$distances1) ## End(Not run)
Get all non-border edges and both faces adjacent to them.
vcgNonBorderEdge(mesh, silent = FALSE)
vcgNonBorderEdge(mesh, silent = FALSE)
mesh |
triangular mesh of class 'mesh3d |
silent |
logical: suppress output of information about number of border edges |
returns a dataframe containing:
vert1 |
integer indicating the position of the first vertex belonging to this edge |
vert2 |
integer indicating the position of the second vertex belonging to this edge |
border |
integer indicating if the edge is at the border of the mesh. 0 = no border, 1 = border |
face1 |
integer pointing to the first face adjacent to the edge |
face2 |
integer pointing to the first face adjacent to the edge |
data(humface) edges <-vcgNonBorderEdge(humface) ## show first edge (not at the border) ## Not run: require(Morpho) require(rgl) lines3d(t(humface$vb[1:3,])[c(edges$vert1[1],edges$vert2[2]),],col=2,lwd=3) ## plot barycenters of adjacent faces bary <- barycenter(humface) points3d(bary[c(edges$face1[1],edges$face2[1]),]) shade3d(humface, col=3) ## now find the edge - hint: it is at the neck. ## End(Not run)
data(humface) edges <-vcgNonBorderEdge(humface) ## show first edge (not at the border) ## Not run: require(Morpho) require(rgl) lines3d(t(humface$vb[1:3,])[c(edges$vert1[1],edges$vert2[2]),],col=2,lwd=3) ## plot barycenters of adjacent faces bary <- barycenter(humface) points3d(bary[c(edges$face1[1],edges$face2[1]),]) shade3d(humface, col=3) ## now find the edge - hint: it is at the neck. ## End(Not run)
Export meshes to OBJ-files
vcgObjWrite(mesh, filename = dataname, writeNormals = TRUE)
vcgObjWrite(mesh, filename = dataname, writeNormals = TRUE)
mesh |
triangular mesh of class 'mesh3d' or a numeric matrix with 3-columns |
filename |
character: filename (file extension '.obj' will be added automatically. |
writeNormals |
write existing normals to file |
data(humface) vcgObjWrite(humface,filename = "humface") unlink("humface.obj")
data(humface) vcgObjWrite(humface,filename = "humface") unlink("humface.obj")
Export meshes to OFF-files
vcgOffWrite(mesh, filename = dataname)
vcgOffWrite(mesh, filename = dataname)
mesh |
triangular mesh of class 'mesh3d' or a numeric matrix with 3-columns |
filename |
character: filename (file extension '.off' will be added automatically. |
data(humface) vcgOffWrite(humface,filename = "humface") unlink("humface.off")
data(humface) vcgOffWrite(humface,filename = "humface") unlink("humface.off")
Reads Polygon File Format (PLY) files and stores the results in an object of class "mesh3d" - momentarily only triangular meshes are supported.
vcgPlyRead(file, updateNormals = TRUE, clean = TRUE)
vcgPlyRead(file, updateNormals = TRUE, clean = TRUE)
file |
character: file to be read. |
updateNormals |
logical: if TRUE and the imported file contais faces, vertex normals will be (re)calculated. |
clean |
logical: if TRUE, duplicated and unreference vertices will be removed. |
Object of class "mesh3d"
with:
vb |
3 x n matrix containing n vertices as homolougous coordinates |
normals |
3 x n matrix containing vertex normals |
it |
3 x m integer matrix containing vertex indices forming triangular faces |
material$color |
Per vertex colors if specified in the imported file |
from version 0.8 on this is only a wrapper for vcgImport (to avoid API breaking).
Stefan Schlager
Export meshes to PLY-files (binary or ascii)
vcgPlyWrite(mesh, filename, binary = TRUE, ...) ## S3 method for class 'mesh3d' vcgPlyWrite( mesh, filename = dataname, binary = TRUE, addNormals = FALSE, writeCol = TRUE, writeNormals = TRUE, ... ) ## S3 method for class 'matrix' vcgPlyWrite(mesh, filename = dataname, binary = TRUE, addNormals = FALSE, ...)
vcgPlyWrite(mesh, filename, binary = TRUE, ...) ## S3 method for class 'mesh3d' vcgPlyWrite( mesh, filename = dataname, binary = TRUE, addNormals = FALSE, writeCol = TRUE, writeNormals = TRUE, ... ) ## S3 method for class 'matrix' vcgPlyWrite(mesh, filename = dataname, binary = TRUE, addNormals = FALSE, ...)
mesh |
triangular mesh of class 'mesh3d' or a numeric matrix with 3-columns |
filename |
character: filename (file extension '.ply' will be added automatically, if missing. |
binary |
logical: write binary file |
... |
additional arguments, currently not used. |
addNormals |
logical: compute per-vertex normals and add to file |
writeCol |
logical: export existing per-vertex color stored in mesh$material$color |
writeNormals |
write existing normals to file |
data(humface) vcgPlyWrite(humface,filename = "humface") ## remove it unlink("humface.ply")
data(humface) vcgPlyWrite(humface,filename = "humface") ## remove it unlink("humface.ply")
Decimates a mesh by adapting the faces of a mesh either to a target face number, a percentage or an approximate mesh resolution (a.k.a. mean edge length
vcgQEdecim( mesh, tarface = NULL, percent = NULL, edgeLength = NULL, topo = FALSE, quality = TRUE, bound = FALSE, optiplace = FALSE, scaleindi = TRUE, normcheck = FALSE, qweightFactor = 100, qthresh = 0.3, boundweight = 1, normalthr = pi/2, silent = FALSE )
vcgQEdecim( mesh, tarface = NULL, percent = NULL, edgeLength = NULL, topo = FALSE, quality = TRUE, bound = FALSE, optiplace = FALSE, scaleindi = TRUE, normcheck = FALSE, qweightFactor = 100, qthresh = 0.3, boundweight = 1, normalthr = pi/2, silent = FALSE )
mesh |
Triangular mesh of class "mesh3d" |
tarface |
Integer: set number of target faces. |
percent |
Numeric: between 0 and 1. Set amount of reduction relative to existing face number. Overrides tarface argument. |
edgeLength |
Numeric: tries to decimate according to a target mean edge length. Under the assumption of regular triangles, the edges are half as long by dividing the triangle into 4 regular smaller triangles. |
topo |
logical: if TRUE, mesh topology is preserved. |
quality |
logical: if TRUE, vertex quality is considered. |
bound |
logical: if TRUE, mesh boundary is preserved. |
optiplace |
logical: if TRUE, mesh boundary is preserved (may lead to unwanted distortions in some cases). |
scaleindi |
logical: if TRUE, decimatiion is scale independent. |
normcheck |
logical: if TRUE, normal directions are considered. |
qweightFactor |
numeric: >= 1. Quality range is mapped into a squared 01 and than into the 1 - QualityWeightFactor range. |
qthresh |
numeric: Quality threshold for decimation process. |
boundweight |
numeric: Weight assigned to mesh boundaries. |
normalthr |
numeric: threshold for normal check in radians. |
silent |
logical, if TRUE no console output is issued. |
This is basically an adaption of the cli tridecimator from vcglib
Returns a reduced mesh of class mesh3d.
Stefan Schlager
data(humface) ##reduce faces to 50% decimface <- vcgQEdecim(humface, percent=0.5) ## view ## Not run: require(rgl) shade3d(decimface, col=3) ## some light smoothing decimface <- vcgSmooth(decimface,iteration = 1) ## End(Not run)
data(humface) ##reduce faces to 50% decimface <- vcgQEdecim(humface, percent=0.5) ## view ## Not run: require(rgl) shade3d(decimface, col=3) ## some light smoothing decimface <- vcgSmooth(decimface,iteration = 1) ## End(Not run)
check if a mesh is intersected by a set of rays (stored as normals)
vcgRaySearch(x, mesh, mintol = 0, maxtol = 1e+15, mindist = FALSE, threads = 1)
vcgRaySearch(x, mesh, mintol = 0, maxtol = 1e+15, mindist = FALSE, threads = 1)
x |
a triangular mesh of class 'mesh3d' or a list containing vertices and vertex normals (fitting the naming conventions of 'mesh3d'). In the second case x must contain x$vb = 3 x n matrix containing 3D-coordinates and x$normals = 3 x n matrix containing normals associated with x$vb. |
mesh |
triangular mesh to be intersected. |
mintol |
minimum distance to target mesh |
maxtol |
maximum distance to search along ray |
mindist |
search both ways (ray and -ray) and select closest point. |
threads |
number of threads used during search. |
vcgRaySearch
projects a mesh (or set of 3D-coordinates) along a set of given rays (stored as normals) onto a target and return the hit points as well as information if the target mesh was hit at all. If nothing is hit along the ray(within the given thresholds), the ordinary closest point's value will be returned and the corresponding entry in quality
will be zero.
list with following items:
vb |
4 x n matrix containing intersection points |
normals |
4 x n matrix containing homogenous coordinates of normals at intersection points |
quality |
integer vector containing a value for each vertex of |
distance |
numeric vector: distances to intersection |
data(humface) #get normals of landmarks lms <- vcgClost(humface.lm, humface) # offset landmarks along their normals for a negative amount of -5mm lms$vb[1:3,] <- lms$vb[1:3,]+lms$normals[1:3,]*-5 intersect <- vcgRaySearch(lms, humface) ## Not run: require(Morpho) require(rgl) spheres3d(vert2points(lms),radius=0.5,col=3) plotNormals(lms,long=5) spheres3d(vert2points(intersect),col=2) #plot intersections wire3d(humface,col="white")#' ## End(Not run)
data(humface) #get normals of landmarks lms <- vcgClost(humface.lm, humface) # offset landmarks along their normals for a negative amount of -5mm lms$vb[1:3,] <- lms$vb[1:3,]+lms$normals[1:3,]*-5 intersect <- vcgRaySearch(lms, humface) ## Not run: require(Morpho) require(rgl) spheres3d(vert2points(lms),radius=0.5,col=3) plotNormals(lms,long=5) spheres3d(vert2points(intersect),col=2) #plot intersections wire3d(humface,col="white")#' ## End(Not run)
Subsamples surface of a triangular mesh and returns a set of points located on that mesh
vcgSample( mesh, SampleNum = 100, type = c("km", "pd", "mc"), MCsamp = 20, geodes = TRUE, strict = FALSE, iter.max = 100, threads = 0 )
vcgSample( mesh, SampleNum = 100, type = c("km", "pd", "mc"), MCsamp = 20, geodes = TRUE, strict = FALSE, iter.max = 100, threads = 0 )
mesh |
triangular mesh of class 'mesh3d' |
SampleNum |
integer: number of sampled points (see |
type |
character: seclect sampling type ("mc"=MonteCarlo Sampling, "pd"=PoissonDisk Sampling,"km"=kmean clustering) |
MCsamp |
integer: MonteCarlo sample iterations used in PoissonDisk sampling. |
geodes |
logical: maximise geodesic distance between sample points (only for Poisson Disk sampling) |
strict |
logical: if |
iter.max |
integer: maximum iterations to use in k-means clustering. |
threads |
integer number of threads to use for k-means clustering |
Poisson disk subsampling will not generate the exact amount of coordinates specified in SampleNum
, depending on MCsamp
the result will contain more or less coordinates.
sampled points
data(humface) ss <- vcgSample(humface,SampleNum = 500, type="km",threads=1) ## Not run: require(rgl) points3d(ss) ## End(Not run)
data(humface) ss <- vcgSample(humface,SampleNum = 500, type="km",threads=1) ## Not run: require(rgl) points3d(ss) ## End(Not run)
search an existing KD-tree
vcgSearchKDtree(kdtree, query, k, threads = 0)
vcgSearchKDtree(kdtree, query, k, threads = 0)
kdtree |
object of class vcgKDtree |
query |
atrix or triangular mesh containing coordinates |
k |
number of k-closest neighbours to query |
threads |
integer: number of threads to use |
a list with
index |
integer matrices with indeces of closest points |
distances |
corresponding distances |
## Not run: data(humface);data(dummyhead) mytree <- vcgCreateKDtree(humface) ## get indices and distances for 10 closest points. closest <- vcgSearchKDtree(mytree,dummyhead.mesh,k=10,threads=1) ## End(Not run)
## Not run: data(humface);data(dummyhead) mytree <- vcgCreateKDtree(humface) ## get indices and distances for 10 closest points. closest <- vcgSearchKDtree(mytree,dummyhead.mesh,k=10,threads=1) ## End(Not run)
Applies different smoothing algorithms on a triangular mesh.
vcgSmooth( mesh, type = c("taubin", "laplace", "HClaplace", "fujiLaplace", "angWeight", "surfPreserveLaplace"), iteration = 10, lambda = 0.5, mu = -0.53, delta = 0.1 )
vcgSmooth( mesh, type = c("taubin", "laplace", "HClaplace", "fujiLaplace", "angWeight", "surfPreserveLaplace"), iteration = 10, lambda = 0.5, mu = -0.53, delta = 0.1 )
mesh |
triangular mesh stored as object of class "mesh3d". |
type |
character: select smoothing algorithm. Available are "taubin", "laplace", "HClaplace", "fujiLaplace", "angWeight" (and any sensible abbreviations). |
iteration |
integer: number of iterations to run. |
lambda |
numeric: parameter for Taubin smooth (see reference below). |
mu |
numeric:parameter for Taubin smooth (see reference below). |
delta |
numeric: parameter for Scale dependent laplacian smoothing (see reference below).and maximum allowed angle (in radians) for deviation between normals Laplacian (surface preserving). |
The algorithms available are Taubin smoothing, Laplacian smoothing and an improved version of Laplacian smoothing ("HClaplace"). Also available are Scale dependent laplacian smoothing ("fujiLaplace") and Laplacian angle weighted smoothing ("angWeight")
returns an object of class "mesh3d" with:
vb |
4xn matrix containing n vertices as homolougous coordinates. |
normals |
4xn matrix containing vertex normals. |
quality |
vector: containing distances to target. |
it |
4xm matrix containing vertex indices forming triangular faces. |
The additional parameters for taubin smooth are hardcoded to the default values of meshlab, as they appear to be the least distorting
Stefan Schlager
Taubin G. 1995. Curve and surface smoothing without shrinkage. In Computer Vision, 1995. Proceedings., Fifth International Conference on, pages 852 - 857.
Vollmer J., Mencl R. and Mueller H. 1999. Improved Laplacian Smoothing of Noisy Surface Meshes. Computer Graphics Forum, 18(3):131 - 138.
Schroeder, P. and Barr, A. H. (1999). Implicit fairing of irregular meshes using diffusion and curvature flow: 317-324.
data(humface) smoothface <- vcgSmooth(humface) ## view ## Not run: require(rgl) shade3d(smoothface, col=3) ## End(Not run)
data(humface) smoothface <- vcgSmooth(humface) ## view ## Not run: require(rgl) shade3d(smoothface, col=3) ## End(Not run)
Applies implicit smoothing algorithms on a triangular mesh.
vcgSmoothImplicit( mesh, lambda = 0.2, useMassMatrix = TRUE, fixBorder = FALSE, useCotWeight = FALSE, degree = 1L, lapWeight = 1, SmoothQ = FALSE )
vcgSmoothImplicit( mesh, lambda = 0.2, useMassMatrix = TRUE, fixBorder = FALSE, useCotWeight = FALSE, degree = 1L, lapWeight = 1, SmoothQ = FALSE )
mesh |
triangular mesh stored as object of class "mesh3d". |
lambda |
numeric: the amount of smoothness, useful only if
|
useMassMatrix |
logical: whether to use mass matrix to keep the mesh
close to its original position (weighted per area distributed on vertices);
default is |
fixBorder |
logical: whether to fix the border vertices of the mesh;
default is |
useCotWeight |
logical: whether to use cotangent weight; default is
|
degree |
integer: degrees of 'Laplacian'; default is |
lapWeight |
numeric: weight when |
SmoothQ |
logical: whether to smooth the quality (distances to target). |
returns an object of class "mesh3d" with:
vb |
4xn matrix containing n vertices as homolougous coordinates. |
normals |
4xn matrix containing vertex normals. |
it |
4xm matrix containing vertex indices forming triangular faces. |
Zhengjia Wang
data(humface) smoothface <- vcgSmoothImplicit(humface) ## view ## Not run: require(rgl) shade3d(smoothface, col=3) ## End(Not run)
data(humface) smoothface <- vcgSmoothImplicit(humface) ## view ## Not run: require(rgl) shade3d(smoothface, col=3) ## End(Not run)
create platonic objects as triangular meshes
vcgSphere(subdivision = 3, normals = TRUE) vcgSphericalCap(angleRad = pi/2, subdivision = 3, normals = TRUE) vcgTetrahedron(normals = TRUE) vcgDodecahedron(normals = TRUE) vcgOctahedron(normals = TRUE) vcgIcosahedron(normals = TRUE) vcgHexahedron(normals = TRUE) vcgSquare(normals = TRUE) vcgBox(mesh = vcgSphere(), normals = TRUE) vcgCone(r1, r2, h, normals = TRUE)
vcgSphere(subdivision = 3, normals = TRUE) vcgSphericalCap(angleRad = pi/2, subdivision = 3, normals = TRUE) vcgTetrahedron(normals = TRUE) vcgDodecahedron(normals = TRUE) vcgOctahedron(normals = TRUE) vcgIcosahedron(normals = TRUE) vcgHexahedron(normals = TRUE) vcgSquare(normals = TRUE) vcgBox(mesh = vcgSphere(), normals = TRUE) vcgCone(r1, r2, h, normals = TRUE)
subdivision |
subdivision level for sphere (the larger the denser the mesh will be) |
normals |
if TRUE vertex normals are calculated |
angleRad |
angle of the spherical cap |
mesh |
mesh to take the bounding box from |
r1 |
radius1 of the cone |
r2 |
radius2 of the cone |
h |
height of the cone |
Export meshes to STL-files (binary or ascii)
vcgStlWrite(mesh, filename = dataname, binary = FALSE)
vcgStlWrite(mesh, filename = dataname, binary = FALSE)
mesh |
triangular mesh of class 'mesh3d' or a numeric matrix with 3-columns |
filename |
character: filename (file extension '.stl' will be added automatically. |
binary |
logical: write binary file |
data(humface) vcgStlWrite(humface,filename = "humface") unlink("humface.stl")
data(humface) vcgStlWrite(humface,filename = "humface") unlink("humface.stl")
subdivide the triangles of a mesh
vcgSubdivide( x, threshold = NULL, type = c("Butterfly", "Loop"), looptype = c("loop", "regularity", "continuity"), iterations = 3, silent = FALSE )
vcgSubdivide( x, threshold = NULL, type = c("Butterfly", "Loop"), looptype = c("loop", "regularity", "continuity"), iterations = 3, silent = FALSE )
x |
triangular mesh of class "mesh3d" |
threshold |
minimum edge length to subdivide |
type |
character: algorithm used. Options are Butterfly and Loop (see notes) |
looptype |
character: method for type = loop options are "loop","regularity","continuity" (see notes) |
iterations |
integer: number of iterations |
silent |
logical: suppress output. |
returns subdivided mesh
The different algorithms are (from meshlab description):
Butterfly Subdivision: Apply Butterfly Subdivision Surface algorithm. It is an interpolated method, defined on arbitrary triangular meshes. The scheme is known to be C1 but not C2 on regular meshes
Loop Subdivision: Apply Loop's Subdivision Surface algorithm. It is an approximant subdivision method and it works for every triangle and has rules for extraordinary vertices. Options are "loop" a simple subdivision, "regularity" to enhance the meshe's regularity and "continuity" to enhance the mesh's continuity.
data(humface) subdivide <- vcgSubdivide(humface,type="Loop",looptype="regularity")
data(humface) subdivide <- vcgSubdivide(humface,type="Loop",looptype="regularity")
Resample a mesh uniformly
vcgUniformRemesh( x, voxelSize = NULL, offset = 0, discretize = FALSE, multiSample = FALSE, absDist = FALSE, mergeClost = FALSE, silent = FALSE )
vcgUniformRemesh( x, voxelSize = NULL, offset = 0, discretize = FALSE, multiSample = FALSE, absDist = FALSE, mergeClost = FALSE, silent = FALSE )
x |
triangular mesh |
voxelSize |
voxel size for space discretization |
offset |
Offset of the created surface (i.e. distance of the created surface from the original one). |
discretize |
If TRUE, the position of the intersected edge of the marching cube grid is not computed by linear interpolation, but it is placed in fixed middle position. As a consequence the resampled object will look severely aliased by a stairstep appearance. |
multiSample |
If TRUE, the distance field is more accurately compute by multisampling the volume (7 sample for each voxel). Much slower but less artifacts. |
absDist |
If TRUE, an unsigned distance field is computed. In this case you have to choose a not zero Offset and a double surface is built around the original surface, inside and outside. |
mergeClost |
logical: merge close vertices |
silent |
logical: suppress messages |
resampled mesh
## Not run: data(humface) humresample <- vcgUniformRemesh(humface,voxelSize=1,multiSample = TRUE) require(rgl) shade3d(humresample,col=3) ## End(Not run)
## Not run: data(humface) humresample <- vcgUniformRemesh(humface,voxelSize=1,multiSample = TRUE) require(rgl) shade3d(humresample,col=3) ## End(Not run)
update vertex normals of a triangular meshes or point clouds
vcgUpdateNormals(mesh, type = 0, pointcloud = c(10, 0), silent = FALSE)
vcgUpdateNormals(mesh, type = 0, pointcloud = c(10, 0), silent = FALSE)
mesh |
triangular mesh of class 'mesh3d' or a n x 3 matrix containing 3D-coordinates. |
type |
select the method to compute per-vertex normals: 0=area weighted average of surrounding face normals; 1 = angle weighted vertex normals. |
pointcloud |
integer vector of length 2: containing optional parameters for normal calculation of point clouds. The first enty specifies the number of neighbouring points to consider. The second entry specifies the amount of smoothing iterations to be performed. |
silent |
logical, if TRUE no console output is issued. |
mesh with updated/created normals, or in case mesh
is a matrix, a list of class "mesh3d" with
vb |
4 x n matrix containing coordinates (as homologous coordinates |
normals |
4 x n matrix containing normals (as homologous coordinates |
data(humface) humface$normals <- NULL # remove normals humface <- vcgUpdateNormals(humface) ## Not run: pointcloud <- t(humface$vb[1:3,]) #get vertex coordinates pointcloud <- vcgUpdateNormals(pointcloud) require(Morpho) plotNormals(pointcloud)#plot normals ## End(Not run)
data(humface) humface$normals <- NULL # remove normals humface <- vcgUpdateNormals(humface) ## Not run: pointcloud <- t(humface$vb[1:3,]) #get vertex coordinates pointcloud <- vcgUpdateNormals(pointcloud) require(Morpho) plotNormals(pointcloud)#plot normals ## End(Not run)
Compute the k
-ring vertex neighborhood for all query vertex indices vi
. If only a mesh is passed (parameter x
) and the other parameters are left at their default values, this compute the adjacency list representation of the mesh.
vcgVertexNeighbors(x, vi = NULL, numstep = 1L, include_self = FALSE)
vcgVertexNeighbors(x, vi = NULL, numstep = 1L, include_self = FALSE)
x |
tmesh3d instance from the |
vi |
optional, vector of positive vertex indices for which to compute the neighborhoods. All vertices are used if left at the default value |
numstep |
positive integer, the number of times to extend the neighborhood from the source vertices (the |
include_self |
logical, whether the returned neighborhood for a vertex |
list of positive integer vectors, the neighborhoods.
data(humface) adjacency_list <- vcgVertexNeighbors(humface) v500_5ring = vcgVertexNeighbors(humface, vi=c(500), numstep = 5)
data(humface) adjacency_list <- vcgVertexNeighbors(humface) v500_5ring = vcgVertexNeighbors(humface, vi=c(500), numstep = 5)
find all faces belonging to each vertex in a mesh and report their indices
vcgVFadj(mesh)
vcgVFadj(mesh)
mesh |
triangular mesh of class "mesh3d" |
list containing one vector per vertex containgin the indices of the adjacent faces
Compute volume for manifold meshes
vcgVolume(x)
vcgVolume(x)
x |
triangular mesh of class mesh3d |
returns volume
Please note, that this function only works reliably on watertight, coherently oriented meshes that constitute a manifold. In case your mesh has some issues regarding non-manifoldness or there are isolated pieces flying around, you can use vcgIsolated and vcgClean to remove those.
mysphere <- vcgSphere() vcgVolume(mysphere) ## Not run: ## here is an example where the mesh has some non-manifold vertices mysphere <- vcgSphere(normals=FALSE) ## add a degenerate face mysphere$it <- cbind(mysphere$it,c(1,2,1)) try(vcgVolume(mysphere)) ## fix the error using vcgClean(): vcgVolume(vcgClean(mysphere,sel=0:6,iterate=TRUE)) ## End(Not run)
mysphere <- vcgSphere() vcgVolume(mysphere) ## Not run: ## here is an example where the mesh has some non-manifold vertices mysphere <- vcgSphere(normals=FALSE) ## add a degenerate face mysphere$it <- cbind(mysphere$it,c(1,2,1)) try(vcgVolume(mysphere)) ## fix the error using vcgClean(): vcgVolume(vcgClean(mysphere,sel=0:6,iterate=TRUE)) ## End(Not run)
Export meshes to WRL-files
vcgWrlWrite(mesh, filename = dataname, writeCol = TRUE, writeNormals = TRUE)
vcgWrlWrite(mesh, filename = dataname, writeCol = TRUE, writeNormals = TRUE)
mesh |
triangular mesh of class 'mesh3d' or a numeric matrix with 3-columns |
filename |
character: filename (file extension '.wrl' will be added automatically. |
writeCol |
logical: export existing per-vertex color stored in mesh$material$color |
writeNormals |
write existing normals to file |
data(humface) vcgWrlWrite(humface,filename = "humface") unlink("humface.wrl")
data(humface) vcgWrlWrite(humface,filename = "humface") unlink("humface.wrl")