DGtal 2.1.0
Loading...
Searching...
No Matches
exampleGenericLatticeConvexHull.cpp File Reference
#include <iostream>
#include <format>
#include <vector>
#include "DGtal/base/Common.h"
#include "DGtal/kernel/SpaceND.h"
#include "DGtal/geometry/tools/GenericLatticeConvexHull.h"
Include dependency graph for exampleGenericLatticeConvexHull.cpp:

Go to the source code of this file.

Functions

int main ()
 [generic-qhull-example]
 

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/08

An example file named exampleGenericLatticeConvexHull

This file is part of the DGtal library.

Definition in file exampleGenericLatticeConvexHull.cpp.

Function Documentation

◆ main()

int main ( void  )

[generic-qhull-example]

Definition at line 53 of file exampleGenericLatticeConvexHull.cpp.

54{
55 using namespace DGtal;
56 // point coordinates are int (32 bits), computation with int64_t
59 typedef Space::Point Point;
60
61 std::size_t nb_ok [ 4 ] = { 0, 0, 0, 0 };
62 std::size_t nb_per_dim [ 4 ] = { 0, 0, 0, 0 };
63 std::size_t nb_vertices[ 4 ] = { 0, 0, 0, 0 };
64 std::size_t nb_facets [ 4 ] = { 0, 0, 0, 0 };
65 const std::size_t nb = 100000;
66 for ( std::size_t n = 0; n < nb; ++n )
67 {
68 // Create a random set of points
69 std::vector< Point > X;
70 int m = 2 + rand() % 9; //< number of points
71 for ( int i = 0; i < m; i++ )
72 X.push_back( Point{ rand() % 8, rand() % 8, rand() % 8 } );
73 // Compute convex hull
74 QHull hull;
75 bool ok = hull.compute( X );
76 auto k = hull.affine_dimension; // affine dimension of X
77 // Compute statistics
78 nb_per_dim [ k ] += 1;
79 nb_vertices[ k ] += hull.positions.size(); // positions of vertices
80 nb_facets [ k ] += hull.facets.size(); // facets
81 nb_ok [ k ] += ok ? 1 : 0;
82 }
83 for ( auto k = 0; k < 4; k++ )
84 std::cout << std::setprecision(3) << ( 100.0 * nb_per_dim[ k ] ) / nb
85 << "% are " << k << "-dimensional"
86 << " #V=" << double( nb_vertices[ k ] ) / nb_per_dim[ k ]
87 << " #F=" << double( nb_facets[ k ] ) / nb_per_dim[ k ]
88 << " (" << nb_ok[ k ] << "/" << nb_per_dim[ k ] << ")\n";
89 return 0;
90}
DGtal is the top-level namespace which contains all DGtal functions and types.
Aim: Implements the quickhull algorithm by Barber et al. , a famous arbitrary dimensional convex hull...
MyPointD Point