Skip to content

shapeworks::Mesh

Module: Mesh Classes

More...

#include <Mesh.h>

Public Types

Name
enum FieldType
enum AlignmentType
enum DistanceMethod
enum CurvatureType
enum SubdivisionType
using vtkSmartPointer< vtkPolyData > MeshType
using vtkSmartPointer< vtkPoints > MeshPoints

Public Functions

Name
Mesh(const std::string & pathname)
Mesh(MeshType meshPtr)
Mesh(const Mesh & orig)
Mesh(Mesh && orig)
Mesh & operator=(const Mesh & orig)
Mesh(const Eigen::MatrixXd & points, const Eigen::MatrixXi & faces)
Mesh & operator=(Mesh && orig)
Mesh & operator+=(const Mesh & otherMesh)
append two meshes
MeshType getVTKMesh() const
return the current mesh
Mesh & write(const std::string & pathname, bool binaryFile =false)
writes mesh, format specified by filename extension
Mesh & coverage(const Mesh & otherMesh, bool allowBackIntersections =true, double angleThreshold =0, double backSearchRadius =0)
determines coverage between current mesh and another mesh (e.g. acetabular cup / femoral head)
Mesh & smooth(int iterations =0, double relaxation =0.0)
applies laplacian smoothing
Mesh & smoothSinc(int iterations =0, double passband =0.0)
applies vtk windowed sinc smoothing
Mesh & remesh(int numVertices, double adaptivity =1.0)
applies remeshing using approximated centroidal voronoi diagrams for a given number of vertices and adaptivity
Mesh & remeshPercent(double percentage, double adaptivity =1.0)
applies remeshing using approximated centroidal voronoi diagrams for a given percentage of vertices and adaptivity
Mesh & invertNormals()
handle flipping normals
Mesh & reflect(const Axis & axis, const Vector3 & origin =makeVector({0.0, 0.0, 0.0}))
reflect meshes with respect to a specified center and specific axis
MeshTransform createTransform(const Mesh & target, AlignmentType align =Similarity, unsigned iterations =10)
Mesh & applyTransform(const MeshTransform transform)
applies the given transformation to the mesh
Mesh & fillHoles()
finds holes in a mesh and closes them
Mesh & clean()
clean mesh
Mesh & probeVolume(const Image & image)
samples image data values at point locations specified by image
Mesh & clip(const Plane plane)
clips a mesh using a cutting plane
Mesh & translate(const Vector3 & v)
helper to translate mesh
Mesh & scale(const Vector3 & v)
helper to scale mesh
PhysicalRegion boundingBox() const
computes bounding box of current mesh
Mesh & fixElement()
fix element winding of mesh
std::vector< Field > distance(const Mesh & target, const DistanceMethod method =PointToCell) const
Mesh & clipClosedSurface(const Plane plane)
clips a mesh using a cutting plane resulting in a closed surface
Mesh & computeNormals()
computes and adds oriented point and cell normals
Point3 closestPoint(const Point3 point, bool & outside, double & distance, vtkIdType & face_id) const
int closestPointId(const Point3 point) const
returns closest point id in this mesh to the given point in space
double geodesicDistance(int source, int target) const
computes geodesic distance between two vertices (specified by their indices) on mesh
Field geodesicDistance(const Point3 landmark) const
computes geodesic distance between a point (landmark) and each vertex on mesh
Field geodesicDistance(const std::vector< Point3 > curve) const
computes geodesic distance between a set of points (curve) and each vertex on mesh
Field curvature(const CurvatureType type =Principal) const
computes curvature using principal (default) or gaussian or mean algorithms
Mesh & applySubdivisionFilter(const SubdivisionType type =Butterfly, int subdivision =1)
applies subdivision filter (butterfly (default) or loop)
Image toImage(PhysicalRegion region =PhysicalRegion(), Point3 spacing =Point3({1., 1., 1.})) const
rasterizes specified region to create binary image of desired dims (default: unit spacing)
Image toDistanceTransform(PhysicalRegion region =PhysicalRegion(), const Point3 spacing =Point3({1., 1., 1.}), const Dims padding =Dims({1, 1, 1})) const
converts specified region to distance transform image (default: unit spacing) with (logical) padding
Point3 center() const
center of mesh
Point3 centerOfMass() const
center of mass of mesh
int numPoints() const
number of points
int numFaces() const
number of faces
Eigen::MatrixXd points() const
matrix with number of points with (x,y,z) coordinates of each point
Eigen::MatrixXi faces() const
matrix with number of faces with indices of the three points from which each face is composed
Point3 getPoint(int id) const
(x,y,z) coordinates of vertex at given index
IPoint3 getFace(int id) const
return indices of the three points with which the face at the given index is composed
std::vector< std::string > getFieldNames() const
print all field names in mesh
Mesh & setField(const std::string name, Array array, const FieldType type)
sets the given field for points or faces with array (*does not copy array's values)
Field getField(const std::string & name, const FieldType type) const
gets a pointer to the requested field of points or faces, null if field doesn't exist
void setFieldValue(const std::string & name, int idx, double value)
sets the given index of field to value
double getFieldValue(const std::string & name, int idx) const
gets the value at the given index of field (NOTE: returns first component of vector fields)
Eigen::VectorXd getMultiFieldValue(const std::string & name, int idx) const
gets the multi value at the given index of [vertex] field
bool compareAllPoints(const Mesh & other_mesh) const
compare if values of the points in two (corresponding) meshes are (eps)equal
bool compareAllFaces(const Mesh & other_mesh) const
compare if face indices in two (corresponding) meshes are equal
bool compareAllFields(const Mesh & other_mesh, const double eps =-1.0) const
compare if all fields in two meshes are (eps)equal
bool compareField(const Mesh & other_mesh, const std::string & name1, const std::string & name2 ="", const double eps =-1.0) const
compare field of meshes to be (eps)equal (same field for both if only one specified)
bool compare(const Mesh & other_mesh, const double eps =-1.0) const
compare meshes
bool operator==(const Mesh & other) const
compare meshes
bool prepareFFCFields(std::vector< std::vector< Eigen::Vector3d >> boundaries, Eigen::Vector3d query, bool onlyGenerateInOut =false)
Prepares the mesh for FFCs by setting scalar and vector fields.
double getFFCValue(Eigen::Vector3d query) const
Gets values for FFCs.
Eigen::Vector3d getFFCGradient(Eigen::Vector3d query) const
Gets gradients for FFCs.
MeshPoints getIGLMesh(Eigen::MatrixXd & V, Eigen::MatrixXi & F) const
Formats mesh into an IGL format.
vtkSmartPointer< vtkPolyData > clipByField(const std::string & name, double value)
Clips the mesh according to a field value.
std::vector< std::string > getSupportedTypes()
getSupportedTypes

Friends

Name
struct SharedCommandData

Detailed Description

class shapeworks::Mesh;

This class encapsulates a Mesh and operations that can be performed on meshes

Public Types Documentation

enum FieldType

Enumerator Value Description
Point
Face

enum AlignmentType

Enumerator Value Description
Rigid
Similarity
Affine

enum DistanceMethod

Enumerator Value Description
PointToPoint
PointToCell

enum CurvatureType

Enumerator Value Description
Principal
Gaussian
Mean

enum SubdivisionType

Enumerator Value Description
Butterfly
Loop

using MeshType

using shapeworks::Mesh::MeshType =  vtkSmartPointer<vtkPolyData>;

using MeshPoints

using shapeworks::Mesh::MeshPoints =  vtkSmartPointer<vtkPoints>;

Public Functions Documentation

function Mesh

Mesh(
    const std::string & pathname
)

function Mesh

inline Mesh(
    MeshType meshPtr
)

function Mesh

inline Mesh(
    const Mesh & orig
)

function Mesh

inline Mesh(
    Mesh && orig
)

function operator=

inline Mesh & operator=(
    const Mesh & orig
)

function Mesh

Mesh(
    const Eigen::MatrixXd & points,
    const Eigen::MatrixXi & faces
)

function operator=

inline Mesh & operator=(
    Mesh && orig
)

function operator+=

Mesh & operator+=(
    const Mesh & otherMesh
)

append two meshes

function getVTKMesh

inline MeshType getVTKMesh() const

return the current mesh

function write

Mesh & write(
    const std::string & pathname,
    bool binaryFile =false
)

writes mesh, format specified by filename extension

function coverage

Mesh & coverage(
    const Mesh & otherMesh,
    bool allowBackIntersections =true,
    double angleThreshold =0,
    double backSearchRadius =0
)

determines coverage between current mesh and another mesh (e.g. acetabular cup / femoral head)

function smooth

Mesh & smooth(
    int iterations =0,
    double relaxation =0.0
)

applies laplacian smoothing

function smoothSinc

Mesh & smoothSinc(
    int iterations =0,
    double passband =0.0
)

applies vtk windowed sinc smoothing

function remesh

Mesh & remesh(
    int numVertices,
    double adaptivity =1.0
)

applies remeshing using approximated centroidal voronoi diagrams for a given number of vertices and adaptivity

function remeshPercent

Mesh & remeshPercent(
    double percentage,
    double adaptivity =1.0
)

applies remeshing using approximated centroidal voronoi diagrams for a given percentage of vertices and adaptivity

function invertNormals

Mesh & invertNormals()

handle flipping normals

function reflect

Mesh & reflect(
    const Axis & axis,
    const Vector3 & origin =makeVector({0.0, 0.0, 0.0})
)

reflect meshes with respect to a specified center and specific axis

function createTransform

MeshTransform createTransform(
    const Mesh & target,
    AlignmentType align =Similarity,
    unsigned iterations =10
)

creates transform to target mesh using specified AlignmentType (Mesh::Rigid, Mesh::Similarity, Mesh::Affine) for specified number of iterations

function applyTransform

Mesh & applyTransform(
    const MeshTransform transform
)

applies the given transformation to the mesh

function fillHoles

Mesh & fillHoles()

finds holes in a mesh and closes them

function clean

Mesh & clean()

clean mesh

function probeVolume

Mesh & probeVolume(
    const Image & image
)

samples image data values at point locations specified by image

function clip

Mesh & clip(
    const Plane plane
)

clips a mesh using a cutting plane

function translate

Mesh & translate(
    const Vector3 & v
)

helper to translate mesh

function scale

Mesh & scale(
    const Vector3 & v
)

helper to scale mesh

function boundingBox

PhysicalRegion boundingBox() const

computes bounding box of current mesh

function fixElement

Mesh & fixElement()

fix element winding of mesh

function distance

std::vector< Field > distance(
    const Mesh & target,
    const DistanceMethod method =PointToCell
) const

Computes distance from each vertex to closest cell or point in target mesh, specified as PointToCell (default) or PointToPoint. Returns Fields containing distance to target and ids of the associated cells or points.

function clipClosedSurface

Mesh & clipClosedSurface(
    const Plane plane
)

clips a mesh using a cutting plane resulting in a closed surface

function computeNormals

Mesh & computeNormals()

computes and adds oriented point and cell normals

function closestPoint

Point3 closestPoint(
    const Point3 point,
    bool & outside,
    double & distance,
    vtkIdType & face_id
) const

Returns closest point on this mesh to the given point in space. In addition, returns by reference:

  • whether the point in space is outside the mesh or not
  • the distance of the point in space from this mesh
  • the face_id containing the closest point

function closestPointId

int closestPointId(
    const Point3 point
) const

returns closest point id in this mesh to the given point in space

function geodesicDistance

double geodesicDistance(
    int source,
    int target
) const

computes geodesic distance between two vertices (specified by their indices) on mesh

function geodesicDistance

Field geodesicDistance(
    const Point3 landmark
) const

computes geodesic distance between a point (landmark) and each vertex on mesh

function geodesicDistance

Field geodesicDistance(
    const std::vector< Point3 > curve
) const

computes geodesic distance between a set of points (curve) and each vertex on mesh

function curvature

Field curvature(
    const CurvatureType type =Principal
) const

computes curvature using principal (default) or gaussian or mean algorithms

function applySubdivisionFilter

Mesh & applySubdivisionFilter(
    const SubdivisionType type =Butterfly,
    int subdivision =1
)

applies subdivision filter (butterfly (default) or loop)

function toImage

Image toImage(
    PhysicalRegion region =PhysicalRegion(),
    Point3 spacing =Point3({1., 1., 1.})
) const

rasterizes specified region to create binary image of desired dims (default: unit spacing)

function toDistanceTransform

Image toDistanceTransform(
    PhysicalRegion region =PhysicalRegion(),
    const Point3 spacing =Point3({1., 1., 1.}),
    const Dims padding =Dims({1, 1, 1})
) const

converts specified region to distance transform image (default: unit spacing) with (logical) padding

function center

Point3 center() const

center of mesh

function centerOfMass

Point3 centerOfMass() const

center of mass of mesh

function numPoints

inline int numPoints() const

number of points

function numFaces

inline int numFaces() const

number of faces

function points

Eigen::MatrixXd points() const

matrix with number of points with (x,y,z) coordinates of each point

function faces

Eigen::MatrixXi faces() const

matrix with number of faces with indices of the three points from which each face is composed

function getPoint

Point3 getPoint(
    int id
) const

(x,y,z) coordinates of vertex at given index

function getFace

IPoint3 getFace(
    int id
) const

return indices of the three points with which the face at the given index is composed

function getFieldNames

std::vector< std::string > getFieldNames() const

print all field names in mesh

function setField

Mesh & setField(
    const std::string name,
    Array array,
    const FieldType type
)

sets the given field for points or faces with array (*does not copy array's values)

function getField

Field getField(
    const std::string & name,
    const FieldType type
) const

gets a pointer to the requested field of points or faces, null if field doesn't exist

function setFieldValue

void setFieldValue(
    const std::string & name,
    int idx,
    double value
)

sets the given index of field to value

function getFieldValue

double getFieldValue(
    const std::string & name,
    int idx
) const

gets the value at the given index of field (NOTE: returns first component of vector fields)

function getMultiFieldValue

Eigen::VectorXd getMultiFieldValue(
    const std::string & name,
    int idx
) const

gets the multi value at the given index of [vertex] field

function compareAllPoints

bool compareAllPoints(
    const Mesh & other_mesh
) const

compare if values of the points in two (corresponding) meshes are (eps)equal

function compareAllFaces

bool compareAllFaces(
    const Mesh & other_mesh
) const

compare if face indices in two (corresponding) meshes are equal

function compareAllFields

bool compareAllFields(
    const Mesh & other_mesh,
    const double eps =-1.0
) const

compare if all fields in two meshes are (eps)equal

function compareField

bool compareField(
    const Mesh & other_mesh,
    const std::string & name1,
    const std::string & name2 ="",
    const double eps =-1.0
) const

compare field of meshes to be (eps)equal (same field for both if only one specified)

function compare

bool compare(
    const Mesh & other_mesh,
    const double eps =-1.0
) const

compare meshes

function operator==

inline bool operator==(
    const Mesh & other
) const

compare meshes

function prepareFFCFields

bool prepareFFCFields(
    std::vector< std::vector< Eigen::Vector3d >> boundaries,
    Eigen::Vector3d query,
    bool onlyGenerateInOut =false
)

Prepares the mesh for FFCs by setting scalar and vector fields.

function getFFCValue

double getFFCValue(
    Eigen::Vector3d query
) const

Gets values for FFCs.

function getFFCGradient

Eigen::Vector3d getFFCGradient(
    Eigen::Vector3d query
) const

Gets gradients for FFCs.

function getIGLMesh

MeshPoints getIGLMesh(
    Eigen::MatrixXd & V,
    Eigen::MatrixXi & F
) const

Formats mesh into an IGL format.

function clipByField

vtkSmartPointer< vtkPolyData > clipByField(
    const std::string & name,
    double value
)

Clips the mesh according to a field value.

function getSupportedTypes

static inline std::vector< std::string > getSupportedTypes()

getSupportedTypes

Friends

friend SharedCommandData

friend struct SharedCommandData();

Updated on 2022-07-23 at 16:40:07 -0600