Euclid
Geometry Processing and Shape Analysis in C++
Euclid::GeodesicsInHeat< Mesh > Class Template Reference

Single source approximate geodesic distance using the heat method. More...

#include <GeodesicsInHeat.h>

Public Types

using FT = FT_t< Mesh >
 
using Vector_3 = Vector_3_t< Mesh >
 
using SpMat = Eigen::SparseMatrix< FT >
 

Public Member Functions

void build (const Mesh &mesh, float scale=1.0f, FT resolution=0, const SpMat *cot_mat=nullptr, const SpMat *mass_mat=nullptr)
 Build up the necesssary computational components. More...
 
void scale (float scale)
 Reset the time scale. More...
 
template<typename T >
void compute (const typename boost::graph_traits< const Mesh >::vertex_descriptor &v, std::vector< T > &geodesics)
 Compute geodesics distance from a vertex. More...
 

Public Attributes

const Mesh * mesh = nullptr
 The target mesh.
 
FT resolution = 0.0
 The resolution of mesh. More...
 
ProPtr< const SpMat > cot_mat = nullptr
 Cotangent matrix. More...
 
ProPtr< const SpMat > mass_mat = nullptr
 Mass matrix. More...
 
Eigen::SimplicialLDLT< SpMat > heat_solver
 The heat equation solver.
 
Eigen::SimplicialLDLT< SpMat > poisson_solver
 The poisson equation solver.
 

Detailed Description

template<typename Mesh>
class Euclid::GeodesicsInHeat< Mesh >

Reference

Crane K, Weischedel C, Wardetzky M. Geodesics in heat: A new approach to computing distance based on heat flow. ACM Transactions on Graphics (TOG), 2013.

Member Function Documentation

template<typename Mesh >
void Euclid::GeodesicsInHeat< Mesh >::build ( const Mesh &  mesh,
float  scale = 1.0f,
FT  resolution = 0,
const SpMat *  cot_mat = nullptr,
const SpMat *  mass_mat = nullptr 
)

Compute the Laplacian matrix and pre-factor the linear system.

Parameters
meshTarget triangle mesh.
scaleThe time scale of the heat diffusion, relative to the average edge length of the mesh.
resolutionThe average edge length of the mesh, provide a value if you have already computed it, otherwise it'll be computed internally.
cot_matThe cotangent matrix of the mesh, provide a value if you have already computed it, otherwise it'll be computed internally.
mass_matThe mass matrix of the mesh, provide a value if you have already computed it, otherwise it'll be computed internally.
template<typename Mesh >
template<typename T >
void Euclid::GeodesicsInHeat< Mesh >::compute ( const typename boost::graph_traits< const Mesh >::vertex_descriptor &  v,
std::vector< T > &  geodesics 
)
Parameters
vThe vertex descriptor.
geodesicsThe output geodesics distances from all the mesh vertices to the target vertex v.
template<typename Mesh >
void Euclid::GeodesicsInHeat< Mesh >::scale ( float  scale)
Parameters
scaleThe time scale of the heat diffusion, relative to the average edge length of the mesh.

Member Data Documentation

template<typename Mesh >
ProPtr<const SpMat> Euclid::GeodesicsInHeat< Mesh >::cot_mat = nullptr

If you directly set this member without calling the build, you need to update heat_solver and poisson_solver manually.

template<typename Mesh >
ProPtr<const SpMat> Euclid::GeodesicsInHeat< Mesh >::mass_mat = nullptr

If you directly set this member without calling the build, you need to update heat_solver and poisson_solver manually.

template<typename Mesh >
FT Euclid::GeodesicsInHeat< Mesh >::resolution = 0.0

Average edge length.


The documentation for this class was generated from the following files: