DGtal 1.4.0
Loading...
Searching...
No Matches
exampleQuickHull3D.cpp
Go to the documentation of this file.
1
75#include "DGtal/base/Common.h"
77#include "DGtal/geometry/volumes/ConvexityHelper.h"
79#include "DGtal/shapes/SurfaceMesh.h"
80#include "DGtal/io/writers/SurfaceMeshWriter.h"
81#include "ConfigExamples.h"
82
83double rand01() { return (double) rand() / (double) RAND_MAX; }
84
85using namespace DGtal;
86using namespace DGtal::Z3i;
87int main( int argc, char* argv[] )
88{
89 int precision = argc > 1 ? atoi( argv[ 1 ] ) : 100; // precision
90 std::string inputFilename = argc > 2
91 ? std::string( argv[ 2 ] )
92 : examplesPath + "samples/bunny.dat" ;
93 std::vector< RealPoint > points;
94 std::ifstream finput( inputFilename.c_str() );
95 std::string linestr;
96 while ( std::getline( finput, linestr ) )
97 {
98 std::istringstream iss( linestr );
99 double a, b, c;
100 if ( ! (iss >> a >> b >> c) ) break;
101 points.push_back( RealPoint( a, b, c ) );
102 }
103 trace.info() << "Read " << points.size() << " 3D points." << std::endl;
104
105 // Build rational polytope
106 typedef ConvexityHelper< 3 > Helper;
107 const auto polytope
108 = Helper::computeRationalPolytope( points, precision );
109 trace.info() << polytope << std::endl;
110
111 // Build the boundary of the convex hull as a surface mesh
113 bool ok = Helper::computeConvexHullBoundary( mesh, points, precision );
114 trace.info() << mesh << std::endl;
115 std::ofstream out( "qhull-mesh.obj" );
117 out.close();
118
119 // Build the boundary of the convex hull as a polygonal surface
121 bool ok2 = Helper::computeConvexHullBoundary( polysurf, points, precision );
122 trace.info() << polysurf << std::endl;
123
124 // Build the convex hull as a convex cell complex.
126 bool ok3 = Helper::computeConvexHullCellComplex( cvx_complex, points, precision );
127 trace.info() << cvx_complex << std::endl;
128
129 return ( ok && ok2 && ok3 ) ? 0 : 1;
130}
131
Aim: Represents a polygon mesh, i.e. a 2-dimensional combinatorial surface whose faces are (topologic...
std::ostream & info()
double rand01()
[QuickHull3D-Includes]
Z3i this namespace gathers the standard of types for 3D imagery.
Space::RealPoint RealPoint
Definition StdDefs.h:170
DGtal is the top-level namespace which contains all DGtal functions and types.
Trace trace
Definition Common.h:153
Aim: represents a d-dimensional complex in a d-dimensional space with the following properties and re...
Aim: Provides a set of functions to facilitate the computation of convex hulls and polytopes,...
static bool writeOBJ(std::ostream &output, const SurfaceMesh &smesh)
Aim: Represents an embedded mesh as faces and a list of vertices. Vertices may be shared among faces ...
Definition SurfaceMesh.h:92
int main()
Definition testBits.cpp:56