DGtal 2.1.0
Loading...
Searching...
No Matches
exampleGenericLatticeConvexHull3D.cpp File Reference
#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"
Include dependency graph for exampleGenericLatticeConvexHull3D.cpp:

Go to the source code of this file.

Functions

std::mt19937 g (rd())
 
template<typename Point >
static std::vector< PointmakeRandomLatticePointsFromDirVectors (Point A, const std::vector< Point > &V, int nb, double radius, int amplitude, int aff_dim)
 
int main (int argc, char *argv[])
 

Variables

polyscope::PointCloud * psPoints
 
polyscope::PointCloud * psVertices
 
polyscope::PointCloud * psBoundary0
 
polyscope::CurveNetwork * psBoundary1
 
polyscope::SurfaceMesh * psBoundary2
 
std::random_device rd
 

Detailed Description

This program is free software: you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

You should have received a copy of the GNU General Public License along with this program. If not, see http://www.gnu.org/licenses/.

Author
Jacques-Olivier Lachaud (jacqu.nosp@m.es-o.nosp@m.livie.nosp@m.r.la.nosp@m.chaud.nosp@m.@uni.nosp@m.v-sav.nosp@m.oie..nosp@m.fr ) Laboratory of Mathematics (CNRS, UMR 5127), University of Savoie, France
Date
2025/09/07

This file is part of the DGtal library.

Definition in file exampleGenericLatticeConvexHull3D.cpp.

Function Documentation

◆ g()

std::mt19937 g ( rd()  )

◆ main()

int main ( int  argc,
char *  argv[] 
)

Definition at line 109 of file exampleGenericLatticeConvexHull3D.cpp.

110{
112 typedef SpaceND< 3, int > Space;
113 typedef Space::Point Point;
114
115 std::cout << "Usage: " << argv[ 0 ] << " [R=30] [N=30] [D=2]\n";
116 std::cout << "Computes the convex hull of N points within a ball of radius R, these points belonging to a lattice of chosen dimension D.\n";
117 double radius = argc > 1 ? atof( argv[ 1 ] ) : 30.0;
118 int nb = argc > 2 ? atoi( argv[ 2 ] ) : 30;
119 int adim = argc > 3 ? atoi( argv[ 3 ] ) : 2;
120 if ( nb < 0 ) return 1;
121 if ( adim < 0 || adim > 3 ) return 1;
122
123 // Create points
124 std::vector< Point > L = { Point{ 4, 1, -3 }, Point{ 0, 2, 5 }, Point{ -1, -3, 5 } };
125 std::vector< Point > X
127 L, nb,
128 radius,
129 int( round( radius+0.5 ) ),
130 adim );
131
132 // Compute convex hull
133 QHull hull;
134 bool ok = hull.compute( X );
135 std:: cout << ( ok ? "[PASSED]" : "[FAILED]" ) << " hull=" << hull << "\n";
136 // Initialize polyscope
137 polyscope::init();
138 psPoints = polyscope::registerPointCloud( "Points", X );
139 psVertices = polyscope::registerPointCloud( "Vertices", hull.positions );
140 if ( hull.affine_dimension <= 1 ) // 1 or 2 points
141 {
142 // no facets
143 psBoundary0 = polyscope::registerPointCloud( "Convex hull bdy dim=0",
144 hull.positions );
145 }
146 else if ( hull.affine_dimension == 2 ) // 2D
147 {
148 // facets are edges (and implicitly converted by polyscope)
149 psBoundary1 = polyscope::registerCurveNetwork( "Convex hull bdy dim=1",
150 hull.positions, hull.facets );
151 psBoundary0 = polyscope::registerPointCloud( "Projected points",
152 hull.projected_points );
153 std::set<Point> S( hull.projected_points.cbegin(),
154 hull.projected_points.cend() );
155 std::cout << "Projection basis=[ " << hull.affine_basis.basis()[0]
156 << "," << hull.affine_basis.basis()[1] << " ]"
157 << " d=" << hull.projected_dilation << "\n";
158 }
159 else if ( hull.affine_dimension == 3 ) // 3D
160 {
161 // facets are polygons
162 psBoundary2 = polyscope::registerSurfaceMesh("Convex hull bdy dim=2",
163 hull.positions, hull.facets );
164 }
165 std::cout << " dilation=" << hull.projected_dilation
166 << " => counting of lattice points is "
167 << (hull.projected_dilation == 1 ? "correct" : "INCORRECT") << ".\n";
168 std::cout << " #(P ∩ Z3)=" << hull.count() << "\n";
169 std::cout << "#(Int(P) ∩ Z3)=" << hull.countInterior() << "\n";
170 std::cout << " #(Bd(P) ∩ Z3)=" << hull.countBoundary() << "\n";
171 polyscope::show();
172 return EXIT_SUCCESS;
173
174}
polyscope::PointCloud * psBoundary0
polyscope::SurfaceMesh * psBoundary2
polyscope::PointCloud * psPoints
polyscope::CurveNetwork * psBoundary1
static std::vector< Point > makeRandomLatticePointsFromDirVectors(Point A, const std::vector< Point > &V, int nb, double radius, int amplitude, int aff_dim)
polyscope::PointCloud * psVertices
Aim: Implements the quickhull algorithm by Barber et al. , a famous arbitrary dimensional convex hull...
MyPointD Point

References DGtal::L, makeRandomLatticePointsFromDirVectors(), psBoundary0, psBoundary1, psBoundary2, psPoints, and psVertices.

◆ makeRandomLatticePointsFromDirVectors()

template<typename Point >
static std::vector< Point > makeRandomLatticePointsFromDirVectors ( Point  A,
const std::vector< Point > &  V,
int  nb,
double  radius,
int  amplitude,
int  aff_dim 
)
static

Definition at line 88 of file exampleGenericLatticeConvexHull3D.cpp.

90{
91 std::uniform_int_distribution<int> U(-amplitude, amplitude);
92 std::vector< Point > P;
93 int m = std::min( aff_dim, (int) V.size() );
94 for ( auto k = 0; P.size() < nb && k < 100000; k++ )
95 {
96 Point B = A;
97 for ( auto i = 0; i < m; i++ )
98 {
99 int l = U( g );
100 B += l * V[ i ];
101 }
102 if ( (B-A).norm() <= radius )
103 P.push_back( B );
104 }
105 std::shuffle( P.begin(), P.end(), g );
106 return P;
107}
std::mt19937 g(rd())

References g().

Referenced by main().

Variable Documentation

◆ psBoundary0

◆ psBoundary1

◆ psBoundary2

◆ psPoints

◆ psVertices

polyscope::PointCloud* psVertices

Definition at line 77 of file exampleGenericLatticeConvexHull3D.cpp.

Referenced by main().

◆ rd

std::random_device rd

Definition at line 82 of file exampleGenericLatticeConvexHull3D.cpp.