Geometry
The geometry of a surface mesh is encoded in its Attributes. This package provides an interface for accessing attributes like vertex position, face normals and vertex normals. Various algorithms use this interface for querying geometry.
Basic representation
AbstractGeometry
provides a common interface to attributes, which encode geometry of a surface mesh. The main idea is collecting the Mesh
and geometry Attributes into some custom subtype such that availability of attributes can be decided "statically" from this type and algorithms that query geometry may dispatch on availability (e.g., whether to read normal vectors stored in an attribute or compute them on-the-fly).
Use geometry
to create a geometry object with (a subset of) vertex positions, face normals and vertex normals. Use triangulated
to "tag" a geometry object as a triangle mesh, i.e., the user asserts that all faces are (and will be) triangles. Note that this assertion is not checked.
Example
using PMesh
using StaticArrays
mesh = createmesh((:x => SVector{3, Float64}, :nrm => SVector{3, Float64}),
(:nrm => SVector{3, Float64}), (), ())
# ...
gp = geometry(mesh, Val(:x)) # vertex positions only
μlen = meanedgelength(g)
gpnn = geometry(mesh, postition=Val(:x), fnormal=Val(:nrm), vnormal=Val(:nrm))
updatenormals!(gpnn, weights=Val(:area))
for v in vertices(mesh)
@show normal(g, v) # compute vertex normal on the fly
@show normal(gpnn, v) # provide stored ("updated") normal
end
gtpnn = triangulated(gpnn) # tag as triangle mesh
updatenormals!(gtpnn, weights=Val(:area))
# calls trianglenormal instead of facenormal
Some methods are available as method(g, ...)
for g=geoemtry(mesh, x, ...)
and attribute x
as well as method(mesh, x, ...)
.
Internals
The internal implementation uses a common interface to access attributes and type traits to encode availability of attributes and the assertion of a triangulation. (See _has_vx(), ..., _istriangulation()
, which return types.)
Any custom defined subtype of AbstractGeometry
that supports this interface can be used with any method/algorithm that takes an geometry
instance.
Reference
PMesh.has_fnormal
— FunctionPMesh.has_position
— FunctionPMesh.has_vnormal
— FunctionPMesh.istriangulation
— Functionistriangulation(geometry(...))
Query type traits: Is triangulation or general polygonal mesh?
This method relies only on static type information. It does not check face degree
s!
See also geometry
PMesh.AbstractGeometry
— TypeBind a mesh to an object that additionally provides attributes with vertex positions and optionally face normals and vertex normals for evaluation and storage of geometry queries.
See also geometry
, triangulated
PMesh.geometry
— Functiong = geometry(mesh, x[, fnrm, vnrm])
g = geometry(mesh, x[, fnormal=fnrm, vnormal=nvrm])
g = geometry(mesh, x, nothing, vnormal=nvrm]
g = geometry(mesh, position=x[, fnormal=fnrm, vnormal=nvrm])
Bind mesh
in an AbstractGeometry
instance that provides attributes with vertex positions and optionally face normals and vertex normals. x
, fnrm
, and vnrm
denote attributes of the mesh.
Use the returned g
for computing geometric quantities which can be stored in and queried from attributes if available.
Here is a typical pattern: facenormal
computes a normal vector, normal
queries a normal vector if an attribute is available, otherwise falls back to computing facenormal
. Finally, updatenormal!
stores a computed normal in an attribute if available (no effect else).
You can pass attributes (inexpensive) or attribute names. The latter is efficient and type-stable if you provide types such as Val(:x)
(instead of :x
), see also Access and performance
geometry
always requires vertex positions. You can but don't have to provide face normals and/or vertex normals.- Providing or not providingg an attribute effects the use of stored data versus computing data on demand.
The current implememtation requires the same element types for positions and normals, e.g., SVector{3, Float64}
.
See also triangulated
PMesh.triangulated
— Functiong = triangulated(geometry(mesh, ...))
Tags g
as a triagulated surface, i.e., all faces of mesh
are assumed to be triangles. For instance, calling facearea
becomes equivalent to calling trianglearea
.
Requires that mesh
is and remains a triangle mesh. This is a contract with the user, which is not checked!
See also geometry
Access mesh and attributes
Base.position
— Functiong = geometry(mesh, x)
p = position(g) # get attribute
pv = position(g, v) # get attribte element
@assert p[v] == pv
Get position of vertex v
given as vertex attribute x
.
PMesh.fnormal
— Methodg = geometry(mesh, ..., fnormal=...)
fnrm = fnormal(g)
Get face normals attribute.
PMesh.fnormals
— Methodg = geometry(mesh, ..., fnormal=...)
fn = fnormals(g)
nrm = fn(f) # map FHnd to face normal
Create a functor that maps f::FHnd
to stored fnormal
.
See also geometry
, maphandles
, fnormal
PMesh.positions
— Functionfp = positions(mesh, x)
g = geometry(mesh, x)
fp = positions(g)
x = fp(v) # map VHnd to position
Create a functor that maps v::VHnd
to stored position
.
xs = positions(g, f)
Get positions of vertices of face f
as 3-tuple for a triangulation or vector.
See also geometry
, maphandles
, position
PMesh.tessellation
— MethodPMesh.vnormal
— Methodg = geometry(mesh, ..., vnormal=...)
vnrm = vnormal(g)
Get vertex normals attribute.
PMesh.vnormals
— Methodg = geometry(mesh, ..., fnormal=...)
vn = vnormals(g)
nrm = vn(v) # map VHnd to vertex normal
Create a functor that maps v::VHnd
to stored vnormal
.
See also geometry
, maphandles
, vnormal
Basic queries
- Compute
edgevector
, their lengths, etc. - Compute (non-orthogonal) bases for face planes
Base.diff
— Methodg = geometry(mesh, x)
d = diff(g, ...)
Synonym for d = edgevector(g, ...)
.
See also geometry
, edgevector
PMesh.barycenter
— Functionc = barycenter(mesh, x, ...)
g = geometry(mesh, x)
c = barycenter(g, h)
c = barycenter(g, e)
c = barycenter(g, (v0, v1))
c = barycenter(g, f)
c = barycenter(g)
Get barycenter of edge or face or all vertex positions.
See also geometry
PMesh.circumference
— Functionlen = circumference(mesh, x, f)
m = geometry(mesh, x)
len = circumferemce(g, f)
Compute circumference of face f
See also geometry
, edgevectors
PMesh.edgelength
— Functionlen = edgelength(mesh, x, e)
len = edgelength(mesh, x, h)
g = geometry(mesh, x)
len = edgelength(g, e)
len = edgelength(g, h)
Compute length of an edge.
See also geometry
, edgevector
PMesh.edgelengths
— Functionlens = collect(edgelengths(mesh, x))
g = geometry(mesh, x)
for len in edgelengths(g)
# ...
end
Generate all edge lengths.
Examples:
μlen = Statistics.mean(edgelengths(g))
minlen, maxlex = extrema(edgelengths(g))
See also geometry
, edgelength
, meanedgelength
PMesh.edgevector
— Functiond = edgevector(mesh, x, h)
d = edgevector(mesh, x, e)
d = edgevector(mesh, x, (v0, v1))
g = geometry(mesh, x)
d = edgevector(g, h)
d = edgevector(g, e)
d = edgevector(g, (v0, v1)
Compute vector spanning an edge.
See also geometry
, edgelength
PMesh.edgevectors
— Functionds = collect(edgevectors(mesh, x, f)
g = geometry(mesh, x)
for d in edgevectors(g, f)
# ...
end
Compute vectors spanning edges of fac f
See also geometry
, edgevector
, circumference
PMesh.meanedgelength
— Functionμlen = meanedgelength(mesh, x)
g = meaure(mesh, x)
μlen = meanedgelength(g)
Compute mean edge length.barycenter(m::Mesh, ax, hs) = barycenter(geometry(mesh, ax), hs)
See also See also geometry
, edgelengths
PMesh.rawfacebase
— Functionx0, u, v = rawfacebase(mesh, x, f)
g = geometry(mesh, x)
x0, u, v = rawfacebase(g, f)
Compute edge vectors that span a planar face with origin shifted to x0
. u, v
provide a basis for the face's plane but they are neither unit length nor orthogonal.
!! warning The caller must not make assumptions on the particular choice of vectors!
The current implementation uses an arbitrary pair of edgevectors
and disregards conditioning/angles.
See also geometry
, rawtrianglebase
, facenormal
, trianglenormal
, tangentspace
](@ref)
PMesh.rawtrianglebase
— Functionx0, u, v = rawtrianglebase(mesh, x, f)
g = geometry(mesh, x)
x0, u, v = rawtrianglebase(g, f)
Compute edge vectors that span a planar triangle with origin shifted to x0
. u, v
provide a basis for the triangle's plane but they are neither unit length nor orthogonal.
See also geometry
, rawfacebase
, rawtrianglenormal
, tangentspace
Areas
- Compute areas of faces
trianglearea
,facearea
- Compute
surfacearea
scalesurfacearea!
scales a surface withcentroid
as origin to match a prescribed area
PMesh.facearea
— Functionar = facearea(mesh, x, f)
g = geometry(mesh, x)
ar = facearea(g, f)
Compute area of face f
.
See also geometry
, trianglearea
, surfacearea
PMesh.facenormalarea
— FunctionSame as trianglenormalarea
but more geneal.
There is no benefit for non-triangular faces. The purpose of this method is providing a unifornm interface.
See also geometry
, trianglenormalarea
PMesh.scalesurfacearea!
— Methodg = geometry(mesh, x)
scalesurfacearea!(g[; area=4π, maxiter=4, rtol=0.001])
scalesurfacearea!(m, x[; area=4π, maxiter=4, rtol=0.001])
@assert surfacearea(g) ≈ area
Scale mesh with centroid
as origin such that surfacearea
becomes area
.
The optional arguments specify the maximum number of iterations and a relative tolerance.
See also geometry
, surfacearea
PMesh.surfacearea
— Functionar = surfacearea(mesh, x, f)
g = geometry(mesh, x)
ar = surfacearea(g, f)
Compute area of face f
.
PMesh.trianglearea
— Functionar = trianglearea(mesh, x, f)
g = geometry(mesh, x)
ar = trianglearea(g, x)
Compute area of triangle f
.
ar = trianglearea(u, v)
Compute area of triangle spanned by vectors u
and v
in 2d or 3d.
See also geometry
, rawtrianglebase
, trianglenormal
Angles
angle
,anglenext
,angleprev
compute angles closed by edges of a facedihedralangle
,computedihedralangle
compute the dihedral angle between two faces sharing an edge
Base.angle
— Functiong = geometry(mesh, position=x)
α = angle(g, h::Hnd)
α = angle(mesh, x, h::Hnd)
Synonym for anglenext
.
Base.angle
— Methodα = angle(v1, v2)
Compute angle (rad) enclosed by vectors v1
and v2
.
PMesh.anglenext
— Methodg = geometry(mesh, position=x)
α = anglenext(g, h::Hnd)
α = anglenext(mesh, x, h::Hnd)
Compute angle α
enclosed by half-edge h
and its successor next(mesh, h)
within their common face.
PMesh.angleprev
— Methodg = geometry(mesh, position=x)
α = angleprev(g, h::Hnd)
α = angleprev(mesh, x, h::Hnd)
Compute angle α
enclosed by half-edge h
and its predecessor pred(mesh, h)
within their common face.
PMesh.computedihedralangle
— Methodg = geometry(mesh, position=x, ...)
α = computedihedralangle(g, h::HHnd)
α = computedihedralangle(g, e::EHnd)
Same as dihedralangle
but enforce computation of face normals, i.e., do not lookup face attribute if present.
See also dihedralangle
PMesh.dihedralangle
— Functiong = geometry(mesh, position=x, ...)
α = dihedralangle(g, h::HHnd)
α = dihedralangle(g, e::EHnd)
α = dihedralangle(mesh, x, ...)
Compute dihedral angle enclosed by the faces that share half-edge h
(or edge e
).
The face normals of g::AbstractGeometry
are determined by normal
, i.e., the normal vector is read from the face attribute if present (and computed otherwise).
See also angle
, computedihedralangle
Statistics
- Compute
centroid
of a single face and area-weighted centroid of the surface mesh - Define
Statistics.mean
,Statistics.std
, andStatistics.cov
forAbstractGeometry
objects to compute mean, standard deviation and covariance of vertex positions - Center surface mesh by shifting
mean
orcentroid
usingcentermean!
orcentercentroid!
pca
computes the principal component analysis of vertex positions or facecentroid
sextrema
computes an axis aligned bounding box spanned by the returned vectorsradius
computes the radius w.r.t to a prescribed center (or themean
)
The centroid
of a surface mesh is an area-weighted average of face centroids, i.e., isolated vertices are ignored. In contrast, Statistics.mean
(and std
and cov
) consider plain vertex positions and include isolated vertices.
Base.extrema
— Methodg = geometry(mesh, x)
xmin, xmax = extrema(g)
Compute axis aligned bounding box spanned by extrema xmin, xmax
.
See also radius
PMesh.centercentroid!
— Methodg = geometry(mesh, x)
centercentroid!(g)
@assert centroid(g) ≈ zeros(eltype(position(g)))
centercentroid!(g; origin = z)
@assert centroid(g) ≈ z
Shift vertex positions such that [centroid
]@(ref) is located in the origin
.
See also [centroid
]@(ref), centermean!
PMesh.centermean!
— Methodg = geometry(mesh, x)
centermean!(g)
@assert Statistics.mean(g) ≈ zeros(eltype(position(g)))
centermean!(g; origin = z)
@assert Statistics.mean(g) ≈ z
Shift vertex positions such that Statistics.mean
is located in the origin
.
See also [mean
]@(ref), centercentroid!
PMesh.centroid
— Functiong = geometry(mesh, x)
c = centroid(g, f)
c = centroid(mesh, x, f)
c = centroid(g)
c = centroid(mesh, x)
Compute centroid of a single face f
or of the entire surface g
using area-weighted centroids of faces.
PMesh.pca
— Methodg = geometry(mesh, x)
U, λ, center = pca(mesh[; centroids=false])
Compute Principle Component Analysis of either position(g)
or centroid
s (centroids=true
).
Compute the covariance matrix C = (X .- center)'*(X .- center) / n
, such that
C = U * Diagonal(λ) * U'
with nonnegative eigenvalues λ
. Here, n
equals the number of vertices or the number of faces (centroids=true
), and center
equals the mean
or the centroid
(centroids=true
).
The additionally returned center
is the mean
or centroid
.
See also centroid
PMesh.radius
— Methodg = geometry(mesh, x)
r = radius(mesh[; center=mean(g)])
Compute (maximum) radius of surface mesh
w.r.t center
.
See also Base.extrema
, Statistics.mean
Normals
- Compute unit-length face normals:
trianglenormal
,facenormal
- Compute vertex normals from weighted face normals:
vertexnormal
- Handle read/write face and vertex normals from/to their attributes if available: read instead of compute (
normal
), explicitlyupdatenormal!
. - Update all normals:
updatefacenormals!
,updatevertexnormals!
,updatenormals!
PMesh.facenormal
— Functionnrm = facenormal(mesh, x, f)
g = geometry(mesh, x)
nrm = facenormal(g, f)
Compute face normal from rawfacebase
.
See also geometry
, rawfacebase
, trianglenormal
PMesh.rawtrianglenormal
— Functionnrm = rawtrianglenormal(mesh, x, f)
g = geometry(mesh, x)
nrm = rawtrianglenormal(g, x)
Compute "unnormalized" normal vector to triangle f
.
See also geometry
, rawtrianglebase
, trianglenormal
PMesh.trianglenormal
— Functionnrm = trianglenormal(mesh, x, f)
g = geometry(mesh, x)
nrm = trianglenormal(g, f)
Compute unit-length normal vector to triangle f
.
The trianglenormal
always returns a finite vector: For a degenerated triangle with zero area, its normal vector is undefined, and trianglenormal
returns zeros(..)
(*not= [NaN,...]
).
See also geometry
, rawtrianglenormal
PMesh.trianglenormalarea
— Functionnrm, ar = trianglenormalarea(mesh, x, f)
g = geometry(mesh, x)
nrm, ar = trianglenormalarea(g, f)
Compute unit-length trianglenormal
and trianglearea
.
See also geometry
, trianglenormal
and trianglearea
PMesh.normal
— Functiong = geometry(mesh, x, ...)
nrm = normal(g, f)
nrm = normal(g, v)
Get or compute face or vertex normal depending on availability of attribute.
See also geometry
, updatenormal!
, facenormal
, vertexnormal
PMesh.updatefacenormals!
— Functiong = geometry(mesh, x, ...)
updatefacenormals!(g)
Recompute and update all face normals.
See also geometry
, updatenormal!
, facenormal
, updatevertexnormals!
, updatenormals!
PMesh.updatenormal!
— Functiong = geometry(mesh, x, ...)
updatenormal!(g, f[, nrm])
updatenormal!(f, v[, nrm])
Update face normal or vertex normal from nrm
or [facenormal
] or vertexnormal
. No effect if g
has no attribute for face normals.
See also geometry
, trianglenormal
, facenormal
, normal
PMesh.updatenormals!
— Functiong = geometry(mesh, x, ...)
updatenormals!(g[, weights=AngleWeights()])
Short for
updatefacenormals!(g)
updatevertexnormals!(g[, weights=AngleWeights()])
PMesh.updatevertexnormals!
— Functiong = geometry(mesh, x, ...)
updatevertexnormals!(g[, weights=AngleWeights()])
Recompute and update vertex normals using weights
(see vertexnormal
.
This method does not update vertex normals. If g
provides an attribute for vertex normals, they are assumed to be up to date and may be used w/o recomputation!
See also geometry
, updatenormal!
, vertexnormal
, updatevertexnormals!
, updatenormals!
PMesh.vertexnormal
— Functionnrm = vertexnormal(mesh, x, v[, weights=AngleWeights()])
g = geometry(mesh, x)
nrm = vertexnormal(g, v[, weights=AngleWeights()])
Compute vertex normal from weighted sum of facenormal
s.
This method uses normal
to query face normals: If face normals are provided as an attribute, use their values assuming they are up to date! Otherwise, compute them explicitly.
Possible values of for weights
are
AngleWeights()
(default)AreaWeights()
InverseAreaWeights()
UniformWeights()
See geometry
, normal
, `updatenormal!, facenormal
, updatevertexnormals!
Barycentric coordinates and interpolation
Get simplex (line segment or triangle) or linear interpolant, likewise:
f = baryint(a, b) # line segment
f = baryint(a, b, c) # triangle
f = baryint(mesh, x, h) # half-edge
f = baryint(mesh, x, f) # triangle
y = f(λ) # barycentric interpolation
Get affine map that computes barycentric coordinates w.r.t. simplex:
b = mapbarycoords(f) # get map as callable object
λ = b(x) # compute barycentric coordinates w.r.t. f
λ = barycoords(f, x) # short cut for single computation
PMesh.BarycentricCoordinates.baryint
— Functionb = baryint(a, b)
b = baryint(a, b, c)
b = baryint((a, b, c))
x = vattr(mesh, Val(:x))
b = baryint(x, v1, v2) # vertex handles
b = baryint(x, v1, v2, v3)
b = baryint(Val(:x), ...)
b = baryint(mesh, x, h) # vertices of half-edge h
b = baryint(mesh, x, f) # vertices of triangular face f
g = triangulated(geometry(mesh, ...))
b = bayrint(g, x, h)
b = bayrint(g, x, f)
Creates a callable object b
for barycentric interpolation between 2 (e.g., algong half-edge) or 3 (e.g., in triangle) values, typically given by a vertex attribute.
b
can be evaluated at barycentric coordinates λ
with sum(λ) == 1
λ = SVector(λ₁, λ₂)
y = b(λ)
y = b(λ₁, λ₂)
λ = SVector(λ₁, λ₂, λ₃)
y = b(λ)
y = b(λ₁, λ₂, λ₃)
where λ
is an SVector
or NTuple
of length 2 or 3.
See also geometry
PMesh.BarycentricCoordinates.barycoords
— Methodf = baryint(x1, x2, ...) :: Simplex λ = barycoords(f, x)
@assert sum(λ) ≈ 1
Evaluate barycentric coordinates λ
for x
in f
.
mapbarycoords
constructs an affine map from x
to λ
. This is more efficient, if multiple x
are mapped for the same simplex.
See also mapbarycoords
PMesh.BarycentricCoordinates.mapbarycoords
— Functionb = mapbarycoords(x1, x2, ...)
f = baryint(x1, x2, ...) :: Simplex b = mapbarycoords(f) :: MapBaryCoords λ = b(x) y = f(λ)
λ = barycoords(f, x)
Get affine map MapBaryCoords
that maps x
to barycentric coordinates λ
.
For computing barycentric coordinates for one single point x
, barycoords
may be more convenient.
In case of a n
-simplex f
in m>n
dimensions, the map "assumes" that is x
is in the affine subspace defined by f
. If this is not the case, there exist no barycentric coordinates, and the map computes barycentric coordinates for a projection into the affine subspace. Note that in this case, we have f(λ) != x
.
In case the simplex f
is degenerated, there exist no barycentric coordinates. In this case, the map computes a "best fit":
Mapping x
to barycentric coordinates λ
never fails, i.e., yields always finite barycentric coordinates λ
. Consequently, one may have f(λ) != x
.
See MapBaryCoords
, barycoords
, baryint
PMesh.BarycentricCoordinates.MapBaryCoords
— TypeAffine map that computes barycentric coordinates.
b = mapbarycoords(x1, x2, ...)
f = baryint(x1, x2, ...) :: Simplex b = mapbarycoords(f) :: MapBaryCoords λ = b(x) y = f(λ)
λ = barycoords(f, x)
@assert sum(λ) ≈ 1
In case of a n
-simplex f
in m>n
dimensions, the map "assumes" that is x
is in the affine subspace defined by f
. If this is not the case, there exist no barycentric coordinates, and the map computes barycentric coordinates for a projection into the affine subspace. Note that in this case, we have f(λ) != x
.
In case the simplex f
is degenerated, there exist no barycentric coordinates. In this case, the map computes a "best fit":
Mapping x
to barycentric coordinates λ
never fails, i.e., yields always finite barycentric coordinates λ
. Consequently, one may have f(λ) != x
.
See also mapbarycoords
, baryint
, Simplex
, barycoords
PMesh.BarycentricCoordinates.Simplex
— TypeSimplex for barycentric interpolation as a callable object.
f = baryint(...) :: Simplex
y = f(λ)
y = f(λ₁, λ₂, ...)
See also baryint
, MapBaryCoords
, mapbarycoords
Transformations
- Represent
AffineMap
s andLinearMap
(possibly followed by nonlinear normalization) as callable objects - Apply transformations on positions, tangent vectors and normal vectors.
mapvattr!
applies transformations on a vertex attribute for all vertices
fx = affinemap(A, b)
fv = linearmap(fx) # == mapvector(fx)
fuv = mapunitvector(fx) # == mapunitvector(fv)
fn = mapnormal(fx)
fx'
inv(fv)
vtolocal = fx ∘ positions(geometry(mesh, ...))
# aply transformation f to all vertices for vertex attribute x
mapvattr!(f, mesh, x)
PMesh.affinemap
— Methodf = affinemap(A, b)
f(x) == A * x - b
affinemap(f) == f
Construct AffineMap
PMesh.linearmap
— Functionf = linearmap(A, b)
f = linearmap(g::AffineMap)
linearmap(f) == f
f(x) == A * x
Construct LinearMap
PMesh.mapnormal
— Functionf = affinemap(A, b)
g = mapnormal(f)
g = mapnormal(linearmap(f))
Nonlinear map that maps normal vectors n
with inv(A)' * n
followed by normalization.
See also LinearMap
, mapvector
, mapunitvector
PMesh.maptangentunitvector
— FunctionSame as mapunitvector
but use only dimensions 1:2
.
Use case: tolocal
maps to local u,v,w
coordinates, map u, v
to 3d ("normalized") unit vector in tangent space.
See also mapvector
PMesh.maptangentvector
— FunctionSame as mapvector
but use only dimensions 1:2
.
Use case: tolocal
maps to local u,v,w
coordinates, map u, v
to 3d vector in tangent space.
See also mapvector
PMesh.mapunitvector
— Functionf = affinemap(A, b)
g = mapunitvector(f)
g = mapunitvector(linearmap(f))
Nonlinear map that composes linearmap
and normalization of result.
PMesh.mapvattr!
— Methodt = affinemap(A, b)
mapvattr(t, mesh, Val(:x))
mapvattr(mapnormal(t), mesh, Val(:n))
mapvattr(x -> 2x, mesh, val(:myattr))
Apply transformation on vertex attribute for all vertices.
See also vattr
, affinemap
, linearmap
, mapunitvector
, normal
, maptangentvector
, maptangentunitvector
PMesh.mapvector
— Functionf = affinemap(A, b)
g = linearmap(f)
g = mapvector(f)
Synomym for linearmap
See also LinearMap
, mapunitvector
, mapnormal
PMesh.AffineMap
— TypeAffine map f(x) = A * x - b
f = affinemap(A, b)
f(x) == A * x - b
g = linearmap(f) # linear part A * x
g = mapvector(f) # same as linearmap(f)
g = mapunitvector(f) # plus normalization
g = mapnormal(f) # map unit normals inv(A)' * x
PMesh.LinearMap
— TypeLinear map f(x) = A * x
f = linearmap(A)
f(x) == A * x
linearmap(f) == f
g = f' # adjoint
g = inv(f) # inverse for maps (for quadratic A)
g = inv(f')
g = mapunitvector(f) # plus normalization
g = mapnormal(f) # map unit normals inv(A)' * x
Maps on triangles
Provide (linear) maps from a standard triangle to a triangle that is congruent to a given triangle.
A = maplinear2d(g, f)
fA = maplinear2d(g)
A = fA(f)
PMesh.linearmap2d
— Functiong = triangulated(geometry(mesh, x, ...))
A = linearmap2d(g, f)
fA = linearmap2d(g)
A = fa(f)
Construct linear map A
as the 2×2 matrix that maps the standard triangle [0,0], [1,0], [0,1]
to a 2d triangle that is congruent (with same orientation) to f
.
See also localframe
Tangent spaces
- Provide orthonormal basis for tangent space (face, vertex):
tangentspace
- Project point into tangent space:
projectintotangentspace
The methods use normal
, i.e., read from attributes if available.
PMesh.localframe
— Functiong = geometry(mesh, x, ...)
x0, u, v = localframe(g, f)
Compute orthonormal basis u, v, n
such that u, v
span plane of face f
and n
is the normal vector. The origin is shifted to x0
.
x0, u, v, n = localframe(g, v)
Compute orthonormal basis u, v, n
such that u, v
span tangent plane defined by vertex normal n
. Here, x0
is the position of vertex v
.
This method uses stored normal
if attribute is available.
This method is similar to tangentspace
but returns also the normal vector.
See also geometry
, rawfacebase
, normal
, tangentspace
, tolocal
PMesh.projectintotangentspace
— Functiong = geometry(mesh, x, ...)
xhat = projectintotangentspace(g, f, x)
xhat = projectintotangentspace(g, v, x)
Project x
into tangentspace
of face f
or vertex v
. The projection xhat
is a 3d vector.
For a planar face f
, its tangent space is the respective plane. For a vertex v
, its tangent space is defined implicitly by the normal
vector.
This method uses the stored normal
if attribute is available. It does require a basis, i.e., does not call tangentspace
.
See also geometry
, normal
, tangentspace
PMesh.tangentspace
— Methodg = geometry(mesh, x, ...)
x0, u, v = tangentspace(g, f)
Compute orthonormal basis u, v, for plane spanning face f
with origin shifted to x0
.
x0, u, v = tangentspace(g, v)
Compute orthonormal basis for tangent plane defined by vertex normal. Here, x0
is the position of vertex v
.
This method uses stored normal
if attribute is available.
This method is similar to localframe
but does not return normal vector.
See also geometry
, rawfacebase
, normal
, localframe
, projectintotangentspace
, tolocal2d
PMesh.tolocal
— Functiong = geometry(mesh, x, ...)
u = tolocal(g, f, z)
u = tolocal(g, v, z)
m = tolocal(g, f) # Provide callable object
m = tolocal(g, v)
u = m(z)
invm = toworld(m)
z = invm(u)
Transform 3d point x
to 3d local coordinates u
defined by localframe
of face f
or vertex v
.
The local coordinate system for a vertex v
has x[v]
as origin, i.e., x[v]
is mapped to the null vector.
Combine with positions
to map, e.g., a 1-ring:
us = map(tolocal(g, v) ∘ positions(g), ring(mesh, v))
See also localframe
, toworld
, tolocal2d
, positions
PMesh.tolocal2d
— Functiong = geometry(mesh, x, ...)
u = tolocal2d(g, f, z)
u = tolocal2d(g, v, z)
f = tolocal2d(g, f) # Provide callable object
f = tolocal2d(g, v)
u = f(x)
Transform point x
2d to local coordinates u
defined by tangentspace
of face f
or vertex v
.
The local coordinate system for a vertex v
has x[v]
as origin, i.e., x[v]
is mapped to the null vector.
Combine with positions
to map, e.g., a 1-ring:
us = map(tolocal2d(g, v) ∘ positions(g), ring(mesh, v))
See also localframe
, tolocal
, positions
PMesh.toworld
— Methodg = geometry(mesh, x, ...)
m = tolocal(g, ...)
im = toworld(m)
u = m(z)
z = im(u)
Inverse map to affine map obtained from tolocal
or tolocal2d
.
Requires rigid transformation as obtained from tolocal
, i.e., A'A == I
.
See also localframe
, tolocal
, tolocal2d
Distance queries and projection
- Setup kd-tree from for
pointlocation
tests - Query
nearest
face and vertex and projection point in barycentric coordinates - Project onto geometry
The PointLocation
context uses kd-trees from DistanceQueries
.
Handles returned by nearestvertex
and nearestface
become invalid if the mesh changes after calling pointlocation
.
Consider using compact
, i.e., ensure iscontiguous
before calling pointlocation
!
Sequences of connected edges
Treat sequences of connected edges as discrete curves on meshes.
arclength
measures length (no checks)ishalfedgecurve
tests for connected (possibly closed) sequencearclengthknots
,uniformknots
knot sequences for parametrizationpositionhalfedgecurve!
specify vertex positions by parametric curve (e.g., layout boundary on circle).
PMesh.arclength
— Methodg = geometry(mesh, x)
len = arclength(g, hs)
len = arclength(g, es)
Compute total length of (half-)edges. Use to compute arc-length of edge curves.
This method does not check whether input is a 1-manifold connected curve (see ishalfedgecurve
).
See also edgelength
, ishalfedgecurve
, arclengthknots
PMesh.arclengthknots
— Methodg = geometry(mesh, x)
ts = arclengthknots(g, halfedges)
Get knot values for each i
-th destination
vertex of halfedges
such that ts[i] == arclength(g, halfedges[1:i])
is a cumulative sum of arclength
s.
ts
does not contain the knot value 0
associated with the source
vertex of first(halfedges)
!
The knots may be used as input to positionhalfedgecurve!
.
This method does not check whether input is a 1-manifold connected curve (see ishalfedgecurve
).
See also arclength
, uniformknots
, positionhalfedgecurve!
PMesh.ishalfedgecurve
— Methodishalfedgecurve(mesh, hs[; closed=false])
Return true
if hs
is a sequence of connected half-edges in mesh
. If closed==true
, the first and the last half-edge must be connected.
See also arclength
PMesh.positionhalfedgecurve!
— Functiong = geometry(mesh, x)
positionhalfedgecurve!(f, g, halfedges
[; knots=arclengthknots(g, halfedges)])
positionhalfedgecurve!(f, g, halfedges;
knots = uniformknots(halfedges))
Set position
s of vertices that span the curve specified by halfedges
as f(t)
where t = knots[i] / knots[end] ∈ [0,1]
.
!! note The method sets x
("output"). Note that knots
may have to be be determined from a different attribute ("input").
See example below with input gx
and output gu
.
gx = geometry(mesh, x)
gu = geometry(mesh, u)
positionhalfedgecurve!(f, gu, halfedges;
knots = arclengthknots(gx, halfedges))
# short for using knots spaced by arc-length w.r.t gx
positionhalfedgecurve!(f, gu, halfedges, gx)
knots
is defined as returned by, e.g., arclengthknots
, i.e., the first knot t = 0
is given implicitly and not contained in knots
.
Example:
circle(t) = (cos(2π*t), sin*(2π*t), zero(t))
g = geometry(mesh, Val(:x))
loop = first(findboundaryloops(mesh))
positionhalfedgecurve!(circle, g, loop)
This method does not check whether input is a 1-manifold connected curve (see ishalfedgecurve
).
See also arclengthknots
, uniformknots
PMesh.uniformknots
— Functionts = uniformknots(halfedges)
g = geometry(mesh, x)
ts = uniformknots(g, halfedges)
Get uniformly spaced knot sequence similar to arclengthknots
.
This method doesn't require the mesh
geometry
!
ts
does not contain the knot value 0
associated with the source
vertex of first(halfedges)
!
The knots may be used as input to positionhalfedgecurve!
.
See also arclengthknots
, positionhalfedgecurve!