A class that stores any affine transformation of the Freudenthal-Kuhn triangulation.
More...
|
| 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...
|
|
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
-
template<class Permutahedral_representation_ = Permutahedral_representation<std::vector<int>, std::vector<std::vector<std::size_t> > >>
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] | simplex | The query simplex. |
[in] | scale | The scale of the triangulation. |
template<class Permutahedral_representation_ = Permutahedral_representation<std::vector<int>, std::vector<std::vector<std::size_t> > >>
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] | vertex | The query vertex. |
[in] | scale | The scale of the triangulation. |
template<class Permutahedral_representation_ = Permutahedral_representation<std::vector<int>, std::vector<std::vector<std::size_t> > >>
template<class Point_d >
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_d | A class that represents a point in d-dimensional Euclidean space. The coordinates should be random-accessible. Needs to provide the method size(). |
- Parameters
-
[in] | point | The query point. |
[in] | scale | The scale of the triangulation. |
- Exceptions
-
std::invalid_argument | In debug mode, if point dimension is different from triangulation one. |