Loading [MathJax]/extensions/MathZoom.js
DGtal 2.0.0
DGtal::ConvexCellComplex< TPoint > Struct Template Reference

Aim: represents a d-dimensional complex in a d-dimensional space with the following properties and restrictions: More...

#include <DGtal/geometry/volumes/ConvexCellComplex.h>

Public Types

typedef TPoint Point
typedef std::size_t Index
typedef std::size_t Size
typedef Index Cell
typedef std::pair< Index, bool > Face
typedef Index Vertex
typedef std::vector< IndexIndexRange
typedef std::vector< VertexVertexRange
typedef std::vector< FaceFaceRange
typedef Point::Coordinate Scalar
typedef PointVector< dimension, ScalarVector
typedef PointVector< dimension, double > RealPoint
typedef PointVector< dimension, double > RealVector

Public Member Functions

Standard services
 ConvexCellComplex ()
 Defaut constructor.
void clear ()
 Clears the complex (as if it was just default constructed).
Size nbCells () const
Size nbFaces () const
Size nbVertices () const
Primal structure topology services
Cell infiniteCell () const
bool isInfinite (const Cell c) const
Face opposite (const Face f) const
const FaceRangecellFaces (const Cell c) const
const VertexRangecellVertices (const Cell c) const
const std::vector< VertexRange > & allCellVertices () const
VertexRange faceVertices (const Face f) const
Cell faceCell (const Face f) const
VertexRange faceComplementVertices (const Face f) const
Primal structure geometry services
std::vector< PointcellVertexPositions (const Cell c) const
std::vector< PointfaceVertexPositions (const Face f) const
Point position (const Vertex v) const
RealPoint toReal (const Point p) const
template<typename LatticePoint>
LatticePoint toLattice (const RealPoint p, double factor=1.0) const
RealPoint cellBarycenter (const Cell c) const
template<typename LatticePolytope>
LatticePolytope cellLatticePolytope (const Cell c, const double factor=1.0) const
const VectorfaceNormal (const Face f) const
Scalar faceIntercept (const Face f) const
bool hasFaceGeometry () const
void requireFaceGeometry ()
 Forces the computation of face geometry.
void unrequireFaceGeometry ()
 Release ressources allocated to face geometry.
Interface
void selfDisplay (std::ostream &out) const
bool isValid () const

Data Fields

Data for primal structure
std::vector< FaceRangecell_faces
 Tells if the face geometry has been computed.
std::vector< VertexRangecell_vertices
 Tells if the face geometry has been computed.
std::vector< Celltrue_face_cell
 Tells if the face geometry has been computed.
std::vector< Cellfalse_face_cell
 Tells if the face geometry has been computed.
std::vector< VertexRangetrue_face_vertices
 Tells if the face geometry has been computed.
std::vector< Pointvertex_position
 Tells if the face geometry has been computed.
bool has_face_geometry
 Tells if the face geometry has been computed.
std::vector< Vectortrue_face_normal
 Contains the outward oriented normal of each 'true' face.
std::vector< Scalartrue_face_intercept
 Contains the intercept of each 'true' face.

Static Public Attributes

static const Dimension dimension = TPoint::dimension
static const Index INFINITE_CELL = (Index) -1

Protected Member Functions

Protected services
void computeCellVertices (const Cell c) const
void computeFaceGeometry ()
 Computes for each face its outward oriented normal vector.
std::pair< Vector, ScalarcomputeHalfSpace (const VertexRange &v) const
void reorderFaceVertices (Index f)

Detailed Description

template<typename TPoint>
struct DGtal::ConvexCellComplex< TPoint >

Aim: represents a d-dimensional complex in a d-dimensional space with the following properties and restrictions:

Description of template class 'ConvexCellComplex'

- all cells (maximal and below) are convex

  • k-cells span a k-dimensional affine space
  • cells are convex, their geometry is described by their vertices
  • the boundary of each d-cell is composed of d-1-cells, called faces
  • each face is shared by two d-cells (one can be infinite), but are viewed with opposite orientation from each d-cell
  • each vertex hold its position
  • cells may be spherical (their vertices are cospherical) if the complex comes from a Delaunay triangulation.

All cells are indexed per dimension, starting from 0, hence:

  • vertices are numbered from 0 to nbVertices()-1
  • faces are numbered from 0 to nbFaces()-1
  • cells are numbered from 0 to nbCells()-1 (the infinite cell is not in this range).

The ConvexCellComplex is used to represent the Delaunay decomposition into convex cells of a range of points. It is built with ConvexityHelper functions like ConvexityHelper::computeDelaunayCellComplex.

Template Parameters
TPointan arbitrary model of Point.
See also
moduleQuickHull
Examples
geometry/tools/exampleLatticeBallDelaunay3D.cpp, geometry/tools/exampleQuickHull3D.cpp, and geometry/tools/exampleRationalBallDelaunay3D.cpp.

Definition at line 85 of file ConvexCellComplex.h.

Member Typedef Documentation

◆ Cell

template<typename TPoint>
typedef Index DGtal::ConvexCellComplex< TPoint >::Cell

Definition at line 93 of file ConvexCellComplex.h.

◆ Face

template<typename TPoint>
typedef std::pair< Index, bool > DGtal::ConvexCellComplex< TPoint >::Face

Definition at line 94 of file ConvexCellComplex.h.

◆ FaceRange

template<typename TPoint>
typedef std::vector< Face > DGtal::ConvexCellComplex< TPoint >::FaceRange

Definition at line 98 of file ConvexCellComplex.h.

◆ Index

template<typename TPoint>
typedef std::size_t DGtal::ConvexCellComplex< TPoint >::Index

◆ IndexRange

template<typename TPoint>
typedef std::vector< Index > DGtal::ConvexCellComplex< TPoint >::IndexRange

Definition at line 96 of file ConvexCellComplex.h.

◆ Point

template<typename TPoint>
typedef TPoint DGtal::ConvexCellComplex< TPoint >::Point

Definition at line 88 of file ConvexCellComplex.h.

◆ RealPoint

template<typename TPoint>
typedef PointVector< dimension, double > DGtal::ConvexCellComplex< TPoint >::RealPoint

Definition at line 102 of file ConvexCellComplex.h.

◆ RealVector

template<typename TPoint>
typedef PointVector< dimension, double > DGtal::ConvexCellComplex< TPoint >::RealVector

Definition at line 103 of file ConvexCellComplex.h.

◆ Scalar

template<typename TPoint>
typedef Point::Coordinate DGtal::ConvexCellComplex< TPoint >::Scalar

Definition at line 100 of file ConvexCellComplex.h.

◆ Size

template<typename TPoint>
typedef std::size_t DGtal::ConvexCellComplex< TPoint >::Size

Definition at line 90 of file ConvexCellComplex.h.

◆ Vector

template<typename TPoint>
typedef PointVector< dimension, Scalar > DGtal::ConvexCellComplex< TPoint >::Vector

Definition at line 101 of file ConvexCellComplex.h.

◆ Vertex

template<typename TPoint>
typedef Index DGtal::ConvexCellComplex< TPoint >::Vertex

Definition at line 95 of file ConvexCellComplex.h.

◆ VertexRange

template<typename TPoint>
typedef std::vector< Vertex > DGtal::ConvexCellComplex< TPoint >::VertexRange

Definition at line 97 of file ConvexCellComplex.h.

Constructor & Destructor Documentation

◆ ConvexCellComplex()

template<typename TPoint>
DGtal::ConvexCellComplex< TPoint >::ConvexCellComplex ( )
inline

Defaut constructor.

Definition at line 111 of file ConvexCellComplex.h.

112 { clear(); }
void clear()
Clears the complex (as if it was just default constructed).

References clear().

Member Function Documentation

◆ allCellVertices()

template<typename TPoint>
const std::vector< VertexRange > & DGtal::ConvexCellComplex< TPoint >::allCellVertices ( ) const
inline
Returns
the range of vertex ranges giving for each cell its vertices.

Definition at line 181 of file ConvexCellComplex.h.

182 {
183 if ( cell_vertices.empty() )
184 cell_vertices.resize( nbCells() );
185 for ( Index c = 0; c < nbCells(); ++c )
186 if ( cell_vertices[ c ].empty() )
188 return cell_vertices;
189 }
Aim: represents a d-dimensional complex in a d-dimensional space with the following properties and re...
void computeCellVertices(const Cell c) const
std::vector< VertexRange > cell_vertices
Tells if the face geometry has been computed.

References cell_vertices, computeCellVertices(), and nbCells().

◆ cellBarycenter()

template<typename TPoint>
RealPoint DGtal::ConvexCellComplex< TPoint >::cellBarycenter ( const Cell c) const
inline
Parameters
[in]cany valid cell
Returns
the cell barycenter.
Examples
geometry/tools/exampleLatticeBallDelaunay3D.cpp, and geometry/tools/exampleRationalBallDelaunay3D.cpp.

Definition at line 286 of file ConvexCellComplex.h.

287 {
288 ASSERT( ! isInfinite( c ) );
289 RealPoint b;
290 const auto& vtcs = cellVertices( c );
291 for ( auto v : vtcs )
292 b += toReal( position( v ) );
293 return b / vtcs.size();
294 }
const VertexRange & cellVertices(const Cell c) const
Point position(const Vertex v) const
PointVector< dimension, double > RealPoint
RealPoint toReal(const Point p) const
bool isInfinite(const Cell c) const

References cellVertices(), isInfinite(), position(), DGtal::PointVector< dim, TEuclideanRing, TContainer >::size(), and toReal().

Referenced by main().

◆ cellFaces()

template<typename TPoint>
const FaceRange & DGtal::ConvexCellComplex< TPoint >::cellFaces ( const Cell c) const
inline
Parameters
[in]cany valid cell
Returns
its range of faces.
Examples
geometry/tools/exampleLatticeBallDelaunay3D.cpp, and geometry/tools/exampleRationalBallDelaunay3D.cpp.

Definition at line 161 of file ConvexCellComplex.h.

162 { return cell_faces[ c ]; }
std::vector< FaceRange > cell_faces
Tells if the face geometry has been computed.

References cell_faces.

Referenced by cellLatticePolytope(), computeFaceGeometry(), and main().

◆ cellLatticePolytope()

template<typename TPoint>
template<typename LatticePolytope>
LatticePolytope DGtal::ConvexCellComplex< TPoint >::cellLatticePolytope ( const Cell c,
const double factor = 1.0 ) const
inline
Precondition
hasFaceGeometry() == true
Parameters
[in]cany valid cell
[in]factorthe dilation applied to all points before integer conversion
Returns
the corresponding polytope

Definition at line 302 of file ConvexCellComplex.h.

303 {
305 ASSERT( ! isInfinite( c ) );
306 typedef typename LatticePolytope::Point LatticePoint;
307 typedef typename LatticePolytope::HalfSpace HalfSpace;
308 typedef typename LatticePolytope::Domain Domain;
309 const auto vtcs = cellVertices( c );
311 for ( auto v : vtcs )
312 pts.push_back( toLattice< LatticePoint >( position ( v ), factor ) );
313 LatticePoint l = pts[ 0 ];
314 LatticePoint u = pts[ 0 ];
315 for ( const auto& p : pts ) {
316 l = l.inf( p );
317 u = u.sup( p );
318 }
319 Domain domain( l, u );
321 for ( auto f : cellFaces( c ) ) {
322 HS.emplace_back( faceNormal( f ), faceIntercept( f ) );
323 }
324 return LatticePolytope( domain, HS.cbegin(), HS.cend(), false );
325 }
const Vector & faceNormal(const Face f) const
Scalar faceIntercept(const Face f) const
LatticePoint toLattice(const RealPoint p, double factor=1.0) const
const FaceRange & cellFaces(const Cell c) const
Domain domain

References cellFaces(), cellVertices(), domain, faceIntercept(), faceNormal(), hasFaceGeometry(), isInfinite(), position(), and toLattice().

◆ cellVertexPositions()

template<typename TPoint>
std::vector< Point > DGtal::ConvexCellComplex< TPoint >::cellVertexPositions ( const Cell c) const
inline
Parameters
[in]cany valid cell (non infinite)
Returns
the range of positions of all cell vertices

Definition at line 233 of file ConvexCellComplex.h.

234 {
235 ASSERT( ! isInfinite( c ) );
236 const auto vtcs = cellVertices( c );
237 std::vector< Point > pts( vtcs.size() );
238 std::transform( vtcs.cbegin(), vtcs.cend(), pts.begin(),
239 [&] ( Vertex v ) { return this->position( v ); } );
240 return pts;
241 }

References cellVertices(), and isInfinite().

◆ cellVertices()

template<typename TPoint>
const VertexRange & DGtal::ConvexCellComplex< TPoint >::cellVertices ( const Cell c) const
inline

Lazy computation of cell vertices from face vertices if they were not given.

Parameters
[in]cany (finite) cell
Returns
its vertices
Examples
geometry/tools/exampleLatticeBallDelaunay3D.cpp, and geometry/tools/exampleRationalBallDelaunay3D.cpp.

Definition at line 169 of file ConvexCellComplex.h.

170 {
171 ASSERT( c != INFINITE_CELL && c < nbCells() );
172 ASSERT( ! cell_vertices.empty() );
173 if ( cell_vertices.empty() )
174 cell_vertices.resize( nbCells() );
175 if ( cell_vertices[ c ].empty() )
177 return cell_vertices[ c ];
178 }
static const Index INFINITE_CELL

References cell_vertices, computeCellVertices(), INFINITE_CELL, and nbCells().

Referenced by cellBarycenter(), cellLatticePolytope(), cellVertexPositions(), faceComplementVertices(), and main().

◆ clear()

template<typename TPoint>
void DGtal::ConvexCellComplex< TPoint >::clear ( )
inline

Clears the complex (as if it was just default constructed).

Definition at line 115 of file ConvexCellComplex.h.

116 {
117 cell_faces.clear();
118 cell_vertices.clear();
119 true_face_cell.clear();
120 false_face_cell.clear();
121 true_face_vertices.clear();
122 vertex_position.clear();
123 has_face_geometry = false;
124 true_face_normal.clear();
125 true_face_intercept.clear();
126 }
std::vector< VertexRange > true_face_vertices
Tells if the face geometry has been computed.
std::vector< Scalar > true_face_intercept
Contains the intercept of each 'true' face.
std::vector< Cell > false_face_cell
Tells if the face geometry has been computed.
std::vector< Cell > true_face_cell
Tells if the face geometry has been computed.
bool has_face_geometry
Tells if the face geometry has been computed.
std::vector< Vector > true_face_normal
Contains the outward oriented normal of each 'true' face.
std::vector< Point > vertex_position
Tells if the face geometry has been computed.

References cell_faces, cell_vertices, false_face_cell, has_face_geometry, true_face_cell, true_face_intercept, true_face_normal, true_face_vertices, and vertex_position.

Referenced by ConvexCellComplex().

◆ computeCellVertices()

template<typename TPoint>
void DGtal::ConvexCellComplex< TPoint >::computeCellVertices ( const Cell c) const
inlineprotected

computes the vertices of a given cell c upon request.

Parameters
cany valid cell

Definition at line 432 of file ConvexCellComplex.h.

433 {
435 for ( auto f : cell_faces[ c ] )
436 vertices.insert( true_face_vertices[ f.first ].cbegin(),
437 true_face_vertices[ f.first ].cend() );
438 cell_vertices[ c ] = VertexRange( vertices.cbegin(), vertices.cend() );
439 }
std::vector< Vertex > VertexRange

References cell_faces, cell_vertices, and true_face_vertices.

Referenced by allCellVertices(), and cellVertices().

◆ computeFaceGeometry()

template<typename TPoint>
void DGtal::ConvexCellComplex< TPoint >::computeFaceGeometry ( )
inlineprotected

Computes for each face its outward oriented normal vector.

Definition at line 442 of file ConvexCellComplex.h.

443 {
444 true_face_normal .resize( nbFaces() );
445 true_face_intercept.resize( nbFaces() );
446 // Compute normals and intercepts
447 for ( Index f = 0; f < nbFaces(); ++f ) {
450 if ( true_face_normal[ f ] == Vector::zero )
451 trace.error() << "[ConvexCellComplex::computeFaceGeometry]"
452 << " null normal vector at face " << f << std::endl;
453 }
454 // Orient them consistently
455 for ( Index c = 0; c < nbCells(); c++ ) {
456 //const auto& c_vtcs = cellVertices( c ); //not used
457 for ( auto f : cellFaces( c ) ) {
458 if ( ! f.second ) continue; // process only 'true' faces
459 const auto ov = faceComplementVertices( f );
460 ASSERT( ! ov.empty() );
461 const Vector N = true_face_normal [ f.first ];
462 const Scalar nu = true_face_intercept[ f.first ];
463 const Scalar iv = N.dot( position( ov.front() ) );
464 if ( ( iv - nu ) > 0 ) {
465 // Should be inside
466 true_face_normal [ f.first ] = -N;
467 true_face_intercept[ f.first ] = -nu;
468 std::reverse( true_face_vertices[ f.first ].begin(),
469 true_face_vertices[ f.first ].end() );
470 }
471 reorderFaceVertices( f.first );
472 }
473 }
474 has_face_geometry = true;
475 }
VertexRange faceComplementVertices(const Face f) const
PointVector< dimension, Scalar > Vector
std::pair< Vector, Scalar > computeHalfSpace(const VertexRange &v) const

References cellFaces(), computeHalfSpace(), DGtal::PointVector< dim, TEuclideanRing, TContainer >::dot(), faceComplementVertices(), has_face_geometry, nbCells(), nbFaces(), position(), reorderFaceVertices(), DGtal::trace, true_face_intercept, true_face_normal, true_face_vertices, and DGtal::PointVector< dimension, Scalar >::zero.

Referenced by requireFaceGeometry().

◆ computeHalfSpace()

template<typename TPoint>
std::pair< Vector, Scalar > DGtal::ConvexCellComplex< TPoint >::computeHalfSpace ( const VertexRange & v) const
inlineprotected
Parameters
va range of vertices of size() >= dimension
Returns
the pair normal/intercept of the face.

Definition at line 480 of file ConvexCellComplex.h.

481 {
483 Matrix A;
484 Vector N;
485 Scalar c;
486 for ( int i = 1; i < dimension; i++ )
487 for ( int j = 0; j < dimension; j++ )
488 A.setComponent( i-1, j, position( v[ i ] )[ j ] - position( v[ 0 ] )[ j ] );
489 for ( int j = 0; j < dimension; j++ )
490 N[ j ] = A.cofactor( dimension-1, j );
491 c = N.dot( position( v[ 0 ] ) );
492 return std::make_pair( N, c );
493 }
static const Dimension dimension

References dimension, DGtal::PointVector< dim, TEuclideanRing, TContainer >::dot(), and position().

Referenced by computeFaceGeometry().

◆ faceCell()

template<typename TPoint>
Cell DGtal::ConvexCellComplex< TPoint >::faceCell ( const Face f) const
inline
Parameters
[in]fany face
Returns
its associated cell (or INFINITE_CELL)

Definition at line 202 of file ConvexCellComplex.h.

203 {
204 return f.second ? true_face_cell[ f.first ] : false_face_cell[ f.first ];
205 }

References false_face_cell, and true_face_cell.

Referenced by faceComplementVertices().

◆ faceComplementVertices()

template<typename TPoint>
VertexRange DGtal::ConvexCellComplex< TPoint >::faceComplementVertices ( const Face f) const
inline
Parameters
[in]fany face
Returns
the vertices within the cell containing f that are not in f.

Definition at line 211 of file ConvexCellComplex.h.

212 {
214 const Cell c = faceCell( f );
215 if ( isInfinite( c ) ) return result;
216 const auto& c_vtcs = cellVertices( c );
217 auto f_vtcs = true_face_vertices[ f.first ];
218 std::sort( f_vtcs.begin(), f_vtcs.end() );
219 for ( auto v : c_vtcs )
220 if ( ! std::binary_search( f_vtcs.cbegin(), f_vtcs.cend(), v ) )
221 result.push_back( v );
222 return result;
223 }
Cell faceCell(const Face f) const

References cellVertices(), faceCell(), isInfinite(), and true_face_vertices.

Referenced by computeFaceGeometry().

◆ faceIntercept()

template<typename TPoint>
Scalar DGtal::ConvexCellComplex< TPoint >::faceIntercept ( const Face f) const
inline
Precondition
hasFaceGeometry() == true
Parameters
[in]fany face
Returns
its intercept.

Definition at line 341 of file ConvexCellComplex.h.

342 {
344 ASSERT( f.first < nbFaces() );
345 if ( f.second ) return true_face_intercept[ f.first ];
346 else return -true_face_intercept[ f.first ];
347 }

References hasFaceGeometry(), nbFaces(), and true_face_intercept.

Referenced by cellLatticePolytope().

◆ faceNormal()

template<typename TPoint>
const Vector & DGtal::ConvexCellComplex< TPoint >::faceNormal ( const Face f) const
inline
Precondition
hasFaceGeometry() == true
Parameters
[in]fany face
Returns
its normal oriented outward its enclosing cell.

Definition at line 330 of file ConvexCellComplex.h.

331 {
333 ASSERT( f.first < nbFaces() );
334 if ( f.second ) return true_face_normal[ f.first ];
335 else return -true_face_normal[ f.first ];
336 }

References hasFaceGeometry(), nbFaces(), and true_face_normal.

Referenced by cellLatticePolytope().

◆ faceVertexPositions()

template<typename TPoint>
std::vector< Point > DGtal::ConvexCellComplex< TPoint >::faceVertexPositions ( const Face f) const
inline
Parameters
[in]fany face
Returns
the range of positions of all face vertices

Definition at line 245 of file ConvexCellComplex.h.

246 {
247 const auto vtcs = faceVertices( f );
248 std::vector< Point > pts( vtcs.size() );
249 std::transform( vtcs.cbegin(), vtcs.cend(), pts.begin(),
250 [&] ( Vertex v ) { return this->position( v ); } );
251 return pts;
252 }
VertexRange faceVertices(const Face f) const

References faceVertices().

◆ faceVertices()

template<typename TPoint>
VertexRange DGtal::ConvexCellComplex< TPoint >::faceVertices ( const Face f) const
inline
Parameters
[in]fany face
Returns
its vertices (in order depending on the side of the face)
Examples
geometry/tools/exampleLatticeBallDelaunay3D.cpp, and geometry/tools/exampleRationalBallDelaunay3D.cpp.

Definition at line 193 of file ConvexCellComplex.h.

194 {
195 if ( f.second ) return true_face_vertices[ f.first ];
196 return VertexRange( true_face_vertices[ f.first ].crbegin(),
197 true_face_vertices[ f.first ].crend() );
198 }

References true_face_vertices.

Referenced by faceVertexPositions(), and main().

◆ hasFaceGeometry()

template<typename TPoint>
bool DGtal::ConvexCellComplex< TPoint >::hasFaceGeometry ( ) const
inline
Returns
'true' iff the face geometry is computed.

Definition at line 350 of file ConvexCellComplex.h.

351 { return has_face_geometry; }

References has_face_geometry.

Referenced by cellLatticePolytope(), faceIntercept(), faceNormal(), requireFaceGeometry(), and selfDisplay().

◆ infiniteCell()

template<typename TPoint>
Cell DGtal::ConvexCellComplex< TPoint >::infiniteCell ( ) const
inline
Returns
the index of the infinite cell.

Definition at line 145 of file ConvexCellComplex.h.

146 { return INFINITE_CELL; }

References INFINITE_CELL.

◆ isInfinite()

template<typename TPoint>
bool DGtal::ConvexCellComplex< TPoint >::isInfinite ( const Cell c) const
inline
Parameters
[in]cany cell
Returns
'true' is the cell is infinite

Definition at line 150 of file ConvexCellComplex.h.

151 { return c == INFINITE_CELL; }

References INFINITE_CELL.

Referenced by cellBarycenter(), cellLatticePolytope(), cellVertexPositions(), and faceComplementVertices().

◆ isValid()

template<typename TPoint>
bool DGtal::ConvexCellComplex< TPoint >::isValid ( ) const
inline

Checks the validity/consistency of the object.

Returns
'true' if the object is valid (at least one cell), 'false' otherwise.

Definition at line 419 of file ConvexCellComplex.h.

420 {
421 return nbCells() != 0;
422 }

References nbCells().

◆ nbCells()

template<typename TPoint>
Size DGtal::ConvexCellComplex< TPoint >::nbCells ( ) const
inline
Returns
the number of finite cells with maximal dimension
Examples
geometry/tools/exampleLatticeBallDelaunay3D.cpp, and geometry/tools/exampleRationalBallDelaunay3D.cpp.

Definition at line 129 of file ConvexCellComplex.h.

130 { return cell_faces.size(); }

References cell_faces.

Referenced by allCellVertices(), cellVertices(), computeFaceGeometry(), isValid(), main(), and selfDisplay().

◆ nbFaces()

template<typename TPoint>
Size DGtal::ConvexCellComplex< TPoint >::nbFaces ( ) const
inline
Returns
the number of unoriented cells of codimension 1

Definition at line 132 of file ConvexCellComplex.h.

133 { return true_face_cell.size(); }

References true_face_cell.

Referenced by computeFaceGeometry(), faceIntercept(), faceNormal(), and selfDisplay().

◆ nbVertices()

template<typename TPoint>
Size DGtal::ConvexCellComplex< TPoint >::nbVertices ( ) const
inline
Returns
the number of vertices

Definition at line 135 of file ConvexCellComplex.h.

136 { return vertex_position.size(); }

References vertex_position.

Referenced by selfDisplay().

◆ opposite()

template<typename TPoint>
Face DGtal::ConvexCellComplex< TPoint >::opposite ( const Face f) const
inline
Parameters
[in]fany face
Returns
its opposite face

Definition at line 156 of file ConvexCellComplex.h.

157 { return std::make_pair( f.first, ! f.second ); }

◆ position()

template<typename TPoint>
Point DGtal::ConvexCellComplex< TPoint >::position ( const Vertex v) const
inline
Parameters
[in]vany vertex
Returns
the vertex position
Examples
geometry/tools/exampleLatticeBallDelaunay3D.cpp, and geometry/tools/exampleRationalBallDelaunay3D.cpp.

Definition at line 256 of file ConvexCellComplex.h.

257 {
258 return vertex_position[ v ];
259 }

References vertex_position.

Referenced by cellBarycenter(), cellLatticePolytope(), computeFaceGeometry(), computeHalfSpace(), and main().

◆ reorderFaceVertices()

template<typename TPoint>
void DGtal::ConvexCellComplex< TPoint >::reorderFaceVertices ( Index f)
inlineprotected

Given a face index f, reorder all its vertices so that they are oriented consistently with respect to the normal.

Parameters
fany valid face index
Note
the order of points is consistent if, picking any d-simplex in order in this range, their associated half-spaces have all the same orientation.

Definition at line 503 of file ConvexCellComplex.h.

504 {
507 // Sort a facet such that its points, taken in order, have
508 // always the same orientation of the facet. More precisely,
509 // facets span a `dimension-1` vector space. There are thus
510 // dimension-2 fixed points, and the last ones (at least two)
511 // may be reordered.
513 for ( Dimension k = 0; k < dimension-2; k++ )
514 splx[ k ] = result[ k ];
515 std::sort( result.begin()+dimension-2, result.end(),
516 [&] ( Index i, Index j )
517 {
518 splx[ dimension-2 ] = i;
519 splx[ dimension-1 ] = j;
520 const auto H = computeHalfSpace( splx );
521 const auto orient = N.dot( H.first );
522 return orient > 0;
523 } );
524 }

References dimension, true_face_normal, and true_face_vertices.

Referenced by computeFaceGeometry().

◆ requireFaceGeometry()

template<typename TPoint>
void DGtal::ConvexCellComplex< TPoint >::requireFaceGeometry ( )
inline

Forces the computation of face geometry.

Examples
geometry/tools/exampleLatticeBallDelaunay3D.cpp, and geometry/tools/exampleRationalBallDelaunay3D.cpp.

Definition at line 354 of file ConvexCellComplex.h.

355 {
356 if ( ! hasFaceGeometry() )
358 }
void computeFaceGeometry()
Computes for each face its outward oriented normal vector.

References computeFaceGeometry(), and hasFaceGeometry().

Referenced by main().

◆ selfDisplay()

template<typename TPoint>
void DGtal::ConvexCellComplex< TPoint >::selfDisplay ( std::ostream & out) const
inline

Writes/Displays the object on an output stream.

Parameters
outthe output stream where the object is written.

Definition at line 405 of file ConvexCellComplex.h.

406 {
407 out << "[ConvexCellComplex<" << dimension << ">"
408 << " #C=" << nbCells()
409 << " #F=" << nbFaces()
410 << " #V=" << nbVertices();
411 if ( hasFaceGeometry() ) out << " hasFaceGeometry";
412 out << " ]";
413 }

References dimension, hasFaceGeometry(), nbCells(), nbFaces(), and nbVertices().

◆ toLattice()

template<typename TPoint>
template<typename LatticePoint>
LatticePoint DGtal::ConvexCellComplex< TPoint >::toLattice ( const RealPoint p,
double factor = 1.0 ) const
inline
Parameters
[in]pany real point
[in]factorthe dilation applied to point p before integer conversion
Returns
its embedding in the real Euclidean space.

Definition at line 275 of file ConvexCellComplex.h.

276 {
278 for ( Dimension k = 0; k < dimension; ++k )
279 x[ k ] = (typename LatticePoint::Coordinate) round( p[ k ] * factor );
280 return x;
281 }

References dimension.

Referenced by cellLatticePolytope().

◆ toReal()

template<typename TPoint>
RealPoint DGtal::ConvexCellComplex< TPoint >::toReal ( const Point p) const
inline
Parameters
[in]pany point
Returns
its embedding in the real Euclidean space.
Examples
geometry/tools/exampleLatticeBallDelaunay3D.cpp.

Definition at line 263 of file ConvexCellComplex.h.

264 {
265 RealPoint x;
266 for ( Dimension k = 0; k < dimension; ++k )
267 x[ k ] = (double) p[ k ];
268 return x;
269 }

References dimension.

Referenced by cellBarycenter(), and main().

◆ unrequireFaceGeometry()

template<typename TPoint>
void DGtal::ConvexCellComplex< TPoint >::unrequireFaceGeometry ( )
inline

Release ressources allocated to face geometry.

Definition at line 361 of file ConvexCellComplex.h.

362 {
363 has_face_geometry = false;
364 true_face_normal.clear();
365 true_face_intercept.clear();
366 }

References has_face_geometry, true_face_intercept, and true_face_normal.

Field Documentation

◆ cell_faces

template<typename TPoint>
std::vector< FaceRange > DGtal::ConvexCellComplex< TPoint >::cell_faces

Tells if the face geometry has been computed.

Definition at line 376 of file ConvexCellComplex.h.

Referenced by cellFaces(), clear(), computeCellVertices(), and nbCells().

◆ cell_vertices

template<typename TPoint>
std::vector< VertexRange > DGtal::ConvexCellComplex< TPoint >::cell_vertices
mutable

Tells if the face geometry has been computed.

Definition at line 378 of file ConvexCellComplex.h.

Referenced by allCellVertices(), cellVertices(), clear(), and computeCellVertices().

◆ dimension

template<typename TPoint>
const Dimension DGtal::ConvexCellComplex< TPoint >::dimension = TPoint::dimension
static

◆ false_face_cell

template<typename TPoint>
std::vector< Cell > DGtal::ConvexCellComplex< TPoint >::false_face_cell

Tells if the face geometry has been computed.

Definition at line 382 of file ConvexCellComplex.h.

Referenced by clear(), and faceCell().

◆ has_face_geometry

template<typename TPoint>
bool DGtal::ConvexCellComplex< TPoint >::has_face_geometry

Tells if the face geometry has been computed.

Definition at line 389 of file ConvexCellComplex.h.

Referenced by clear(), computeFaceGeometry(), hasFaceGeometry(), and unrequireFaceGeometry().

◆ INFINITE_CELL

template<typename TPoint>
const Index DGtal::ConvexCellComplex< TPoint >::INFINITE_CELL = (Index) -1
static

Definition at line 92 of file ConvexCellComplex.h.

Referenced by cellVertices(), infiniteCell(), and isInfinite().

◆ true_face_cell

template<typename TPoint>
std::vector< Cell > DGtal::ConvexCellComplex< TPoint >::true_face_cell

Tells if the face geometry has been computed.

Definition at line 380 of file ConvexCellComplex.h.

Referenced by clear(), faceCell(), and nbFaces().

◆ true_face_intercept

template<typename TPoint>
std::vector< Scalar > DGtal::ConvexCellComplex< TPoint >::true_face_intercept

Contains the intercept of each 'true' face.

Definition at line 393 of file ConvexCellComplex.h.

Referenced by clear(), computeFaceGeometry(), faceIntercept(), and unrequireFaceGeometry().

◆ true_face_normal

template<typename TPoint>
std::vector< Vector > DGtal::ConvexCellComplex< TPoint >::true_face_normal

Contains the outward oriented normal of each 'true' face.

Definition at line 391 of file ConvexCellComplex.h.

Referenced by clear(), computeFaceGeometry(), faceNormal(), reorderFaceVertices(), and unrequireFaceGeometry().

◆ true_face_vertices

template<typename TPoint>
std::vector< VertexRange > DGtal::ConvexCellComplex< TPoint >::true_face_vertices

Tells if the face geometry has been computed.

Definition at line 384 of file ConvexCellComplex.h.

Referenced by clear(), computeCellVertices(), computeFaceGeometry(), faceComplementVertices(), faceVertices(), and reorderFaceVertices().

◆ vertex_position

template<typename TPoint>
std::vector< Point > DGtal::ConvexCellComplex< TPoint >::vertex_position

Tells if the face geometry has been computed.

Definition at line 386 of file ConvexCellComplex.h.

Referenced by clear(), nbVertices(), and position().


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