Lagrange
Loading...
Searching...
No Matches
Mesh utilities

Various mesh processing utilities. More...

Classes

struct  FacetAreaOptions
 Option struct for computing per-facet area. More...
 
struct  FacetVectorAreaOptions
 
struct  MeshAreaOptions
 Option struct for computing mesh area. More...
 
struct  FacetCentroidOptions
 Option struct for computing per-facet centroid. More...
 
struct  MeshCentroidOptions
 Option struct for computing mesh centroid. More...
 
struct  ComponentOptions
 Options to control connected components computation. More...
 
struct  DihedralAngleOptions
 Option struct for computing dihedral angles. More...
 
struct  DijkstraDistanceOptions< Scalar, Index >
 Option struct for compute_dijkstra_distance. More...
 
struct  EdgeLengthOptions
 
struct  FacetCircumcenterOptions
 
struct  FacetNormalOptions
 Option struct for computing per-facet mesh normals. More...
 
struct  GreedyColoringOptions
 Option struct for computing dihedral angles. More...
 
struct  MeshCovarianceOptions
 Options struct for computing mesh covariance. More...
 
struct  NormalOptions
 Option struct for computing indexed mesh normals. More...
 
struct  ComputePointcloudPCAOptions
 
struct  PointcloudPCAOutput< Scalar >
 
struct  SeamEdgesOptions
 Options for computing seam edges. More...
 
struct  TangentBitangentOptions
 Option struct for computing tangent and bitangent vectors. More...
 
struct  TangentBitangentResult
 Result type of the compute_tangent_bitangent function. More...
 
struct  UVChartOptions
 
struct  UVDistortionOptions
 Option struct for compute uv distortion. More...
 
struct  VertexNormalOptions
 Option struct for computing per-vertex mesh normals. More...
 
struct  VertexValenceOptions
 Option struct for computing vertex valence. More...
 
struct  CornerNormalOptions
 Option struct for computing per-corner mesh normals. More...
 
struct  SubmeshOptions
 Options for extract submesh. More...
 
struct  AttributeFilter
 Helper object to filter attributes based on name, id, usage or element type. More...
 
struct  AttributeMatcher
 Helper object to match attributes based on usage, element type, and number of channels. More...
 
struct  IsolineOptions
 Options for isoline extraction/trimming. More...
 
struct  OrientOptions
 Options for orienting the facets of a mesh. More...
 
struct  OrientationOptions
 Option struct for computing if edges are oriented. More...
 
struct  RemapVerticesOptions
 Remap vertices options. More...
 
struct  SelectFacetsByNormalSimilarityOptions
 Option struct for selecting facets based on normal similarity. More...
 
struct  Frustum< Scalar >
 An array of four planes that define a frustum. More...
 
struct  FrustumSelectionOptions
 Option struct for selecting facets. More...
 
struct  SeparateByComponentsOptions
 Option settings for separate_by_components. More...
 
struct  SeparateByFacetGroupsOptions
 Option settings for separate_by_facet_groups. More...
 
struct  ThickenAndCloseOptions
 Options for thicken_and_close_mesh. More...
 
struct  VertexManifoldOptions
 Option struct for computing manifold flags. More...
 
struct  EdgeManifoldOptions
 Option struct for computing edge manifold flags. More...
 
struct  TriangulationOptions
 
struct  TransformOptions
 Options available when applying affine transforms to a mesh. More...
 
struct  UVMeshOptions
 

Enumerations

enum class  DistortionMetric {
  Dirichlet , InverseDirichlet , SymmetricDirichlet , AreaRatio ,
  MIPS
}
 UV distortion metric type. More...
 
enum class  NormalWeightingType : char { Uniform = 0 , CornerTriangleArea = 1 , Angle = 2 }
 Weighting types for averaging corner normals around a vertex. More...
 
enum class  ReorderingMethod { Lexicographic , Morton , Hilbert , None }
 Mesh reordering method to apply before decimation. More...
 

Functions

template<typename ToScalar, typename ToIndex, typename FromScalar, typename FromIndex>
SurfaceMesh< ToScalar, ToIndex > cast (const SurfaceMesh< FromScalar, FromIndex > &source_mesh, const AttributeFilter &convertible_attributes={}, std::vector< std::string > *converted_attributes_names=nullptr)
 Cast a mesh to a mesh of different scalar and/or index type.
 
template<typename ToValueType, typename Scalar, typename Index>
AttributeId cast_attribute (SurfaceMesh< Scalar, Index > &mesh, AttributeId source_id, std::string_view target_name)
 Cast an attribute in place to a different value type.
 
template<typename ToValueType, typename Scalar, typename Index>
AttributeId cast_attribute (SurfaceMesh< Scalar, Index > &mesh, std::string_view source_name, std::string_view target_name)
 Cast an attribute in place to a different value type.
 
template<typename ToValueType, typename Scalar, typename Index>
AttributeId cast_attribute_in_place (SurfaceMesh< Scalar, Index > &mesh, AttributeId attribute_id)
 Cast an attribute in place to a different value type.
 
template<typename ToValueType, typename Scalar, typename Index>
AttributeId cast_attribute_in_place (SurfaceMesh< Scalar, Index > &mesh, std::string_view name)
 Cast an attribute in place to a different value type.
 
template<typename Scalar, typename Index>
SurfaceMesh< Scalar, Index > combine_meshes (std::initializer_list< const SurfaceMesh< Scalar, Index > * > meshes, bool preserve_attributes=true)
 Combine multiple meshes into a single mesh.
 
template<typename Scalar, typename Index>
SurfaceMesh< Scalar, Index > combine_meshes (span< const SurfaceMesh< Scalar, Index > > meshes, bool preserve_attributes=true)
 Combine multiple meshes into a single mesh.
 
template<typename Scalar, typename Index>
SurfaceMesh< Scalar, Index > combine_meshes (size_t num_meshes, function_ref< const SurfaceMesh< Scalar, Index > &(size_t)> get_mesh, bool preserve_attributes=true)
 Combine multiple meshes into a single mesh.
 
template<typename Scalar, typename Index>
AttributeId compute_facet_area (SurfaceMesh< Scalar, Index > &mesh, FacetAreaOptions options={})
 Compute per-facet area.
 
template<typename Scalar, typename Index, int Dimension>
AttributeId compute_facet_area (SurfaceMesh< Scalar, Index > &mesh, const Eigen::Transform< Scalar, Dimension, Eigen::Affine > &transformation, FacetAreaOptions options={})
 Compute per-facet area.
 
template<typename Scalar, typename Index>
AttributeId compute_facet_vector_area (SurfaceMesh< Scalar, Index > &mesh, FacetVectorAreaOptions options={})
 Compute per-facet vector area.
 
template<typename Scalar, typename Index>
Scalar compute_mesh_area (const SurfaceMesh< Scalar, Index > &mesh, MeshAreaOptions options={})
 Compute mesh area.
 
template<typename Scalar, typename Index, int Dimension>
Scalar compute_mesh_area (const SurfaceMesh< Scalar, Index > &mesh, const Eigen::Transform< Scalar, Dimension, Eigen::Affine > &transformation, MeshAreaOptions options={})
 Compute mesh area.
 
template<typename Scalar, typename Index>
Scalar compute_uv_area (const SurfaceMesh< Scalar, Index > &mesh, MeshAreaOptions options={})
 Compute UV mesh area.
 
template<typename Scalar, typename Index>
AttributeId compute_facet_centroid (SurfaceMesh< Scalar, Index > &mesh, FacetCentroidOptions options={})
 Compute per-facet centroid.
 
template<typename Scalar, typename Index>
void compute_mesh_centroid (const SurfaceMesh< Scalar, Index > &mesh, span< Scalar > centroid, MeshCentroidOptions options={})
 Compute mesh centroid, where mesh centroid is defined as the weighted sum of facet centroids.
 
template<typename Scalar, typename Index>
size_t compute_components (SurfaceMesh< Scalar, Index > &mesh, ComponentOptions options={})
 Compute connected components of an input mesh.
 
template<typename Scalar, typename Index>
size_t compute_components (SurfaceMesh< Scalar, Index > &mesh, span< const Index > blocker_elements, ComponentOptions options={})
 Compute connected components of an input mesh.
 
template<typename Scalar, typename Index>
AttributeId compute_dihedral_angles (SurfaceMesh< Scalar, Index > &mesh, const DihedralAngleOptions &options={})
 Computes dihedral angles for each edge in the mesh.
 
template<typename Scalar, typename Index>
std::optional< std::vector< Index > > compute_dijkstra_distance (SurfaceMesh< Scalar, Index > &mesh, const DijkstraDistanceOptions< Scalar, Index > &options={})
 Computes dijkstra distance from a seed facet.
 
template<typename Scalar, typename Index>
AttributeId compute_edge_lengths (SurfaceMesh< Scalar, Index > &mesh, const EdgeLengthOptions &options={})
 Computes edge lengths attribute.
 
template<typename Scalar, typename Index>
AttributeId compute_facet_circumcenter (SurfaceMesh< Scalar, Index > &mesh, FacetCircumcenterOptions options={})
 Compute per-facet circumcenter.
 
template<typename Scalar, typename Index>
AdjacencyList< Index > compute_facet_facet_adjacency (SurfaceMesh< Scalar, Index > &mesh)
 Compute facet-facet adjacency information based on shared edges.
 
template<typename Scalar, typename Index>
AttributeId compute_facet_normal (SurfaceMesh< Scalar, Index > &mesh, FacetNormalOptions options={})
 Compute facet normals.
 
template<typename Scalar, typename Index>
AttributeId compute_greedy_coloring (SurfaceMesh< Scalar, Index > &mesh, const GreedyColoringOptions &options={})
 Compute a greedy graph coloring of the mesh.
 
template<typename Scalar, typename Index>
std::array< std::array< Scalar, 3 >, 3 > compute_mesh_covariance (const SurfaceMesh< Scalar, Index > &mesh, const MeshCovarianceOptions &options={})
 Compute the covariance matrix w.r.t.
 
template<typename Scalar, typename Index>
AttributeId compute_normal (SurfaceMesh< Scalar, Index > &mesh, function_ref< bool(Index)> is_edge_smooth, span< const Index > cone_vertices={}, NormalOptions options={})
 Compute smooth normals based on specified sharp edges and cone vertices.
 
template<typename Scalar, typename Index>
AttributeId compute_normal (SurfaceMesh< Scalar, Index > &mesh, function_ref< bool(Index, Index)> is_edge_smooth, span< const Index > cone_vertices={}, NormalOptions options={})
 Compute smooth normals based on specified sharp edges and cone vertices.
 
template<typename Scalar, typename Index>
AttributeId compute_normal (SurfaceMesh< Scalar, Index > &mesh, Scalar feature_angle_threshold, span< const Index > cone_vertices={}, NormalOptions options={})
 Compute smooth normal based on specified dihedral angle threshold and cone vertices.
 
template<typename Scalar>
PointcloudPCAOutput< Scalarcompute_pointcloud_pca (span< const Scalar > points, ComputePointcloudPCAOptions options={})
 Finds the principal components for a pointcloud.
 
template<typename Scalar, typename Index>
AttributeId compute_seam_edges (SurfaceMesh< Scalar, Index > &mesh, AttributeId indexed_attribute_id, const SeamEdgesOptions &options={})
 Computes the seam edges for a given indexed attribute.
 
template<typename Scalar, typename Index>
TangentBitangentResult compute_tangent_bitangent (SurfaceMesh< Scalar, Index > &mesh, TangentBitangentOptions options={})
 Compute mesh tangent and bitangent vectors orthogonal to the input mesh normals.
 
template<typename Scalar, typename Index>
size_t compute_uv_charts (SurfaceMesh< Scalar, Index > &mesh, const UVChartOptions &options={})
 Compute UV charts of an input mesh.
 
template<typename Scalar, typename Index>
AttributeId compute_uv_distortion (SurfaceMesh< Scalar, Index > &mesh, const UVDistortionOptions &options={})
 Compute uv distortion using the selected distortion measure.
 
template<typename Scalar, typename Index>
std::vector< std::pair< int32_t, int32_t > > compute_uv_tile_list (const SurfaceMesh< Scalar, Index > &mesh)
 Extract the list of all UV tiles that a mesh's parametrization spans.
 
template<typename Scalar, typename Index>
AttributeId compute_vertex_normal (SurfaceMesh< Scalar, Index > &mesh, VertexNormalOptions options={})
 Compute per-vertex normals based on specified weighting type.
 
template<typename Scalar, typename Index>
AttributeId compute_vertex_valence (SurfaceMesh< Scalar, Index > &mesh, VertexValenceOptions options={})
 Compute vertex valence.
 
template<typename Scalar, typename Index>
AdjacencyList< Index > compute_vertex_vertex_adjacency (SurfaceMesh< Scalar, Index > &mesh)
 Compute vertex-vertex adjacency information.
 
template<typename Scalar, typename Index>
AttributeId compute_weighted_corner_normal (SurfaceMesh< Scalar, Index > &mesh, CornerNormalOptions option={})
 Compute corner normals.
 
template<typename Scalar, typename Index, typename DerivedV, typename DerivedF>
SurfaceMesh< Scalar, Index > eigen_to_surface_mesh (const Eigen::MatrixBase< DerivedV > &V, const Eigen::MatrixBase< DerivedF > &F)
 Create a SurfaceMesh from a igl-style pair of matrices (V, F).
 
template<typename Scalar, typename Index>
std::vector< std::vector< Index > > extract_boundary_loops (const SurfaceMesh< Scalar, Index > &mesh)
 Extract boundary loops from a surface mesh.
 
template<typename Scalar, typename Index>
SurfaceMesh< Scalar, Index > extract_submesh (const SurfaceMesh< Scalar, Index > &mesh, span< const Index > selected_facets, const SubmeshOptions &options={})
 Extract a submesh that consists of a subset of the facets of the source mesh.
 
template<typename Scalar, typename Index>
std::vector< AttributeIdfiltered_attribute_ids (const SurfaceMesh< Scalar, Index > &mesh, const AttributeFilter &options)
 Create a list of attribute ids corresponding to the given filter.
 
template<typename Scalar, typename Index>
SurfaceMesh< Scalar, Index > filter_attributes (SurfaceMesh< Scalar, Index > source_mesh, const AttributeFilter &options={})
 Filters the attributes of mesh according to user specifications.
 
template<typename Scalar, typename Index>
std::optional< AttributeIdfind_matching_attribute (const SurfaceMesh< Scalar, Index > &mesh, const AttributeMatcher &options)
 Finds the first attribute with the specified usage/element type/number of channels.
 
template<typename Scalar, typename Index>
std::optional< AttributeIdfind_matching_attribute (const SurfaceMesh< Scalar, Index > &mesh, AttributeUsage usage)
 
template<typename Scalar, typename Index>
std::optional< AttributeIdfind_matching_attribute (const SurfaceMesh< Scalar, Index > &mesh, BitField< AttributeElement > element_types)
 
template<typename Scalar, typename Index>
std::vector< AttributeIdfind_matching_attributes (const SurfaceMesh< Scalar, Index > &mesh, const AttributeMatcher &options)
 Finds all attributes with the specified usage/element type/number of channels.
 
template<typename Scalar, typename Index>
std::vector< AttributeIdfind_matching_attributes (const SurfaceMesh< Scalar, Index > &mesh, AttributeUsage usage)
 Finds all attributes with the specified usage.
 
template<typename Scalar, typename Index>
std::vector< AttributeIdfind_matching_attributes (const SurfaceMesh< Scalar, Index > &mesh, BitField< AttributeElement > element_types)
 Finds all attributes with the specified element types.
 
template<typename Scalar, typename Index>
SurfaceMesh< Scalar, Index > trim_by_isoline (const SurfaceMesh< Scalar, Index > &mesh, const IsolineOptions &options={})
 Trim a mesh by the isoline of an implicit function defined on the mesh vertices/corners.
 
template<typename Scalar, typename Index>
SurfaceMesh< Scalar, Index > extract_isoline (const SurfaceMesh< Scalar, Index > &mesh, const IsolineOptions &options={})
 Extract the isoline of an implicit function defined on the mesh vertices/corners.
 
template<size_t Dimension, typename Scalar, typename Index>
Eigen::AlignedBox< Scalar, static_cast< int >(Dimension)> mesh_bbox (const SurfaceMesh< Scalar, Index > &mesh)
 Compute the axis-aligned bounding box of a mesh.
 
template<typename Scalar, typename Index, typename MeshType>
SurfaceMesh< Scalar, Index > to_surface_mesh_copy (const MeshType &mesh)
 Convert a legacy mesh object to a surface mesh object.
 
template<typename Scalar, typename Index, typename MeshType>
SurfaceMesh< Scalar, Index > to_surface_mesh_wrap (MeshType &&mesh)
 Wrap a legacy mesh object as a surface mesh object.
 
template<typename MeshType, typename Scalar, typename Index>
std::unique_ptr< MeshTypeto_legacy_mesh (const SurfaceMesh< Scalar, Index > &mesh)
 Convert a surface mesh object to a legacy mesh object.
 
template<size_t Dimension = 3, typename Scalar, typename Index>
auto normalize_mesh_with_transform (SurfaceMesh< Scalar, Index > &mesh, const TransformOptions &options={}) -> Eigen::Transform< Scalar, Dimension, Eigen::Affine >
 Normalize a mesh to fit in a unit box centered at the origin.
 
template<typename Scalar, typename Index>
void normalize_mesh (SurfaceMesh< Scalar, Index > &mesh, const TransformOptions &options={})
 Normalize a mesh to fit in a unit box centered at the origin.
 
template<size_t Dimension = 3, typename Scalar, typename Index>
auto normalize_meshes_with_transform (span< SurfaceMesh< Scalar, Index > * > meshes, const TransformOptions &options={}) -> Eigen::Transform< Scalar, Dimension, Eigen::Affine >
 Normalize a list of meshes to fit in a unit box centered at the origin.
 
template<typename Scalar, typename Index>
void normalize_meshes (span< SurfaceMesh< Scalar, Index > * > meshes, const TransformOptions &options={})
 Normalize a list of meshes to fit in a unit box centered at the origin.
 
template<typename Scalar, typename Index>
void orient_outward (lagrange::SurfaceMesh< Scalar, Index > &mesh, const OrientOptions &options={})
 Orient the facets of a mesh so that the signed volume of each connected component is positive or negative.
 
template<typename Scalar, typename Index>
bool is_oriented (const SurfaceMesh< Scalar, Index > &mesh)
 Check if a mesh is oriented.
 
template<typename Scalar, typename Index>
AttributeId compute_edge_is_oriented (SurfaceMesh< Scalar, Index > &mesh, const OrientationOptions &options={})
 Compute a mesh attribute indicating whether an edge is oriented.
 
template<typename Scalar, typename Index>
void permute_facets (SurfaceMesh< Scalar, Index > &mesh, span< const Index > new_to_old)
 Reorder facets of a mesh based on a given permutation.
 
template<typename Scalar, typename Index>
void permute_vertices (SurfaceMesh< Scalar, Index > &mesh, span< const Index > new_to_old)
 Reorder vertices of a mesh based on a given permutation.
 
template<typename Scalar, typename Index>
void remap_vertices (SurfaceMesh< Scalar, Index > &mesh, span< const Index > forward_mapping, RemapVerticesOptions options={})
 Remap vertices of a mesh based on provided forward mapping.
 
template<typename Scalar, typename Index>
void reorder_mesh (SurfaceMesh< Scalar, Index > &mesh, ReorderingMethod method)
 Mesh reordering to improve cache locality.
 
template<typename PointType>
auto segment_segment_squared_distance (const Eigen::MatrixBase< PointType > &U0, const Eigen::MatrixBase< PointType > &U1, const Eigen::MatrixBase< PointType > &V0, const Eigen::MatrixBase< PointType > &V1, Eigen::PlainObjectBase< PointType > &closest_pointU, Eigen::PlainObjectBase< PointType > &closest_pointV, ScalarOf< PointType > &lambdaU, ScalarOf< PointType > &lambdaV) -> ScalarOf< PointType >
 Computes the squared distance between two N-d line segments, and the closest pair of points whose separation is this distance.
 
template<typename Scalar, typename Index>
AttributeId select_facets_by_normal_similarity (SurfaceMesh< Scalar, Index > &mesh, const Index seed_facet_id, const SelectFacetsByNormalSimilarityOptions &options={})
 Given a seed facet, selects facets around it based on the change in triangle normals.
 
template<typename Scalar, typename Index>
bool select_facets_in_frustum (SurfaceMesh< Scalar, Index > &mesh, const Frustum< Scalar > &frustum, const FrustumSelectionOptions &options={})
 Select all facets that intersect the cone/frustrum bounded by 4 planes defined by (n_i, p_i), where n_i is the plane normal and p_i is a point on the plane.
 
template<typename Scalar, typename Index>
std::vector< SurfaceMesh< Scalar, Index > > separate_by_components (const SurfaceMesh< Scalar, Index > &mesh, const SeparateByComponentsOptions &options={})
 Separate a mesh by connected components.
 
template<typename Scalar, typename Index>
std::vector< SurfaceMesh< Scalar, Index > > separate_by_facet_groups (const SurfaceMesh< Scalar, Index > &mesh, size_t num_groups, span< const Index > facet_group_indices, const SeparateByFacetGroupsOptions &options={})
 Extract a set of submeshes based on facet groups.
 
template<typename Scalar, typename Index>
std::vector< SurfaceMesh< Scalar, Index > > separate_by_facet_groups (const SurfaceMesh< Scalar, Index > &mesh, span< const Index > facet_group_indices, const SeparateByFacetGroupsOptions &options={})
 Extract a set of submeshes based on facet groups.
 
template<typename Scalar, typename Index>
std::vector< SurfaceMesh< Scalar, Index > > separate_by_facet_groups (const SurfaceMesh< Scalar, Index > &mesh, size_t num_groups, function_ref< Index(Index)> get_facet_group, const SeparateByFacetGroupsOptions &options={})
 Extract a set of submeshes based on facet groups.
 
template<typename Scalar, typename Index>
void split_facets_by_material (SurfaceMesh< Scalar, Index > &mesh, std::string_view material_attribute_name)
 Split mesh facets based on material labels.
 
template<typename Scalar, typename Index>
SurfaceMesh< Scalar, Index > thicken_and_close_mesh (SurfaceMesh< Scalar, Index > input_mesh, const ThickenAndCloseOptions &options={})
 Thicken a mesh by offsetting it, and close the shape into a thick 3D solid.
 
template<typename Scalar, typename Index>
int compute_euler (const SurfaceMesh< Scalar, Index > &mesh)
 Compute Euler characteristic of a mesh.
 
template<typename Scalar, typename Index>
bool is_closed (const SurfaceMesh< Scalar, Index > &mesh)
 Check if a mesh is closed.
 
template<typename Scalar, typename Index>
bool is_vertex_manifold (const SurfaceMesh< Scalar, Index > &mesh)
 Check if a mesh is vertex-manifold.
 
template<typename Scalar, typename Index>
bool is_edge_manifold (const SurfaceMesh< Scalar, Index > &mesh)
 Check if a mesh is edge-manifold.
 
template<typename Scalar, typename Index>
bool is_manifold (const SurfaceMesh< Scalar, Index > &mesh)
 Check if a mesh is both vertex-manifold and edge-manifold.
 
template<typename Scalar, typename Index>
AttributeId compute_vertex_is_manifold (SurfaceMesh< Scalar, Index > &mesh, const VertexManifoldOptions &options={})
 Compute a mesh attribute of value type uint8_t indicating vertex manifoldness.
 
template<typename Scalar, typename Index>
AttributeId compute_edge_is_manifold (SurfaceMesh< Scalar, Index > &mesh, const EdgeManifoldOptions &options={})
 Compute a mesh attribute of value type uint8_t indicating edge manifoldness.
 
template<typename Scalar, typename Index, int Dimension>
void transform_mesh (SurfaceMesh< Scalar, Index > &mesh, const Eigen::Transform< Scalar, Dimension, Eigen::Affine > &transform, const TransformOptions &options={})
 Apply an affine transform \( M \) to a mesh in-place.
 
template<typename Scalar, typename Index, int Dimension>
SurfaceMesh< Scalar, Index > transformed_mesh (SurfaceMesh< Scalar, Index > mesh, const Eigen::Transform< Scalar, Dimension, Eigen::Affine > &transform, const TransformOptions &options={})
 Apply an affine transform to a mesh and return the transformed mesh.
 
template<typename Scalar, typename Index>
void triangulate_polygonal_facets (SurfaceMesh< Scalar, Index > &mesh, const TriangulationOptions &options={})
 Triangulate polygonal facets of a mesh using a prescribed set of rules.
 
template<typename Scalar, typename Index, typename UVScalar = Scalar>
SurfaceMesh< UVScalar, Index > uv_mesh_ref (SurfaceMesh< Scalar, Index > &mesh, const UVMeshOptions &options={})
 Extract a UV mesh reference from an input mesh.
 
template<typename Scalar, typename Index, typename UVScalar = Scalar>
SurfaceMesh< UVScalar, Index > uv_mesh_view (const SurfaceMesh< Scalar, Index > &mesh, const UVMeshOptions &options={})
 Extract a UV mesh view from an input mesh.
 

Detailed Description

Various mesh processing utilities.

Enumeration Type Documentation

◆ DistortionMetric

enum class DistortionMetric
strong

#include <lagrange/DistortionMetric.h>

UV distortion metric type.

Enumerator
Dirichlet 

Dirichlet energy.

InverseDirichlet 

Inverse Dirichlet energy.

SymmetricDirichlet 

Symmetric Dirichlet energy.

MIPS 

UV triangle area / 3D triangle area.

◆ NormalWeightingType

enum class NormalWeightingType : char
strong

#include <lagrange/NormalWeightingType.h>

Weighting types for averaging corner normals around a vertex.

Enumerator
Uniform 

Incident face normals have uniform influence on vertex normal.

CornerTriangleArea 

Incident face normals are averaged weighted by area of the corner triangle.

Angle 

Incident face normals are averaged weighted by incident angle of vertex.

◆ ReorderingMethod

enum class ReorderingMethod
strong

#include <lagrange/reorder_mesh.h>

Mesh reordering method to apply before decimation.

Enumerator
Lexicographic 

Sort vertices/facets lexicographically.

Morton 

Spatial sort vertices/facets using Morton encoding.

Hilbert 

Spatial sort vertices/facets using Hilbert curve.

None 

Do not reorder mesh vertices/facets.

Function Documentation

◆ cast()

template<typename ToScalar, typename ToIndex, typename FromScalar, typename FromIndex>
SurfaceMesh< ToScalar, ToIndex > cast ( const SurfaceMesh< FromScalar, FromIndex > & source_mesh,
const AttributeFilter & convertible_attributes = {},
std::vector< std::string > * converted_attributes_names = nullptr )

#include <lagrange/cast.h>

Cast a mesh to a mesh of different scalar and/or index type.

Note
To filter only certain attributes prior to casting a mesh, use the filter_attributes function.
Parameters
[in]source_meshInput mesh.
[in]convertible_attributesFilter to determine which attribute are convertible.
[out]converted_attributes_namesOptional output arg storing the list of non-reserved attribute names that were actually converted to a different type.
Template Parameters
ToScalarScalar type of the output mesh.
ToIndexIndex type of the output mesh.
FromScalarScalar type of the input mesh.
FromIndexIndex type of the input mesh.
Returns
Output mesh.
See also
filter_attributes

◆ cast_attribute() [1/2]

template<typename ToValueType, typename Scalar, typename Index>
AttributeId cast_attribute ( SurfaceMesh< Scalar, Index > & mesh,
AttributeId source_id,
std::string_view target_name )

#include <lagrange/cast_attribute.h>

Cast an attribute in place to a different value type.

This effectively replaces the existing attribute with a new one.

Parameters
[in,out]meshInput mesh. Modified to replace the attribute being cast.
[in]source_idId of the source attribute to cast.
[in]target_nameName of the target attribute to be created.
Template Parameters
ToValueTypeTarget value type for the attribute.
ScalarMesh scalar type.
IndexMesh value type.

◆ cast_attribute() [2/2]

template<typename ToValueType, typename Scalar, typename Index>
AttributeId cast_attribute ( SurfaceMesh< Scalar, Index > & mesh,
std::string_view source_name,
std::string_view target_name )

#include <lagrange/cast_attribute.h>

Cast an attribute in place to a different value type.

This effectively replaces the existing attribute with a new one.

Parameters
[in,out]meshInput mesh. Modified to replace the attribute being cast.
[in]source_nameName of the source attribute being cast.
[in]target_nameName of the target attribute to be created.
Template Parameters
ToValueTypeTarget value type for the attribute.
ScalarMesh scalar type.
IndexMesh value type.

◆ cast_attribute_in_place() [1/2]

template<typename ToValueType, typename Scalar, typename Index>
AttributeId cast_attribute_in_place ( SurfaceMesh< Scalar, Index > & mesh,
AttributeId attribute_id )

#include <lagrange/cast_attribute.h>

Cast an attribute in place to a different value type.

This effectively replaces the existing attribute with a new one.

Parameters
[in,out]meshInput mesh. Modified to replace the attribute being cast.
[in]attribute_idId of the attribute to cast.
Template Parameters
ToValueTypeTarget value type for the attribute.
ScalarMesh scalar type.
IndexMesh value type.

◆ cast_attribute_in_place() [2/2]

template<typename ToValueType, typename Scalar, typename Index>
AttributeId cast_attribute_in_place ( SurfaceMesh< Scalar, Index > & mesh,
std::string_view name )

#include <lagrange/cast_attribute.h>

Cast an attribute in place to a different value type.

This effectively replaces the existing attribute with a new one.

Parameters
[in,out]meshInput mesh. Modified to replace the attribute being cast.
[in]nameName of the attribute to cast.
Template Parameters
ToValueTypeTarget value type for the attribute.
ScalarMesh scalar type.
IndexMesh value type.

◆ combine_meshes() [1/3]

template<typename Scalar, typename Index>
SurfaceMesh< Scalar, Index > combine_meshes ( std::initializer_list< const SurfaceMesh< Scalar, Index > * > meshes,
bool preserve_attributes = true )

#include <lagrange/combine_meshes.h>

Combine multiple meshes into a single mesh.

Parameters
[in]meshesThe set of input mesh pointers.
[in]preserve_attributesPreserve shared attributes and map them to the output mesh.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
The combined mesh.

◆ combine_meshes() [2/3]

template<typename Scalar, typename Index>
SurfaceMesh< Scalar, Index > combine_meshes ( span< const SurfaceMesh< Scalar, Index > > meshes,
bool preserve_attributes = true )

#include <lagrange/combine_meshes.h>

Combine multiple meshes into a single mesh.

Parameters
[in]meshesMeshes to combine.
[in]preserve_attributesPreserve shared attributes and map them to the output mesh.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
The combined mesh.

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

◆ combine_meshes() [3/3]

template<typename Scalar, typename Index>
SurfaceMesh< Scalar, Index > combine_meshes ( size_t num_meshes,
function_ref< const SurfaceMesh< Scalar, Index > &(size_t)> get_mesh,
bool preserve_attributes = true )

#include <lagrange/combine_meshes.h>

Combine multiple meshes into a single mesh.

This is the most generic version, where get_mesh(i) provides the ith mesh.

Parameters
[in]num_meshesNumber of meshes to combine.
[in]get_meshRetrieve the i-th mesh to combine.
[in]preserve_attributesPreserve shared attributes and map them to the output mesh.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
The combined mesh.

This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.

◆ compute_facet_area() [1/2]

template<typename Scalar, typename Index>
AttributeId compute_facet_area ( SurfaceMesh< Scalar, Index > & mesh,
FacetAreaOptions options = {} )

#include <lagrange/compute_area.h>

Compute per-facet area.

Parameters
[in,out]meshThe input mesh.
[in]optionsThe options controlling the computation.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
The attribute id of the facet area attribute.
See also
FacetAreaOptions

◆ compute_facet_area() [2/2]

template<typename Scalar, typename Index, int Dimension>
AttributeId compute_facet_area ( SurfaceMesh< Scalar, Index > & mesh,
const Eigen::Transform< Scalar, Dimension, Eigen::Affine > & transformation,
FacetAreaOptions options = {} )

#include <lagrange/compute_area.h>

Compute per-facet area.

Parameters
[in,out]meshThe input mesh.
[in]transformationAffine transformation to apply on mesh geometry.
[in]optionsThe options controlling the computation.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
DimensionMesh dimension.
Returns
The attribute id of the facet area attribute.
See also
FacetAreaOptions

◆ compute_facet_vector_area()

template<typename Scalar, typename Index>
AttributeId compute_facet_vector_area ( SurfaceMesh< Scalar, Index > & mesh,
FacetVectorAreaOptions options = {} )

#include <lagrange/compute_area.h>

Compute per-facet vector area.

The vector area of a facet is defined as the facet's area multiplied by its normal vector. For triangular facets, it is equivalent to the cross product of two edges divided by 2. For non-planar polygonal facets, the vector area offers a more robust measure of both the facet's orientation and area. The magnitude of the vector area corresponds to the largest area of any orthogonal planar projection of the facet, and its direction aligns with the normal of the projection plane [1, 2].

References:

[1] Sullivan, John M. "Curvatures of smooth and discrete surfaces." Discrete differential geometry. Basel: Birkhäuser Basel, 2008. 175-188.

[2] Alexa, Marc, and Max Wardetzky. "Discrete Laplacians on general polygonal meshes." ACM SIGGRAPH 2011 papers. 2011. 1-10.

Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Parameters
[in,out]meshThe input mesh.
[in]optionsThe options controlling the computation.
Returns
The attribute id of the facet vector area attribute.

◆ compute_mesh_area() [1/2]

template<typename Scalar, typename Index>
Scalar compute_mesh_area ( const SurfaceMesh< Scalar, Index > & mesh,
MeshAreaOptions options = {} )

#include <lagrange/compute_area.h>

Compute mesh area.

Parameters
[in]meshThe input mesh.
[in]optionsThe options controlling the computation.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
The computed mesh area.
See also
MeshAreaOptions

◆ compute_mesh_area() [2/2]

template<typename Scalar, typename Index, int Dimension>
Scalar compute_mesh_area ( const SurfaceMesh< Scalar, Index > & mesh,
const Eigen::Transform< Scalar, Dimension, Eigen::Affine > & transformation,
MeshAreaOptions options = {} )

#include <lagrange/compute_area.h>

Compute mesh area.

Parameters
[in]meshThe input mesh.
[in]transformationAffine transformation to apply on mesh geometry.
[in]optionsThe options controlling the computation.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
DimensionMesh dimension.
Returns
The computed mesh area.
See also
MeshAreaOptions

◆ compute_uv_area()

template<typename Scalar, typename Index>
Scalar compute_uv_area ( const SurfaceMesh< Scalar, Index > & mesh,
MeshAreaOptions options = {} )

#include <lagrange/compute_area.h>

Compute UV mesh area.

It is equivalent to calling compute_mesh_area on the UV mesh extracted with uv_mesh_view. However, this function will also work if UV attribute is using a different scalar type than the mesh vertex positions.

Parameters
[in]meshThe input mesh.
[in]optionsThe options controlling the computation.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
The computed mesh area.
See also
MeshAreaOptions, uv_mesh_view

◆ compute_facet_centroid()

template<typename Scalar, typename Index>
AttributeId compute_facet_centroid ( SurfaceMesh< Scalar, Index > & mesh,
FacetCentroidOptions options = {} )

#include <lagrange/compute_centroid.h>

Compute per-facet centroid.

Parameters
[in,out]meshThe input mesh.
[in]optionsOption settings to control the computation.
Template Parameters
ScalarMesh Scalar type.
IndexMesh Index type.
Returns
The id of the facet centroid attribute.
See also
FacetCentroidOptions

◆ compute_mesh_centroid()

template<typename Scalar, typename Index>
void compute_mesh_centroid ( const SurfaceMesh< Scalar, Index > & mesh,
span< Scalar > centroid,
MeshCentroidOptions options = {} )

#include <lagrange/compute_centroid.h>

Compute mesh centroid, where mesh centroid is defined as the weighted sum of facet centroids.

Parameters
[in]meshThe input mesh.
[out]centroidThe buffer to store centroid coordinates.
[in]optionsOption settings to control the computation.
Template Parameters
ScalarMesh Scalar type.
IndexMesh Index type.
See also
MeshCentroidOptions

◆ compute_components() [1/2]

template<typename Scalar, typename Index>
size_t compute_components ( SurfaceMesh< Scalar, Index > & mesh,
ComponentOptions options = {} )

#include <lagrange/compute_components.h>

Compute connected components of an input mesh.

This method will create a per-facet component id in an attribute named ComponentOptions::output_attribute_name. Each component id is in [0, num_components-1].

Parameters
meshInput mesh.
optionsOptions to control component computation.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
The total number of connected components.
See also
ComponentOptions

◆ compute_components() [2/2]

template<typename Scalar, typename Index>
size_t compute_components ( SurfaceMesh< Scalar, Index > & mesh,
span< const Index > blocker_elements,
ComponentOptions options = {} )

#include <lagrange/compute_components.h>

Compute connected components of an input mesh.

This method will create a per-facet component id in an attribute named ComponentOptions::output_attribute_name. Each component id is in [0, num_components-1].

Parameters
meshInput mesh.
blocker_elementsAn array of blocker element indices. The blocker element index is either a vertex index or an edge index depending on options.connectivity_type. If options.connectivity_type is ConnectivityType::Edge, facets adjacent to a blocker edge are not considered as connected through this edge. If options.connectivity_type is ConnectivityType::Vertex, facets sharing a blocker vertex are not considered as connected through this vertex. If empty, no blocker elements are used.
optionsOptions to control component computation.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
The total number of connected components.
See also
ComponentOptions

◆ compute_dihedral_angles()

template<typename Scalar, typename Index>
AttributeId compute_dihedral_angles ( SurfaceMesh< Scalar, Index > & mesh,
const DihedralAngleOptions & options = {} )

#include <lagrange/compute_dihedral_angles.h>

Computes dihedral angles for each edge in the mesh.

The dihedral angle of an edge is defined as the angle between the normals of two facets adjacent to the edge. The dihedral angle is always in the range \([0, \pi]\) for manifold edges. For boundary edges, the dihedral angle defaults to 0. For non-manifold edges, the dihedral angle is not well-defined and will be set to the special value \( 2\pi \).

Template Parameters
ScalarMesh scalar type
IndexMesh index type
Parameters
[in]meshThe input mesh.
[in]optionsOptions for computing dihedral angles.
Returns
The id of the dihedral angle attribute.
See also
DihedralAngleOptions

◆ compute_dijkstra_distance()

template<typename Scalar, typename Index>
std::optional< std::vector< Index > > compute_dijkstra_distance ( SurfaceMesh< Scalar, Index > & mesh,
const DijkstraDistanceOptions< Scalar, Index > & options = {} )

#include <lagrange/compute_dijkstra_distance.h>

Computes dijkstra distance from a seed facet.

Parameters
meshInput mesh.
optionsOptions for computing dijkstra distance.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
Optionally, a vector of indices of vertices involved

◆ compute_edge_lengths()

template<typename Scalar, typename Index>
AttributeId compute_edge_lengths ( SurfaceMesh< Scalar, Index > & mesh,
const EdgeLengthOptions & options = {} )

#include <lagrange/compute_edge_lengths.h>

Computes edge lengths attribute.

Template Parameters
ScalarMesh scalar type
IndexMesh index type
Parameters
meshThe input mesh
optionsOptions for computing edge lengths.
Returns
Attribute ID of the computed edge lengths attribute.

◆ compute_facet_circumcenter()

template<typename Scalar, typename Index>
AttributeId compute_facet_circumcenter ( SurfaceMesh< Scalar, Index > & mesh,
FacetCircumcenterOptions options = {} )

#include <lagrange/compute_facet_circumcenter.h>

Compute per-facet circumcenter.

Parameters
[in,out]meshInput mesh.
[in]optionsOption settings to control the computation.
Returns
Attribute ID of the facet circumcenters.

◆ compute_facet_facet_adjacency()

template<typename Scalar, typename Index>
AdjacencyList< Index > compute_facet_facet_adjacency ( SurfaceMesh< Scalar, Index > & mesh)

#include <lagrange/compute_facet_facet_adjacency.h>

Compute facet-facet adjacency information based on shared edges.

Two facets are considered adjacent if they share an edge. For non-manifold edges with 3 or more incident facets, a complete clique is formed (every pair of facets around the edge is adjacent).

Note
If two facets share multiple edges, the neighboring facet will appear in the adjacency list once for each time the shared edge is referenced by its incident facets (i.e., neighbors are not deduplicated).
This function calls initialize_edges() if edges have not been initialized yet, which mutates the mesh.
Parameters
meshThe input mesh (edges will be initialized if needed).
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
The facet-facet adjacency list.

◆ compute_facet_normal()

template<typename Scalar, typename Index>
AttributeId compute_facet_normal ( SurfaceMesh< Scalar, Index > & mesh,
FacetNormalOptions options = {} )

#include <lagrange/compute_facet_normal.h>

Compute facet normals.

Parameters
[in,out]meshThe input mesh.
[in]optionsOptional arguments to control normal generation.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
AttributeId The attribute id of the facet normal attribute.
Postcondition
The computed facet normals are stored in mesh as a facet attribute named options.output_attribute_name.
Note
Non-planar polygonal facet's normal is not well defined. This method can only compute an approximated normal using a triangle fan.
See also
FacetNormalOptions

◆ compute_greedy_coloring()

template<typename Scalar, typename Index>
AttributeId compute_greedy_coloring ( SurfaceMesh< Scalar, Index > & mesh,
const GreedyColoringOptions & options = {} )

#include <lagrange/compute_greedy_coloring.h>

Compute a greedy graph coloring of the mesh.

The selected mesh element type (either Vertex or Facet) will be colored in such a way that no two adjacent element share the same color.

Parameters
[in,out]meshInput mesh to be colored. Modified to compute edge information and the new color attribute.
[in]optionsColoring options.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
Id of the newly computed color attribute. The value type of the created attribute will be the same as the mesh index type.

◆ compute_mesh_covariance()

template<typename Scalar, typename Index>
std::array< std::array< Scalar, 3 >, 3 > compute_mesh_covariance ( const SurfaceMesh< Scalar, Index > & mesh,
const MeshCovarianceOptions & options = {} )

#include <lagrange/compute_mesh_covariance.h>

Compute the covariance matrix w.r.t.

a given center (default to zeros).

Parameters
[in]meshThe input mesh.
[in]optionsOptional arguments.
Template Parameters
ScalarThe scalar type of the mesh.
IndexThe index type of the mesh.
Returns
The covariance matrix in column-major order (but it should be symmetric).
See also
MeshCovarianceOptions.
Note
Adapted from https://github.com/mkazhdan/ShapeSPH/blob/master/Util/TriangleMesh.h#L101

◆ compute_normal() [1/3]

template<typename Scalar, typename Index>
AttributeId compute_normal ( SurfaceMesh< Scalar, Index > & mesh,
function_ref< bool(Index)> is_edge_smooth,
span< const Index > cone_vertices = {},
NormalOptions options = {} )

#include <lagrange/compute_normal.h>

Compute smooth normals based on specified sharp edges and cone vertices.

Parameters
[in]meshThe input mesh.
[in]is_edge_smoothReturns true on e if the edge is smooth.
[in]cone_verticesA list of cone vertices.
[in]optionsOptional arguments to control normal generation.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
The indexed attribute id of normal attribute.
See also
NormalOptions

◆ compute_normal() [2/3]

template<typename Scalar, typename Index>
AttributeId compute_normal ( SurfaceMesh< Scalar, Index > & mesh,
function_ref< bool(Index, Index)> is_edge_smooth,
span< const Index > cone_vertices = {},
NormalOptions options = {} )

#include <lagrange/compute_normal.h>

Compute smooth normals based on specified sharp edges and cone vertices.

Parameters
[in]meshThe input mesh.
[in]is_edge_smoothReturns true on (fi, fj) if the edge between fi and fj is smooth. Assumes fi and fi are adjacent.
[in]cone_verticesA list of cone vertices.
[in]optionsOptional arguments to control normal generation.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
The indexed attribute id of normal attribute.
See also
NormalOptions

◆ compute_normal() [3/3]

template<typename Scalar, typename Index>
AttributeId compute_normal ( SurfaceMesh< Scalar, Index > & mesh,
Scalar feature_angle_threshold,
span< const Index > cone_vertices = {},
NormalOptions options = {} )

#include <lagrange/compute_normal.h>

Compute smooth normal based on specified dihedral angle threshold and cone vertices.

Parameters
[in]meshThe input mesh.
[in]feature_angle_thresholdAn edge with dihedral angle larger than this threshold is considered as an feature edge. The angle is expressed in radian.
[in]cone_verticesA list of cone vertices.
[in]optionsOptional arguments to control normal generation.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
The indexed attribute id of normal attribute.
See also
NormalOptions

◆ compute_pointcloud_pca()

template<typename Scalar>
PointcloudPCAOutput< Scalar > compute_pointcloud_pca ( span< const Scalar > points,
ComputePointcloudPCAOptions options = {} )

#include <lagrange/compute_pointcloud_pca.h>

Finds the principal components for a pointcloud.

Assumes that the points are supplied in a matrix where each `‘row’' is a point.

This is closely related to the inertia tensor, principal directions and principal moments. But it is not exactly the same.

COVARIANCE_MATRIX (tensor) = (P-eC)^T (P-eC) Where C is the centroid, and e is column vector of ones. eigenvalues and eigenvectors of this matrix would be the principal weights and components.

MOMENT OF INERTIA (tensor) = trace(P^T P)I - P^T P eigenvalues and eigenvectors of this matrix would be the principal moments and directions.

Parameters
[in]pointsThe input points.
[in]optionsOptions struct, see above.
Returns
PointCloudPCAOutput Contains center, eigenvectors (components) and eigenvalues (weights). See above for more details.

◆ compute_seam_edges()

template<typename Scalar, typename Index>
AttributeId compute_seam_edges ( SurfaceMesh< Scalar, Index > & mesh,
AttributeId indexed_attribute_id,
const SeamEdgesOptions & options = {} )

#include <lagrange/compute_seam_edges.h>

Computes the seam edges for a given indexed attribute.

A seam edge is an edge where either of its endpoint vertex has corners whose attribute index differ across the edge.

Parameters
[in,out]meshInput mesh. Modified to add edge information and compute output seam edges.
[in]indexed_attribute_idId of the indexed attribute for which to compute seam edges.
[in]optionsOptional parameters.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
Id of the newly added per-edge attribute. The newly computed attribute will have value type uint8_t, and contain 0 for non-seam edges, and 1 for seam edges.

◆ compute_tangent_bitangent()

template<typename Scalar, typename Index>
TangentBitangentResult compute_tangent_bitangent ( SurfaceMesh< Scalar, Index > & mesh,
TangentBitangentOptions options = {} )

#include <lagrange/compute_tangent_bitangent.h>

Compute mesh tangent and bitangent vectors orthogonal to the input mesh normals.

Note
Unless options.keep_existing_tangent is true, the input mesh must have existing indexed normal and UV attributes. The input UV attribute is used to orient the resulting T/B vectors coherently wrt to the UV mapping.
Parameters
[in]meshThe input mesh.
[in]optionsOptional arguments to control tangent/bitangent generation.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
A struct containing the id of the generated tangent/bitangent attributes.
See also
TangentBitangentOptions

◆ compute_uv_charts()

template<typename Scalar, typename Index>
size_t compute_uv_charts ( SurfaceMesh< Scalar, Index > & mesh,
const UVChartOptions & options = {} )

#include <lagrange/compute_uv_charts.h>

Compute UV charts of an input mesh.

This method will create a per-vertex chart id in an attribute named UVChartOptions::output_attribute_name. Each chart id is in [0, num_charts-1].

Parameters
meshInput mesh.
optionsOptions to control chart computation.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
The total number of UV charts.
See also
UVChartOptions

◆ compute_uv_distortion()

template<typename Scalar, typename Index>
AttributeId compute_uv_distortion ( SurfaceMesh< Scalar, Index > & mesh,
const UVDistortionOptions & options = {} )

#include <lagrange/compute_uv_distortion.h>

Compute uv distortion using the selected distortion measure.

Parameters
[in]meshThe input mesh.
[in]optionsThe computation option settings.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
The attribute id of the distortion measure facet attribute.
See also
UVDistortionOptions

◆ compute_uv_tile_list()

template<typename Scalar, typename Index>
std::vector< std::pair< int32_t, int32_t > > compute_uv_tile_list ( const SurfaceMesh< Scalar, Index > & mesh)

#include <lagrange/compute_uv_tile_list.h>

Extract the list of all UV tiles that a mesh's parametrization spans.

UV tiles are usually understood to be a regular unit grid in UV space. This process thus reads UV for all vertices of the input mesh, and adds an entry to the output for each new integer pair that it finds.

Parameters
meshMesh to be analyzed
Returns
A list of integer coordinates with one entry for each UV tile of mesh.

◆ compute_vertex_normal()

template<typename Scalar, typename Index>
AttributeId compute_vertex_normal ( SurfaceMesh< Scalar, Index > & mesh,
VertexNormalOptions options = {} )

#include <lagrange/compute_vertex_normal.h>

Compute per-vertex normals based on specified weighting type.

Parameters
[in]meshThe input mesh.
[in]optionsOptional arguments to control normal generation.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
The attribute id of vertex normal attribute.
See also
VertexNormalOptions

◆ compute_vertex_valence()

template<typename Scalar, typename Index>
AttributeId compute_vertex_valence ( SurfaceMesh< Scalar, Index > & mesh,
VertexValenceOptions options = {} )

#include <lagrange/compute_vertex_valence.h>

Compute vertex valence.

Parameters
meshThe input mesh.
optionsOptional settings to control valence computation.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
The vertex attribute id containing valence information.
See also
VertexValenceOptions

◆ compute_vertex_vertex_adjacency()

template<typename Scalar, typename Index>
AdjacencyList< Index > compute_vertex_vertex_adjacency ( SurfaceMesh< Scalar, Index > & mesh)

#include <lagrange/compute_vertex_vertex_adjacency.h>

Compute vertex-vertex adjacency information.

Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Parameters
meshThe input mesh.
Returns
The vertex-vertex adjacency data and adjacency indices.

◆ compute_weighted_corner_normal()

template<typename Scalar, typename Index>
AttributeId compute_weighted_corner_normal ( SurfaceMesh< Scalar, Index > & mesh,
CornerNormalOptions option = {} )

#include <lagrange/compute_weighted_corner_normal.h>

Compute corner normals.

Parameters
[in]meshThe input mesh.
[in]optionOptional arguments to control normal generation.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
The corner attribute id of corner normal attribute.
Note
If options.weight_type is not Uniform, the resulting corner normal's lengths encodes the corresponding weight.
Corner normals around a given vertex could be different even when the vertex is at a smooth region. For computing smooth normal, use compute_normal instead.
See also
CornerNormalOptions, compute_normal

◆ eigen_to_surface_mesh()

template<typename Scalar, typename Index, typename DerivedV, typename DerivedF>
SurfaceMesh< Scalar, Index > eigen_to_surface_mesh ( const Eigen::MatrixBase< DerivedV > & V,
const Eigen::MatrixBase< DerivedF > & F )

#include <lagrange/eigen_convert.h>

Create a SurfaceMesh from a igl-style pair of matrices (V, F).

Note
The target Scalar x Index type must be explicitly specified. In the future, we may allow automatic deduction based on the input matrix types.
Parameters
[in]VV x d matrix of vertex positions.
[in]FF x k matrix of facet indices. E.g. k=3 for triangle meshes, k=4 for quad meshes.
Template Parameters
ScalarTarget mesh scalar type (required). Either float or double.
IndexTarget mesh index type (required). Either uint32_t or uint64_t.
DerivedVInput vertex matrix type.
DerivedFInput facet matrix type.
Returns
New mesh object.

◆ extract_boundary_loops()

template<typename Scalar, typename Index>
std::vector< std::vector< Index > > extract_boundary_loops ( const SurfaceMesh< Scalar, Index > & mesh)

#include <lagrange/extract_boundary_loops.h>

Extract boundary loops from a surface mesh.

Template Parameters
ScalarThe scalar type of the mesh.
IndexThe index type of the mesh.
Parameters
meshThe input surface mesh.
Returns
A vector of boundary loops, each represented by a vector of vertex indices.

◆ extract_submesh()

template<typename Scalar, typename Index>
SurfaceMesh< Scalar, Index > extract_submesh ( const SurfaceMesh< Scalar, Index > & mesh,
span< const Index > selected_facets,
const SubmeshOptions & options = {} )

#include <lagrange/extract_submesh.h>

Extract a submesh that consists of a subset of the facets of the source mesh.

Template Parameters
ScalarThe scalar type.
IndexThe index type.
Parameters
[in]meshThe source mesh.
[in]selected_facetsThe set of selected facets to extract.
[in]optionsExtraction options.
Returns
The mesh containing the selected facets.
See also
SubmeshOptions

◆ filtered_attribute_ids()

template<typename Scalar, typename Index>
std::vector< AttributeId > filtered_attribute_ids ( const SurfaceMesh< Scalar, Index > & mesh,
const AttributeFilter & options )

#include <lagrange/filter_attributes.h>

Create a list of attribute ids corresponding to the given filter.

Parameters
[in]meshMesh whose attributes are being filtered.
[in]optionsFilter options.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
A list of attribute ids matching the filter.

◆ filter_attributes()

template<typename Scalar, typename Index>
SurfaceMesh< Scalar, Index > filter_attributes ( SurfaceMesh< Scalar, Index > source_mesh,
const AttributeFilter & options = {} )

#include <lagrange/filter_attributes.h>

Filters the attributes of mesh according to user specifications.

Note
If the filter option does not contains AttributeElement::Edge as one of its element type, mesh edge information will be removed in the output mesh.
To convert a mesh and its attributes to different types, use the cast function.
Parameters
[in]source_meshInput mesh.
[in]optionsFilter options.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
Output mesh.
See also
cast

◆ find_matching_attribute()

template<typename Scalar, typename Index>
std::optional< AttributeId > find_matching_attribute ( const SurfaceMesh< Scalar, Index > & mesh,
const AttributeMatcher & options )

#include <lagrange/find_matching_attributes.h>

Finds the first attribute with the specified usage/element type/number of channels.

Parameters
[in]meshMesh whose attribute to retrieve.
[in]optionsAttribute properties to match against.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
The attribute id of the first matching attribute.

◆ find_matching_attributes() [1/3]

template<typename Scalar, typename Index>
std::vector< AttributeId > find_matching_attributes ( const SurfaceMesh< Scalar, Index > & mesh,
const AttributeMatcher & options )

#include <lagrange/find_matching_attributes.h>

Finds all attributes with the specified usage/element type/number of channels.

Parameters
[in]meshMesh whose attribute to retrieve.
[in]optionsAttribute properties to match against.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
A list of attribute ids of each matching attribute.

◆ find_matching_attributes() [2/3]

template<typename Scalar, typename Index>
std::vector< AttributeId > find_matching_attributes ( const SurfaceMesh< Scalar, Index > & mesh,
AttributeUsage usage )

#include <lagrange/find_matching_attributes.h>

Finds all attributes with the specified usage.

Parameters
[in]meshMesh whose attribute to retrieve.
[in]usageAttribute usage to match against.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
A list of attribute ids of each matching attribute.

◆ find_matching_attributes() [3/3]

template<typename Scalar, typename Index>
std::vector< AttributeId > find_matching_attributes ( const SurfaceMesh< Scalar, Index > & mesh,
BitField< AttributeElement > element_types )

#include <lagrange/find_matching_attributes.h>

Finds all attributes with the specified element types.

Parameters
[in]meshMesh whose attribute to retrieve.
[in]element_typesElement types to match against.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
A list of attribute ids of each matching attribute.

◆ trim_by_isoline()

template<typename Scalar, typename Index>
SurfaceMesh< Scalar, Index > trim_by_isoline ( const SurfaceMesh< Scalar, Index > & mesh,
const IsolineOptions & options = {} )

#include <lagrange/isoline.h>

Trim a mesh by the isoline of an implicit function defined on the mesh vertices/corners.

Note
The input must be a triangle mesh.
Parameters
[in]meshInput triangle mesh to trim.
[in]optionsIsoline options.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
The trimmed mesh.

◆ extract_isoline()

template<typename Scalar, typename Index>
SurfaceMesh< Scalar, Index > extract_isoline ( const SurfaceMesh< Scalar, Index > & mesh,
const IsolineOptions & options = {} )

#include <lagrange/isoline.h>

Extract the isoline of an implicit function defined on the mesh vertices/corners.

Note
The input must be a triangle mesh.
Parameters
[in]meshInput triangle mesh to extract the isoline from.
[in]optionsIsoline options.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
A mesh whose facets is a collection of size 2 elements representing the extracted isoline.

◆ mesh_bbox()

template<size_t Dimension, typename Scalar, typename Index>
Eigen::AlignedBox< Scalar, static_cast< int >(Dimension)> mesh_bbox ( const SurfaceMesh< Scalar, Index > & mesh)

#include <lagrange/mesh_bbox.h>

Compute the axis-aligned bounding box of a mesh.

If the mesh has no vertices, the returned bounding box is empty (default-constructed Eigen::AlignedBox, where isEmpty() returns true).

Parameters
[in]meshInput mesh. Its dimension must match the Dimension template parameter.
Template Parameters
DimensionSpatial dimension of the bounding box. Must be 2 or 3.
ScalarMesh scalar type.
IndexMesh index type.
Returns
The axis-aligned bounding box of the mesh vertices.

◆ to_surface_mesh_copy()

template<typename Scalar, typename Index, typename MeshType>
SurfaceMesh< Scalar, Index > to_surface_mesh_copy ( const MeshType & mesh)

#include <lagrange/mesh_convert.h>

Convert a legacy mesh object to a surface mesh object.

Parameters
[in]meshMesh object to convert.
Template Parameters
ScalarOutput mesh scalar type. Must be either float or double.
IndexOutput mesh index type. Must be either uint32_t or uint64_t.
MeshTypeInput mesh type.
Returns
New mesh object.

◆ to_surface_mesh_wrap()

template<typename Scalar, typename Index, typename MeshType>
SurfaceMesh< Scalar, Index > to_surface_mesh_wrap ( MeshType && mesh)

#include <lagrange/mesh_convert.h>

Wrap a legacy mesh object as a surface mesh object.

The mesh scalar & index types must match.

Parameters
[in]meshMesh object to convert. The mesh object must be a lvalue reference (no temporary).
Template Parameters
ScalarOutput mesh scalar type. Must be either float or double.
IndexOutput mesh index type. Must be either uint32_t or uint64_t.
MeshTypeInput mesh type.
Returns
New mesh object.

◆ to_legacy_mesh()

template<typename MeshType, typename Scalar, typename Index>
std::unique_ptr< MeshType > to_legacy_mesh ( const SurfaceMesh< Scalar, Index > & mesh)

#include <lagrange/mesh_convert.h>

Convert a surface mesh object to a legacy mesh object.

The mesh must be a regular mesh object.

Parameters
[in]meshMesh object to convert.
Template Parameters
MeshTypeOutput mesh type.
ScalarInput mesh scalar type. Must be either float or double.
IndexInput mesh index type. Must be either uint32_t or uint64_t.
Returns
New mesh object.

◆ normalize_mesh_with_transform()

template<size_t Dimension = 3, typename Scalar, typename Index>
auto normalize_mesh_with_transform ( SurfaceMesh< Scalar, Index > & mesh,
const TransformOptions & options = {} ) -> Eigen::Transform<Scalar, Dimension, Eigen::Affine>

#include <lagrange/normalize_meshes.h>

Normalize a mesh to fit in a unit box centered at the origin.

Parameters
[out]meshInput mesh.
[in]optionsTransform options.
Template Parameters
DimensionMesh dimension.
ScalarMesh scalar type.
IndexMesh index type.
Returns
The inverse transform, can be used to undo the normalization process.

◆ normalize_mesh()

template<typename Scalar, typename Index>
void normalize_mesh ( SurfaceMesh< Scalar, Index > & mesh,
const TransformOptions & options = {} )

#include <lagrange/normalize_meshes.h>

Normalize a mesh to fit in a unit box centered at the origin.

Parameters
[out]meshInput mesh.
[in]optionsTransform options.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.

◆ normalize_meshes_with_transform()

template<size_t Dimension = 3, typename Scalar, typename Index>
auto normalize_meshes_with_transform ( span< SurfaceMesh< Scalar, Index > * > meshes,
const TransformOptions & options = {} ) -> Eigen::Transform<Scalar, Dimension, Eigen::Affine>

#include <lagrange/normalize_meshes.h>

Normalize a list of meshes to fit in a unit box centered at the origin.

Parameters
[out]meshesList of pointers to the meshes to modify.
[in]optionsTransform options.
Template Parameters
DimensionMesh dimension.
ScalarMesh scalar type.
IndexMesh index type.
Returns
The inverse transform, can be used to undo the normalization process.

◆ normalize_meshes()

template<typename Scalar, typename Index>
void normalize_meshes ( span< SurfaceMesh< Scalar, Index > * > meshes,
const TransformOptions & options = {} )

#include <lagrange/normalize_meshes.h>

Normalize a list of meshes to fit in a unit box centered at the origin.

Parameters
[out]meshesList of pointers to the meshes to modify.
[in]optionsTransform options.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.

◆ orient_outward()

template<typename Scalar, typename Index>
void orient_outward ( lagrange::SurfaceMesh< Scalar, Index > & mesh,
const OrientOptions & options = {} )

#include <lagrange/orient_outward.h>

Orient the facets of a mesh so that the signed volume of each connected component is positive or negative.

Parameters
[in,out]meshMesh whose facets needs to be re-oriented.
[in]optionsOrientation options.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.

◆ is_oriented()

template<typename Scalar, typename Index>
bool is_oriented ( const SurfaceMesh< Scalar, Index > & mesh)

#include <lagrange/orientation.h>

Check if a mesh is oriented.

A mesh is oriented if interior edges have the same number of half-edges for each of the edge direction.

Parameters
meshThe input mesh.
Template Parameters
ScalarThe scalar type of the mesh.
IndexThe index type of the mesh.
Returns
True if the mesh is oriented.

< initial value of the result.

◆ compute_edge_is_oriented()

template<typename Scalar, typename Index>
AttributeId compute_edge_is_oriented ( SurfaceMesh< Scalar, Index > & mesh,
const OrientationOptions & options = {} )

#include <lagrange/orientation.h>

Compute a mesh attribute indicating whether an edge is oriented.

Parameters
[in,out]meshInput mesh.
[in]optionsOutput attribute options.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
Id of the newly added per-edge attribute, of type uint8_t.

◆ permute_facets()

template<typename Scalar, typename Index>
void permute_facets ( SurfaceMesh< Scalar, Index > & mesh,
span< const Index > new_to_old )

#include <lagrange/permute_facets.h>

Reorder facets of a mesh based on a given permutation.

i.e. rearrangement facets so that they are ordered as specified by the new_to_old index array. The total number of facets is unchanged.

Parameters
[in,out]meshThe target mesh whose facets will be reordered in place.
[in]new_to_oldThe permutation index array specifying the new facet order. This array can often be obtained via index-based sorting of the facets with customized comparison.
See also
extract_submesh to extract a subset of the facets.

◆ permute_vertices()

template<typename Scalar, typename Index>
void permute_vertices ( SurfaceMesh< Scalar, Index > & mesh,
span< const Index > new_to_old )

#include <lagrange/permute_vertices.h>

Reorder vertices of a mesh based on a given permutation.

i.e. rearrangement vertices so that they are ordered as specified by the new_to_old index array. The total number of vertices is unchanged.

Parameters
[in,out]meshThe target mesh whose vertices will be reordered in place.
[in]new_to_oldThe permutation index array specifying the new vertex order. This array can often be obtained via index-based sorting of the vertices with customized comparison.
See also
remap_vertices if two or more vertices may be combined.

◆ remap_vertices()

template<typename Scalar, typename Index>
void remap_vertices ( SurfaceMesh< Scalar, Index > & mesh,
span< const Index > forward_mapping,
RemapVerticesOptions options = {} )

#include <lagrange/remap_vertices.h>

Remap vertices of a mesh based on provided forward mapping.

Parameters
[in,out]meshThe target mesh.
[in]forward_mappingVertex mapping where vertex i will be remapped to vertex forward_mapping[i].
Precondition
  • forward_mapping must be surjective.
  • Edge information cannot be updated, thus its presence will cause an exception.
Postcondition
  • All vertex attributes will be updated.
  • The order of facets are unchanged.
  • If two vertices are mapped to the same index, they will be merged based on the collision policy specified in options.
See also
permute_vertices for simply permuting the vertex order.
RemapVerticesOptions

◆ reorder_mesh()

template<typename Scalar, typename Index>
void reorder_mesh ( SurfaceMesh< Scalar, Index > & mesh,
ReorderingMethod method )

#include <lagrange/reorder_mesh.h>

Mesh reordering to improve cache locality.

Reorder mesh vertices using Morton encoding.

The reordering is done in place.

Parameters
[in,out]meshMesh to reorder in place. For now we only support triangle meshes.
[in]methodReordering method.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Todo
Reorder mesh facets as well.
Parameters
[in,out]meshMesh to reorder.
Template Parameters
MeshTypeMesh type.

◆ segment_segment_squared_distance()

template<typename PointType>
auto segment_segment_squared_distance ( const Eigen::MatrixBase< PointType > & U0,
const Eigen::MatrixBase< PointType > & U1,
const Eigen::MatrixBase< PointType > & V0,
const Eigen::MatrixBase< PointType > & V1,
Eigen::PlainObjectBase< PointType > & closest_pointU,
Eigen::PlainObjectBase< PointType > & closest_pointV,
ScalarOf< PointType > & lambdaU,
ScalarOf< PointType > & lambdaV ) -> ScalarOf<PointType>

#include <lagrange/segment_segment_squared_distance.h>

Computes the squared distance between two N-d line segments, and the closest pair of points whose separation is this distance.

Parameters
[in]U0first extremity of the first segment
[in]U1second extremity of the first segment
[in]V0first extremity of the second segment
[in]V1second extremity of the second segment
[out]closest_pointUthe closest point on segment [U0, U1]
[out]closest_pointVthe closest point on segment [V0, V1]
[out]lambdaUbarycentric coordinate of the closest point relative to [U0, U1]
[out]lambdaVbarycentric coordinate of the closest point relative to [V0, V1]
Template Parameters
PointTypethe class that represents the points.
Returns
the squared distance between the segments [U0, U1] and [V0, V1]
Note
Adapted from Real-Time Collision Detection by Christer Ericson, published by Morgan Kaufmann Publishers, (c) 2005 Elsevier Inc. This function was modified for use by Siddhartha Chaudhuri in the Thea library, from where this version was taken on 28 Sep 2022. All modifications beyond Thea are Copyright 2022 Adobe.

◆ select_facets_by_normal_similarity()

template<typename Scalar, typename Index>
AttributeId select_facets_by_normal_similarity ( SurfaceMesh< Scalar, Index > & mesh,
const Index seed_facet_id,
const SelectFacetsByNormalSimilarityOptions & options = {} )

#include <lagrange/select_facets_by_normal_similarity.h>

Given a seed facet, selects facets around it based on the change in triangle normals.

Parameters
[in,out]meshThe input mesh.
[in]seed_facet_idThe index of the seed facet
[in]optionsOptional arguments (greedy, output_attribute_name).
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
AttributeId The attribute id of the selection, whose attribute name is given by options.output_attribute_name.
Note
Currently only support triangular mesh and will throw error if mesh.is_triangle_mesh() is false! The function will check if the mesh contains facet normal by looking for options.facet_normal_attribute_name, and if not found, will call compute_facet_normal(meshm, options.facet_normal_attribute_name).
See also
SelectFacetsByNormalSimilarityOptions

◆ select_facets_in_frustum()

template<typename Scalar, typename Index>
bool select_facets_in_frustum ( SurfaceMesh< Scalar, Index > & mesh,
const Frustum< Scalar > & frustum,
const FrustumSelectionOptions & options = {} )

#include <lagrange/select_facets_in_frustum.h>

Select all facets that intersect the cone/frustrum bounded by 4 planes defined by (n_i, p_i), where n_i is the plane normal and p_i is a point on the plane.

Parameters
[in,out]meshThe input mesh.
[in]frustumA collection of four planes.
[in]optionsOptional arguments (greedy, output_attribute_name).
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
bool Whether any facet is selected.
Note
When options.greedy is true, this function returns as soon as the first facet is selected.
Postcondition
If options.greedy is false, the computed selection is stored in mesh as a facet attribute named options.output_attribute_name.
See also
Frustum and FrustumSelectionOptions.

◆ separate_by_components()

template<typename Scalar, typename Index>
std::vector< SurfaceMesh< Scalar, Index > > separate_by_components ( const SurfaceMesh< Scalar, Index > & mesh,
const SeparateByComponentsOptions & options = {} )

#include <lagrange/separate_by_components.h>

Separate a mesh by connected components.

Template Parameters
ScalarThe scalar type.
IndexThe index type.
Parameters
[in]meshThe source mesh.
[in]optionsOption settings.
Returns
A list of meshes representing the set of connected components.
See also
SeparateByComponentsOptions

◆ separate_by_facet_groups() [1/3]

template<typename Scalar, typename Index>
std::vector< SurfaceMesh< Scalar, Index > > separate_by_facet_groups ( const SurfaceMesh< Scalar, Index > & mesh,
size_t num_groups,
span< const Index > facet_group_indices,
const SeparateByFacetGroupsOptions & options = {} )

#include <lagrange/separate_by_facet_groups.h>

Extract a set of submeshes based on facet groups.

Facets with the same group index are grouped together in a single submesh.

Template Parameters
ScalarThe scalar type.
IndexThe index type.
Parameters
[in]meshThe source mesh.
[in]num_groupsThe number of face groups.
[in]facet_group_indicesThe group index of each facet. Each group index must be in the range of [0, num_groups - 1].
[in]optionsExtraction options.
Returns
A list of submeshes representing each facet group.

◆ separate_by_facet_groups() [2/3]

template<typename Scalar, typename Index>
std::vector< SurfaceMesh< Scalar, Index > > separate_by_facet_groups ( const SurfaceMesh< Scalar, Index > & mesh,
span< const Index > facet_group_indices,
const SeparateByFacetGroupsOptions & options = {} )

#include <lagrange/separate_by_facet_groups.h>

Extract a set of submeshes based on facet groups.

Facets with the same group index are grouped together in a single submesh.

Template Parameters
ScalarThe scalar type.
IndexThe index type.
Parameters
[in]meshThe source mesh.
[in]facet_group_indicesThe group index of each facet. Each group index must be in the range of [0, max(facet_group_indices)].
[in]optionsExtraction options.
Returns
A list of submeshes representing each facet group.

◆ separate_by_facet_groups() [3/3]

template<typename Scalar, typename Index>
std::vector< SurfaceMesh< Scalar, Index > > separate_by_facet_groups ( const SurfaceMesh< Scalar, Index > & mesh,
size_t num_groups,
function_ref< Index(Index)> get_facet_group,
const SeparateByFacetGroupsOptions & options = {} )

#include <lagrange/separate_by_facet_groups.h>

Extract a set of submeshes based on facet groups.

Facets with the same group index are grouped together in a single submesh.

Template Parameters
ScalarThe scalar type.
IndexThe index type.
Parameters
[in]meshThe source mesh.
[in]num_groupsThe number of face groups.
[in]get_facet_groupFunction that returns the facet group id from facet id. The groud id must be in [0, num_groups - 1].
[in]optionsExtraction options.
Returns
A list of submeshes representing each facet group.

◆ split_facets_by_material()

template<typename Scalar, typename Index>
void split_facets_by_material ( SurfaceMesh< Scalar, Index > & mesh,
std::string_view material_attribute_name )

#include <lagrange/split_facets_by_material.h>

Split mesh facets based on material labels.

Parameters
meshMesh on which material segmentation will be applied in place.
material_attribute_nameName of the material attribute to use for splitting.
Note
The attribute should be a vertex attribute with k channels, each channel encodes the probability of the vertex being part of material i (0 <= i < k). The most probable material will be assigned to each vertex, and facets will be split at the boundaries of different materials.
Template Parameters
ScalarMesh scalar type.
Parameters
IndexMesh index type.

◆ thicken_and_close_mesh()

template<typename Scalar, typename Index>
SurfaceMesh< Scalar, Index > thicken_and_close_mesh ( SurfaceMesh< Scalar, Index > input_mesh,
const ThickenAndCloseOptions & options = {} )

#include <lagrange/thicken_and_close_mesh.h>

Thicken a mesh by offsetting it, and close the shape into a thick 3D solid.

Input mesh vertices are duplicated and can additionally be mirrored wrt to a plane going through the origin.

Parameters
[in]input_meshInput mesh, assumed to have a disk topology. Must have edge information included.
[in]optionsThickening options.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
A mesh thickened into a 3D solid shape.

◆ compute_euler()

template<typename Scalar, typename Index>
int compute_euler ( const SurfaceMesh< Scalar, Index > & mesh)

#include <lagrange/topology.h>

Compute Euler characteristic of a mesh.

Template Parameters
ScalarThe scalar type of the mesh.
IndexThe index type of the mesh.
Parameters
meshThe input mesh.
Returns
The Euler characteristic of the mesh.

◆ is_closed()

template<typename Scalar, typename Index>
bool is_closed ( const SurfaceMesh< Scalar, Index > & mesh)

#include <lagrange/topology.h>

Check if a mesh is closed.

A mesh is closed if and only if it has no boundary edges.

Template Parameters
ScalarThe scalar type of the mesh.
IndexThe index type of the mesh.
Parameters
meshThe input mesh.
Returns
True if the mesh is closed.

◆ is_vertex_manifold()

template<typename Scalar, typename Index>
bool is_vertex_manifold ( const SurfaceMesh< Scalar, Index > & mesh)

#include <lagrange/topology.h>

Check if a mesh is vertex-manifold.

A mesh is vertex-manifold if and only if the one-ring neighborhood of each vertex is of disc topology. I.e. The boundary of the 1-ring neighborhood is a simple loop for interior vertices, and a simple chain for boundary vertices.

Template Parameters
ScalarThe scalar type of the mesh.
IndexThe index type of the mesh.
Parameters
meshThe input mesh.
Returns
True if the mesh is vertex-manifold.

< initial value of the result.

◆ is_edge_manifold()

template<typename Scalar, typename Index>
bool is_edge_manifold ( const SurfaceMesh< Scalar, Index > & mesh)

#include <lagrange/topology.h>

Check if a mesh is edge-manifold.

A mesh is edge-manifold if and only if every interior edge is incident to exactly two facets, and every boundary edge is incident to exactly one facet.

Template Parameters
ScalarThe scalar type of the mesh.
IndexThe index type of the mesh.
Parameters
meshThe input mesh.
Returns
True if the mesh is edge-manifold.

◆ is_manifold()

template<typename Scalar, typename Index>
bool is_manifold ( const SurfaceMesh< Scalar, Index > & mesh)

#include <lagrange/topology.h>

Check if a mesh is both vertex-manifold and edge-manifold.

Template Parameters
ScalarThe scalar type of the mesh.
IndexThe index type of the mesh.
Parameters
meshThe input mesh.
Returns
True if the mesh is both vertex-manifold and edge-manifold.

◆ compute_vertex_is_manifold()

template<typename Scalar, typename Index>
AttributeId compute_vertex_is_manifold ( SurfaceMesh< Scalar, Index > & mesh,
const VertexManifoldOptions & options = {} )

#include <lagrange/topology.h>

Compute a mesh attribute of value type uint8_t indicating vertex manifoldness.

Parameters
[in,out]meshInput mesh.
[in]optionsOutput attribute options.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
Id of the newly added per-vertex attribute.

◆ compute_edge_is_manifold()

template<typename Scalar, typename Index>
AttributeId compute_edge_is_manifold ( SurfaceMesh< Scalar, Index > & mesh,
const EdgeManifoldOptions & options = {} )

#include <lagrange/topology.h>

Compute a mesh attribute of value type uint8_t indicating edge manifoldness.

Parameters
[in,out]meshInput mesh.
[in]optionsOutput attribute options.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
Returns
Id of the newly added per-edge attribute.

◆ transform_mesh()

template<typename Scalar, typename Index, int Dimension>
void transform_mesh ( SurfaceMesh< Scalar, Index > & mesh,
const Eigen::Transform< Scalar, Dimension, Eigen::Affine > & transform,
const TransformOptions & options = {} )

#include <lagrange/transform_mesh.h>

Apply an affine transform \( M \) to a mesh in-place.

All mesh attributes are transformed based on their usage tags:

  • Position: Applies \( P \to M * P \)
  • Normal: Applies \( P \to \det(M) M^{-T} * P \)
  • Tangent: Applies \( P \to normalize(M * P) \)
  • Bitangent: Applies \( P \to normalize(M * P) \)

    Todo
    Add an overload for 2D transforms.
Parameters
[in,out]meshMesh to transform in-place.
[in]transformAffine transform to apply.
[in]optionsTransform options.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
DimensionTransform dimension (either 2 or 3). Must match mesh dimension.

◆ transformed_mesh()

template<typename Scalar, typename Index, int Dimension>
SurfaceMesh< Scalar, Index > transformed_mesh ( SurfaceMesh< Scalar, Index > mesh,
const Eigen::Transform< Scalar, Dimension, Eigen::Affine > & transform,
const TransformOptions & options = {} )

#include <lagrange/transform_mesh.h>

Apply an affine transform to a mesh and return the transformed mesh.

All mesh attributes are transformed based on their usage tags:

  • Position: Applies \( P \to M * P \)
  • Normal: Applies \( P \to \det(M) M^{-T} * P \)
  • Tangent: Applies \( P \to normalize(M * P) \)
  • Bitangent: Applies \( P \to normalize(M * P) \)

    Todo
    Add an overload for 2D transforms.
Parameters
[in]meshMesh to transform in-place.
[in]transformAffine transform to apply.
[in]optionsTransform options.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
DimensionTransform dimension (either 2 or 3). Must match mesh dimension.
Returns
Transformed mesh.

◆ triangulate_polygonal_facets()

template<typename Scalar, typename Index>
void triangulate_polygonal_facets ( SurfaceMesh< Scalar, Index > & mesh,
const TriangulationOptions & options = {} )

#include <lagrange/triangulate_polygonal_facets.h>

Triangulate polygonal facets of a mesh using a prescribed set of rules.

Parameters
[in,out]meshPolygonal mesh to triangulate in place.
[in]optionsOptions for triangulation.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.

◆ uv_mesh_ref()

template<typename Scalar, typename Index, typename UVScalar = Scalar>
SurfaceMesh< UVScalar, Index > uv_mesh_ref ( SurfaceMesh< Scalar, Index > & mesh,
const UVMeshOptions & options = {} )

#include <lagrange/uv_mesh.h>

Extract a UV mesh reference from an input mesh.

This method will create a new mesh with the same topology as the input mesh, but with vertex positions set to the corresponding UV coordinates. Modification of UV mesh vertices will be reflected in the input mesh.

Parameters
meshInput mesh.
optionsOptions to control UV mesh extraction.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
UVScalarTarget UV attribute value type.
Returns
The extracted UV mesh reference.
See also
UVMeshOptions

◆ uv_mesh_view()

template<typename Scalar, typename Index, typename UVScalar = Scalar>
SurfaceMesh< UVScalar, Index > uv_mesh_view ( const SurfaceMesh< Scalar, Index > & mesh,
const UVMeshOptions & options = {} )

#include <lagrange/uv_mesh.h>

Extract a UV mesh view from an input mesh.

This method will create a new mesh with the same topology as the input mesh, but with vertex positions set to the corresponding UV coordinates. The output UV mesh cannot be modified.

Parameters
meshInput mesh.
optionsOptions to control UV mesh extraction.
Template Parameters
ScalarMesh scalar type.
IndexMesh index type.
UVScalarTarget UV attribute value type.
Returns
The extracted UV mesh reference.
See also
UVMeshOptions