Computes the convex hull of N 4D points within a ball of radius R, these points belonging to a lattice of chosen dimension D. Displays their projections along the different main axes.
#include <iostream>
#include <vector>
#include <random>
#include <algorithm>
#include <polyscope/polyscope.h>
#include <polyscope/surface_mesh.h>
#include <polyscope/point_cloud.h>
#include <polyscope/curve_network.h>
#include "DGtal/base/Common.h"
#include "DGtal/helpers/StdDefs.h"
#include "DGtal/geometry/tools/GenericLatticeConvexHull.h"
polyscope::Group *
group[4];
template < Dimension dim, typename TComponent, typename TContainer >
static
{
if ( i != k ) pp[ l++ ] = p[ i ];
return pp;
}
template < Dimension dim, typename TComponent, typename TContainer >
static
{
for ( auto i = 0; i < pV.size(); i++ )
return pV;
}
template < typename Point >
static
std::vector< Point >
int nb, double radius, int amplitude, int aff_dim )
{
std::uniform_int_distribution<int> U(-amplitude, amplitude);
std::vector< Point > P;
int m = std::min( aff_dim, (int) V.size() );
for ( auto k = 0; P.size() < nb && k < 100000; k++ )
{
for ( auto i = 0; i < m; i++ )
{
}
if ( (
B-
A).norm() <= radius )
}
std::shuffle( P.begin(), P.end(),
g );
return P;
}
int main(
int argc,
char* argv[] )
{
typedef Space::Point Point4;
std::cout << "Usage: " << argv[ 0 ] << " [R=30] [N=30] [D=2]\n";
std::cout << "Computes the convex hull of N 4D points within a ball of radius R, these points belonging to a lattice of chosen dimension 0<=D<=3. The output is projected along the 4 canonic projections onto 3D space. You cannot choose D=4 since we cannot diplay the result in 3D.\n";
double radius = argc > 1 ? atof( argv[ 1 ] ) : 30.0;
int nb = argc > 2 ? atoi( argv[ 2 ] ) : 30;
int adim = argc > 3 ? atoi( argv[ 3 ] ) : 2;
if ( nb < 0 ) return 1;
if ( adim < 0 || adim > 3 ) return 1;
std::vector< Point4 >
L = { Point4{ 4, 1, -3, 1 },
Point4{ 0, 2, 5, -1 },
Point4{ -1, -3, 5, 0 },
Point4{ 2, 0, 3, -3 } };
std::vector< Point4 > X
L, nb,
radius,
int( round( radius+0.5 ) ),
adim );
QHull hull;
bool ok = hull.compute( X );
std:: cout << ( ok ? "[PASSED]" : "[FAILED]" ) << " hull=" << hull << "\n";
polyscope::init();
std::string proj[ 4 ] = { "(yzw)", "(xzw)", "(xyw)", "(xyz)" };
{
group [k] = polyscope::createGroup( proj[k] );
psPoints [k] = polyscope::registerPointCloud( proj[k]+
" Points",
psVertices[k] = polyscope::registerPointCloud( proj[k]+
" Vertices",
if ( hull.affine_dimension <= 1 )
{
psBoundary0[k] = polyscope::registerPointCloud( proj[k]+
" Convex hull bdy dim=0",
}
else if ( hull.affine_dimension == 2 )
{
psBoundary1[k] = polyscope::registerCurveNetwork( proj[k]+
" Convex hull bdy dim=1",
hull.facets );
}
else if ( hull.affine_dimension == 3 )
{
psBoundary2[k] = polyscope::registerSurfaceMesh( proj[k]+
" Convex hull bdy dim=2",
hull.facets );
}
}
polyscope::show();
return EXIT_SUCCESS;
}
void project(ProjectedPoint &pp, const Point &p, Dimension a)
Aim: Implements basic operations that will be used in Point and Vector classes.
polyscope::PointCloud * psVertices
polyscope::PointCloud * psBoundary0
polyscope::SurfaceMesh * psBoundary2
polyscope::PointCloud * psPoints
polyscope::CurveNetwork * psBoundary1
polyscope::Group * group[4]
DGtal is the top-level namespace which contains all DGtal functions and types.
DGtal::uint32_t Dimension
Aim: Implements the quickhull algorithm by Barber et al. , a famous arbitrary dimensional convex hull...
std::vector< Point > makeRandomLatticePointsFromDirVectors(int nb, const vector< Point > &V)