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
Note

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.)

Tip

Any custom defined subtype of AbstractGeometry that supports this interface can be used with any method/algorithm that takes an geometry instance.

Reference

PMesh.istriangulationFunction
istriangulation(geometry(...))

Query type traits: Is triangulation or general polygonal mesh?

Note

This method relies only on static type information. It does not check face degrees!

See also geometry

source
PMesh.AbstractGeometryType

Bind 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

source
PMesh.geometryFunction
g = 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).

Tip

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

Note
  1. geometry always requires vertex positions. You can but don't have to provide face normals and/or vertex normals.
  2. Providing or not providingg an attribute effects the use of stored data versus computing data on demand.
Warning

The current implememtation requires the same element types for positions and normals, e.g., SVector{3, Float64}.

See also triangulated

source
PMesh.triangulatedFunction
g = 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.

precondition

Requires that mesh is and remains a triangle mesh. This is a contract with the user, which is not checked!

See also geometry

source

Access mesh and attributes

Base.positionFunction
g = 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.

See also geometry, positions

source
PMesh.fnormalMethod
g = geometry(mesh, ..., fnormal=...)
fnrm = fnormal(g)

Get face normals attribute.

precondition

Requires has_fnormal(g). Otherwise method missing.

See also geometry, fnormals

source
PMesh.fnormalsMethod
g = geometry(mesh, ..., fnormal=...)

fn = fnormals(g)
nrm = fn(f)          # map FHnd to face normal

Create a functor that maps f::FHnd to stored fnormal.

precondition

Requires has_fnormal(g). Otherwise method missing.

See also geometry, maphandles, fnormal

source
PMesh.positionsFunction
fp = 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

source
PMesh.vnormalMethod
g = geometry(mesh, ..., vnormal=...)
vnrm = vnormal(g)

Get vertex normals attribute.

precondition

Requires has_fnormal(g). Otherwise method missing.

See also geometry, vnormals

source
PMesh.vnormalsMethod
g = geometry(mesh, ..., fnormal=...)

vn = vnormals(g)
nrm = vn(v)          # map VHnd to vertex normal

Create a functor that maps v::VHnd to stored vnormal.

precondition

Requires has_fnormal(g). Otherwise method missing.

See also geometry, maphandles, vnormal

source

Basic queries

  • Compute edgevector, their lengths, etc.
  • Compute (non-orthogonal) bases for face planes
PMesh.barycenterFunction
c = 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

source
PMesh.edgelengthFunction
len = 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

source
PMesh.edgelengthsFunction
lens = 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

source
PMesh.edgevectorFunction
d = 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

source
PMesh.meanedgelengthFunction
μ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

source
PMesh.rawfacebaseFunction
x0, 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.

precondition

f must be a planar face.

!! warning The caller must not make assumptions on the particular choice of vectors!

Todo

The current implementation uses an arbitrary pair of edgevectors and disregards conditioning/angles.

See also geometry, rawtrianglebase, facenormal, trianglenormal, tangentspace](@ref)

source
PMesh.rawtrianglebaseFunction
x0, 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.

precondition

Face f refers to a triangle.

See also geometry, rawfacebase, rawtrianglenormal, tangentspace

source

Areas

Todo

There is currently no method for computing the area of non-triangular faces.

PMesh.scalesurfacearea!Method
g = 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

source
PMesh.triangleareaFunction
ar = 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.

precondition

Face f refers to a triangle.

See also geometry, rawtrianglebase, trianglenormal

source

Angles

Base.angleMethod
α = angle(v1, v2)

Compute angle (rad) enclosed by vectors v1 and v2.

precondition

Both vectors v1, v2 must have positive length.

source
PMesh.anglenextMethod
g = 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.

See also angle, next, angleprev

source
PMesh.angleprevMethod
g = 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.

See also prev, anglenext, angle,

source
PMesh.dihedralangleFunction
g = 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).

precondition

h and e must denote an inner (half-)edge, i.e., both faces are well-defined.

Note

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

source

Statistics

  • Compute centroid of a single face and area-weighted centroid of the surface mesh
  • Define Statistics.mean, Statistics.std, and Statistics.cov for AbstractGeometry objects to compute mean, standard deviation and covariance of vertex positions
  • Center surface mesh by shifting mean or centroid using centermean! or centercentroid!
  • pca computes the principal component analysis of vertex positions or face centroids
  • extrema computes an axis aligned bounding box spanned by the returned vectors
  • radius computes the radius w.r.t to a prescribed center (or the mean)
Note

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.extremaMethod
g = geometry(mesh, x)

xmin, xmax = extrema(g)

Compute axis aligned bounding box spanned by extrema xmin, xmax.

See also radius

source
PMesh.centercentroid!Method
g = 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!

source
PMesh.centermean!Method
g = 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!

source
PMesh.centroidFunction
g = 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.

source
PMesh.pcaMethod
g = geometry(mesh, x)
U, λ, center = pca(mesh[; centroids=false])

Compute Principle Component Analysis of either position(g) or centroids (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

source

Normals

PMesh.trianglenormalFunction
nrm = trianglenormal(mesh, x, f)

g = geometry(mesh, x)
nrm = trianglenormal(g, f)

Compute unit-length normal vector to triangle f.

precondition

Face f refers to a triangle.

Note

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

source
PMesh.updatenormals!Function
g = geometry(mesh, x, ...)
updatenormals!(g[, weights=AngleWeights()])

Short for

updatefacenormals!(g)
updatevertexnormals!(g[, weights=AngleWeights()])
source
PMesh.updatevertexnormals!Function
g = geometry(mesh, x, ...)
updatevertexnormals!(g[, weights=AngleWeights()])

Recompute and update vertex normals using weights (see vertexnormal.

Note

This method has no effect if g does not provide an attribute for vertex normals.

Note

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!

source
PMesh.vertexnormalFunction
nrm = vertexnormal(mesh, x, v[, weights=AngleWeights()])

g = geometry(mesh, x)
nrm = vertexnormal(g, v[, weights=AngleWeights()])

Compute vertex normal from weighted sum of facenormals.

Note

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!

source

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.baryintFunction
b = 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

source
PMesh.BarycentricCoordinates.mapbarycoordsFunction

b = 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 λ.

Note

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":

Note

Mapping x to barycentric coordinates λ never fails, i.e., yields always finite barycentric coordinates λ. Consequently, one may have f(λ) != x.

See MapBaryCoords, barycoords, baryint

source
PMesh.BarycentricCoordinates.MapBaryCoordsType

Affine 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":

Note

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

source

Transformations

  • Represent AffineMaps and LinearMap (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.mapnormalFunction
f = affinemap(A, b)
g = mapnormal(f)
g = mapnormal(linearmap(f))

Nonlinear map that maps normal vectors n with inv(A)' * n followed by normalization.

Note

Null vectors are mapped to null vectors. They are not "normalized"!

See also LinearMap, mapvector, mapunitvector

source
PMesh.AffineMapType

Affine 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

See also affinemap, LinearMap

source
PMesh.LinearMapType

Linear 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

See also AffineMap, linearmap

source

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.linearmap2dFunction
g = 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.

precondition

This method requires a triangle mesh.

See also localframe

source

Tangent spaces

The methods use normal, i.e., read from attributes if available.

PMesh.localframeFunction
g = 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.

Note

This method uses stored normal if attribute is available.

Note

This method is similar to tangentspace but returns also the normal vector.

See also geometry, rawfacebase, normal, tangentspace, tolocal

source
PMesh.projectintotangentspaceFunction
g = 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.

Note

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

source
PMesh.tangentspaceMethod
g = 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.

Note

This method uses stored normal if attribute is available.

Note

This method is similar to localframe but does not return normal vector.

See also geometry, rawfacebase, normal, localframe, projectintotangentspace, tolocal2d

source
PMesh.tolocalFunction
g = 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.

Note

The local coordinate system for a vertex v has x[v] as origin, i.e., x[v] is mapped to the null vector.

Tip

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

source
PMesh.tolocal2dFunction
g = 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.

Note

The local coordinate system for a vertex v has x[v] as origin, i.e., x[v] is mapped to the null vector.

Tip

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

source
PMesh.toworldMethod
g = 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.

Note

toworld maps points in local coordinates to points in world coordinates by an affine map.

Use linearmap or mapvector to obtain a linear map to map vectors in local coordinates to vectors in world coordinates.

precondition

Requires rigid transformation as obtained from tolocal, i.e., A'A == I.

See also localframe, tolocal, tolocal2d

source

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.

Warning

Handles returned by nearestvertex and nearestface become invalid if the mesh changes after calling pointlocation.

Note

Consider using compact, i.e., ensure iscontiguous before calling pointlocation!

Sequences of connected edges

Treat sequences of connected edges as discrete curves on meshes.

PMesh.arclengthknotsMethod
g = 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 arclengths.

Note

ts does not contain the knot value 0 associated with the source vertex of first(halfedges)!

Note

The knots may be used as input to positionhalfedgecurve!.

Note

This method does not check whether input is a 1-manifold connected curve (see ishalfedgecurve).

See also arclength, uniformknots, positionhalfedgecurve!

source
PMesh.ishalfedgecurveMethod
ishalfedgecurve(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

source
PMesh.positionhalfedgecurve!Function
g = geometry(mesh, x)
positionhalfedgecurve!(f, g, halfedges
                       [; knots=arclengthknots(g, halfedges)])

positionhalfedgecurve!(f, g, halfedges;
                       knots = uniformknots(halfedges))

Set positions 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)
precondition

If halfedges specifies a closed curve, f must be periodic such that f(0) == f(1).

Note

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)
Note

This method does not check whether input is a 1-manifold connected curve (see ishalfedgecurve).

See also arclengthknots, uniformknots

source