31#if defined(GenericLatticeConvexHull_RECURSES)
32#error Recursive header files inclusion detected in GenericLatticeConvexHull.h
35#define GenericLatticeConvexHull_RECURSES
37#if !defined GenericLatticeConvexHull_h
39#define GenericLatticeConvexHull_h
48#include "DGtal/base/Common.h"
49#include "DGtal/kernel/SpaceND.h"
50#include "DGtal/geometry/tools/AffineGeometry.h"
51#include "DGtal/geometry/tools/AffineBasis.h"
52#include "DGtal/geometry/tools/QuickHull.h"
53#include "DGtal/geometry/tools/QuickHullKernels.h"
54#include "DGtal/geometry/volumes/BoundedLatticePolytope.h"
60 typename TCoordinateInteger,
61 typename TInternalInteger >
66 typename TCoordinateInteger,
67 typename TInternalInteger,
72 <
K,TCoordinateInteger,TInternalInteger >
Kernel;
108 template <
typename TInputPo
int>
110 const std::vector< TInputPoint >& X )
112 typedef TInputPoint InputPoint;
138 basis = Basis( X[0], normal, Basis::Type::ECHELON_REDUCED );
142 basis = Basis( X, Basis::Type::SHORTEST_ECHELON_REDUCED );
150 bool ok_hull =
hull.computeConvexHull( QHull::Status::VerticesCompleted );
151 if ( ! ok_hull || ! ok_input )
153 trace.error() <<
"[GenericLatticeConvexHullComputers::compute]"
154 <<
" Error in quick hull computation.\n"
155 <<
"qhull=" <<
hull <<
"\n";
159 points.resize( X.size() );
160 for (
Size i = 0; i < points.size(); i++ )
161 points[ i ] = Affine::transform( X[ i ] );
162 hull.getVertex2Point( v2p );
163 hull.getFacetVertices( facets );
164 positions.resize( v2p.size() );
165 for (
Size i = 0; i < positions.size(); i++ )
166 positions[ i ] = X[ v2p[ i ] ];
168 for (
Size i = 0; i < ppoints.size(); i++ )
176 Basis::Type::SHORTEST_ECHELON_REDUCED,
194 if (
polytope.nbHalfSpaces() > 0 )
return true;
198 for ( std::size_t i = 1; i <
proj_points.size(); i++ )
207 std::vector< ConvexHullHalfSpace > HS;
208 std::vector< PolytopeHalfSpace > PHS;
209 hull.getFacetHalfSpaces( HS );
210 PHS.reserve( HS.size() );
211 for (
auto&
H : HS ) {
216 ::cast(
H.internalNormal()[ i ] );
217 Point Ns = Affine::simplifiedVector( N );
220 PHS.emplace_back( Ns, nu /
g );
240 if (
polytope.nbHalfSpaces() == 0 )
return -1;
260 if (
polytope.nbHalfSpaces() == 0 )
return -1;
280 if (
polytope.nbHalfSpaces() == 0 )
return -1;
306 if (
polytope.nbHalfSpaces() == 0 )
return -1;
320 typename TCoordinateInteger,
321 typename TInternalInteger >
325 < 1,TCoordinateInteger,TInternalInteger >
Kernel;
354 template <
typename TInputPo
int>
356 const std::vector< TInputPoint >& X )
358 typedef TInputPoint InputPoint;
361 typedef std::size_t
Index;
379 points.resize( X.size() );
380 for (
Size i = 0; i < points.size(); i++ )
381 points[ i ] = Affine::transform( X[ i ] );
383 positions.resize( 1 );
384 positions[ 0 ] = X[ 0 ];
406 basis = Basis( X, Basis::Type::SHORTEST_ECHELON_REDUCED );
420 points.resize( X.size() );
421 for (
Size i = 0; i < points.size(); i++ )
422 points[ i ] = Affine::transform( X[ i ] );
424 for (
Size i = 0; i < ppoints.size(); i++ )
432 positions.resize( 2 );
433 positions[ 0 ] = X[ v2p[ 0 ] ];
434 positions[ 1 ] = X[ v2p[ 1 ] ];
437 Basis::Type::SHORTEST_ECHELON_REDUCED,
439 auto dx = Affine::transform( points[ v2p[ 1 ] ] - points[ v2p[ 0 ] ] );
440 auto sdx = Affine::simplifiedVector( dx );
441 Integer n = dx.normInfinity() / sdx.normInfinity();
585 template <
typename InputPo
int >
586 bool compute(
const std::vector< InputPoint >& input_points )
590 std::vector< Size > indices = Affine::affineSubset( input_points );
594 std::cout <<
"Generic Convex hull #V=" <<
positions.size()
595 <<
" #F=" <<
facets.size() <<
"\n";
598 std::cout <<
"F_" << i <<
" = (";
599 for (
auto v :
facets[ i ] ) std::cout <<
" " << v;
603 std::cout <<
"V_" << i
718 out <<
"[GenericLatticeConvexHull"
720 <<
" #in=" <<
points.size()
723 <<
" #F=" <<
facets.size()
789 typename TCoordinateInteger,
790 typename TInternalInteger >
795 object.selfDisplay( out );
808#undef GenericLatticeConvexHull_RECURSES
ClosedIntegerHalfPlane< Space > HalfSpace
HyperRectDomain< Space > Domain
UnsignedComponent norm1() const
auto inf(const PointVector< dim, OtherComponent, OtherStorage > &aPoint) const -> decltype(DGtal::inf(*this, aPoint))
Implements the infimum (or greatest lower bound).
auto sup(const PointVector< dim, OtherComponent, OtherStorage > &aPoint) const -> decltype(DGtal::sup(*this, aPoint))
Implements the supremum (or least upper bound).
static const Dimension dimension
detail namespace gathers internal classes and functions.
static void getOrthogonalVector(TPoint &w, const std::vector< TInputPoint > &X, const TIndexRange &I, const double tolerance=1e-12)
DGtal is the top-level namespace which contains all DGtal functions and types.
std::int64_t int64_t
signed 94-bit integer.
std::ostream & operator<<(std::ostream &out, const ClosedIntegerHalfPlane< TSpace > &object)
DGtal::uint32_t Dimension
Aim: Utility class to determine the affine geometry of an input set of points. It provides exact resu...
Aim: Utility class to determine the affine geometry of an input set of points. It provides exact resu...
InternalInteger InternalScalar
CoordinateInteger CoordinateScalar
DGtal::PointVector< dim, CoordinateInteger > CoordinatePoint
Aim: a geometric kernel to compute the convex hull of digital points with integer-only arithmetic.
Aim: Implements the quickhull algorithm by Barber et al. barber1996, a famous arbitrary dimensional c...
Kernel::InternalScalar InternalInteger
std::vector< Index > IndexRange
bool compute(const std::vector< InputPoint > &input_points)
AffineBasis< OutputPoint > affine_basis
std::vector< OutputPoint > points
GenericLatticeConvexHull(const Kernel &K_=Kernel(), int dbg=0)
detail::GenericLatticeConvexHullComputers< dim, TCoordinateInteger, TInternalInteger, dim > GenericComputers
GenericComputers generic_computers
void selfDisplay(std::ostream &out) const
void clear()
Clears the object as if no computations have been made.
std::vector< OutputPoint > positions
Kernel::CoordinatePoint Point
static const Size dimension
Integer projected_dilation
std::vector< OutputPoint > projected_points
Integer countUpTo(Integer max)
std::vector< IndexRange > facets
Kernel::CoordinateScalar Integer
ConvexHullIntegralKernel< dim, TCoordinateInteger, TInternalInteger > Kernel
----------— INTEGER/POINT CONVERSION SERVICES -----------------—
static Integer cast(Integer i)
Kernel::HalfSpace HalfSpace
std::vector< Point > proj_points
the projected points, as points in lower dimension
static const Dimension dimension
Point::Coordinate Integer
ConvexHullIntegralKernel< 1, TCoordinateInteger, TInternalInteger > Kernel
GenericLatticeConvexHullComputers(Computer *ptrGenQHull)
bool compute(const std::vector< Size > &I, const std::vector< TInputPoint > &X)
Integer countUpTo(Integer max)
void clear()
Clears the object as if no computations have been made.
bool makePolytope()
Does nothing for this 1-dimensional specialization.
Computer::OutputPoint OutputPoint
Integer nb_in_hull
the number of lattice points in the convex hull.
Computer * ptr_gen_qhull
the pointer on the parent computer
GenericLatticeConvexHull< dim, TCoordinateInteger, TInternalInteger > Computer
Kernel::CoordinatePoint Point
SpaceND< 1, Integer > Space
detail::GenericLatticeConvexHullComputers< dim, TCoordinateInteger, TInternalInteger, K-1 > LowerKernels
GenericLatticeConvexHullComputers(Computer *ptrGenQHull)
GenericLatticeConvexHull< dim, TCoordinateInteger, TInternalInteger > Computer
LowerKernels lower_kernels
void clear()
Clears the object as if no computations have been made.
Point::Coordinate Integer
BoundedLatticePolytope< Space > LatticePolytope
static const Dimension dimension
Kernel::CoordinatePoint Point
SpaceND< K, Integer > Space
bool compute(const std::vector< Size > &I, const std::vector< TInputPoint > &X)
ConvexHullIntegralKernel< K, TCoordinateInteger, TInternalInteger > Kernel
Computer::OutputPoint OutputPoint
QuickHull< Kernel > QHull
Integer countUpTo(Integer max)
std::vector< Point > proj_points