Loading...
Searching...
No Matches
Gudhi::coxeter_triangulation::Coxeter_triangulation< Permutahedral_representation_ > Class Template Reference

A class that stores Coxeter triangulation of type \(\tilde{A}_d\). This triangulation has the greatest simplex quality out of all linear transformations of the Freudenthal-Kuhn triangulation. More...

#include <include/gudhi/Coxeter_triangulation.h>

Public Member Functions

 Coxeter_triangulation (std::size_t dimension)
 Constructor of Coxeter triangulation of a given dimension. More...
 
- Public Member Functions inherited from Gudhi::coxeter_triangulation::Freudenthal_triangulation< Permutahedral_representation_ >
 Freudenthal_triangulation (std::size_t dimension)
 Constructor of the Freudenthal-Kuhn triangulation of a given dimension. More...
 
 Freudenthal_triangulation (std::size_t dimension, const Matrix &matrix)
 Constructor of the Freudenthal-Kuhn triangulation of a given dimension under a linear transformation by a given matrix. More...
 
 Freudenthal_triangulation (unsigned dimension, const Matrix &matrix, const Vector &offset)
 Constructor of the Freudenthal-Kuhn triangulation of a given dimension under an affine transformation by a given matrix and a translation vector. More...
 
unsigned dimension () const
 Dimension of the triangulation.
 
const Matrix & matrix () const
 Matrix that defines the linear transformation of the triangulation.
 
const Vector & offset () const
 Vector that defines the offset of the triangulation.
 
void change_matrix (const Eigen::MatrixXd &matrix)
 Change the linear transformation matrix to a given value. More...
 
void change_offset (const Eigen::VectorXd &offset)
 Change the offset vector to a given value. More...
 
template<class Point_d >
Simplex_handle locate_point (const Point_d &point, double scale=1) const
 Returns the permutahedral representation of the simplex in the triangulation that contains a given query point. More...
 
Eigen::VectorXd cartesian_coordinates (const Vertex_handle &vertex, double scale=1) const
 Returns the Cartesian coordinates of the given vertex. More...
 
Eigen::VectorXd barycenter (const Simplex_handle &simplex, double scale=1) const
 Returns the Cartesian coordinates of the barycenter of a given simplex. More...
 

Additional Inherited Members

- Public Types inherited from Gudhi::coxeter_triangulation::Freudenthal_triangulation< Permutahedral_representation_ >
using Simplex_handle = Permutahedral_representation_
 Type of the simplices in the triangulation.
 
using Vertex_handle = typename Permutahedral_representation_::Vertex
 Type of the vertices in the triangulation.
 

Detailed Description

template<class Permutahedral_representation_ = Permutahedral_representation<std::vector<int>, std::vector<std::vector<std::size_t> > >>
class Gudhi::coxeter_triangulation::Coxeter_triangulation< Permutahedral_representation_ >

A class that stores Coxeter triangulation of type \(\tilde{A}_d\). This triangulation has the greatest simplex quality out of all linear transformations of the Freudenthal-Kuhn triangulation.

Template Parameters
Permutahedral_representation_Type of a simplex given by a permutahedral representation. Needs to be a model of SimplexInCoxeterTriangulation.
Examples
cell_complex_from_basic_circle_manifold.cpp, manifold_tracing_custom_function.cpp, and manifold_tracing_flat_torus_with_boundary.cpp.

Constructor & Destructor Documentation

◆ Coxeter_triangulation()

template<class Permutahedral_representation_ = Permutahedral_representation<std::vector<int>, std::vector<std::vector<std::size_t> > >>
Gudhi::coxeter_triangulation::Coxeter_triangulation< Permutahedral_representation_ >::Coxeter_triangulation ( std::size_t  dimension)
inline

Constructor of Coxeter triangulation of a given dimension.

Parameters
[in]dimensionThe dimension of the triangulation.

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