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

A class that stores any affine transformation of the Freudenthal-Kuhn triangulation. More...

#include <include/gudhi/Freudenthal_triangulation.h>

Public Types

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.
 

Public Member Functions

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

Detailed Description

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

A class that stores any affine transformation of the Freudenthal-Kuhn triangulation.

The data structure is a record that consists of a matrix that represents the linear transformation of the Freudenthal-Kuhn triangulation and a vector that represents the offset.

Template Parameters
Permutahedral_representation_Type of a simplex given by a permutahedral representation. Needs to be a model of SimplexInCoxeterTriangulation.

Constructor & Destructor Documentation

◆ Freudenthal_triangulation() [1/3]

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

Constructor of the Freudenthal-Kuhn triangulation of a given dimension.

Parameters
[in]dimensionThe dimension of the triangulation.

◆ Freudenthal_triangulation() [2/3]

template<class Permutahedral_representation_ = Permutahedral_representation<std::vector<int>, std::vector<std::vector<std::size_t> > >>
Gudhi::coxeter_triangulation::Freudenthal_triangulation< Permutahedral_representation_ >::Freudenthal_triangulation ( std::size_t  dimension,
const Matrix &  matrix 
)
inline

Constructor of the Freudenthal-Kuhn triangulation of a given dimension under a linear transformation by a given matrix.

Parameters
[in]dimensionThe dimension of the triangulation.
[in]matrixThe matrix that defines the linear transformation. Needs to be invertible.

◆ Freudenthal_triangulation() [3/3]

template<class Permutahedral_representation_ = Permutahedral_representation<std::vector<int>, std::vector<std::vector<std::size_t> > >>
Gudhi::coxeter_triangulation::Freudenthal_triangulation< Permutahedral_representation_ >::Freudenthal_triangulation ( unsigned  dimension,
const Matrix &  matrix,
const Vector &  offset 
)
inline

Constructor of the Freudenthal-Kuhn triangulation of a given dimension under an affine transformation by a given matrix and a translation vector.

Parameters
[in]dimensionThe dimension of the triangulation.
[in]matrixThe matrix that defines the linear transformation. Needs to be invertible.
[in]offsetThe offset vector.
Exceptions
std::invalid_argumentIn debug mode, if offset size is different from dimension.

Member Function Documentation

◆ barycenter()

template<class Permutahedral_representation_ = Permutahedral_representation<std::vector<int>, std::vector<std::vector<std::size_t> > >>
Eigen::VectorXd Gudhi::coxeter_triangulation::Freudenthal_triangulation< Permutahedral_representation_ >::barycenter ( const Simplex_handle simplex,
double  scale = 1 
) const
inline

Returns the Cartesian coordinates of the barycenter of a given simplex.

Using the additional parameter scale, the search can be done in a triangulation that shares the origin, but is scaled by a given factor. This parameter can be useful to simulate the computation of Cartesian coordinates of the barycenter of a simplex in a subdivided triangulation.

Parameters
[in]simplexThe query simplex.
[in]scaleThe scale of the triangulation.

◆ cartesian_coordinates()

template<class Permutahedral_representation_ = Permutahedral_representation<std::vector<int>, std::vector<std::vector<std::size_t> > >>
Eigen::VectorXd Gudhi::coxeter_triangulation::Freudenthal_triangulation< Permutahedral_representation_ >::cartesian_coordinates ( const Vertex_handle vertex,
double  scale = 1 
) const
inline

Returns the Cartesian coordinates of the given vertex.

Using the additional parameter scale, the search can be done in a triangulation that shares the origin, but is scaled by a given factor. This parameter can be useful to simulate the computation of Cartesian coordinates of a vertex in a subdivided triangulation.

Parameters
[in]vertexThe query vertex.
[in]scaleThe scale of the triangulation.

◆ change_matrix()

template<class Permutahedral_representation_ = Permutahedral_representation<std::vector<int>, std::vector<std::vector<std::size_t> > >>
void Gudhi::coxeter_triangulation::Freudenthal_triangulation< Permutahedral_representation_ >::change_matrix ( const Eigen::MatrixXd &  matrix)
inline

Change the linear transformation matrix to a given value.

Parameters
[in]matrixNew value of the linear transformation matrix.
Examples
manifold_tracing_flat_torus_with_boundary.cpp.

◆ change_offset()

template<class Permutahedral_representation_ = Permutahedral_representation<std::vector<int>, std::vector<std::vector<std::size_t> > >>
void Gudhi::coxeter_triangulation::Freudenthal_triangulation< Permutahedral_representation_ >::change_offset ( const Eigen::VectorXd &  offset)
inline

Change the offset vector to a given value.

Parameters
[in]offsetNew value of the offset vector.
Examples
cell_complex_from_basic_circle_manifold.cpp, and manifold_tracing_flat_torus_with_boundary.cpp.

◆ locate_point()

template<class Permutahedral_representation_ = Permutahedral_representation<std::vector<int>, std::vector<std::vector<std::size_t> > >>
template<class Point_d >
Simplex_handle Gudhi::coxeter_triangulation::Freudenthal_triangulation< Permutahedral_representation_ >::locate_point ( const Point_d &  point,
double  scale = 1 
) const
inline

Returns the permutahedral representation of the simplex in the triangulation that contains a given query point.

Using the additional parameter scale, the search can be done in a triangulation that shares the origin, but is scaled by a given factor. This parameter can be useful to simulate the point location in a subdivided triangulation. The returned simplex is always minimal by inclusion.

Template Parameters
Point_dA class that represents a point in d-dimensional Euclidean space. The coordinates should be random-accessible. Needs to provide the method size().
Parameters
[in]pointThe query point.
[in]scaleThe scale of the triangulation.
Exceptions
std::invalid_argumentIn debug mode, if point dimension is different from triangulation one.

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