DGtal 1.4.0
Loading...
Searching...
No Matches
exampleQuickHull3D.cpp File Reference
#include "DGtal/base/Common.h"
#include "DGtal/geometry/volumes/ConvexityHelper.h"
#include "DGtal/shapes/SurfaceMesh.h"
#include "DGtal/io/writers/SurfaceMeshWriter.h"
#include "ConfigExamples.h"
Include dependency graph for exampleQuickHull3D.cpp:

Go to the source code of this file.

Functions

double rand01 ()
 [QuickHull3D-Includes]
 
int main (int argc, char *argv[])
 

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
2021/01/01

An example file named exampleQuickHull3D.

This file is part of the DGtal library.

Definition in file exampleQuickHull3D.cpp.

Function Documentation

◆ main()

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

Definition at line 87 of file exampleQuickHull3D.cpp.

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}
Aim: Represents a polygon mesh, i.e. a 2-dimensional combinatorial surface whose faces are (topologic...
std::ostream & info()
Space::RealPoint RealPoint
Definition StdDefs.h:170
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

References DGtal::Trace::info(), DGtal::trace, and DGtal::SurfaceMeshWriter< TRealPoint, TRealVector >::writeOBJ().

◆ rand01()

double rand01 ( )

[QuickHull3D-Includes]

[QuickHull3D-Includes]

Examples
geometry/tools/exampleQuickHull3D.cpp, and geometry/tools/exampleRationalBallQuickHull3D.cpp.

Definition at line 83 of file exampleQuickHull3D.cpp.

83{ return (double) rand() / (double) RAND_MAX; }