DGtal  1.5.beta
tangency-reconstruction.cpp
1 
25 #include <iostream>
26 #include <string>
27 #include <DGtal/base/Common.h>
28 #include <DGtal/helpers/StdDefs.h>
29 #include <DGtal/helpers/Shortcuts.h>
30 #include <DGtal/helpers/ShortcutsGeometry.h>
31 #include <DGtal/shapes/SurfaceMesh.h>
32 #include <DGtal/geometry/volumes/TangencyComputer.h>
33 #include <DGtal/geometry/volumes/ConvexityHelper.h>
34 #include <DGtal/geometry/tools/QuickHull.h>
35 #include <DGtal/geometry/surfaces/estimation/PlaneProbingTetrahedronEstimator.h>
36 
37 #include <polyscope/pick.h>
38 #include <polyscope/polyscope.h>
39 #include <polyscope/surface_mesh.h>
40 #include <polyscope/point_cloud.h>
41 
42 
43 using namespace DGtal;
44 using namespace Z3i;
45 
46 // Using standard 3D digital space.
49 // The following typedefs are useful
52 typedef SurfMesh::Face Face;
53 typedef SurfMesh::Vertex Vertex;
56 typedef std::size_t Size;
57 //typedef Z3i::KSpace::SCell SCell;
58 //Polyscope global
59 polyscope::SurfaceMesh *psMesh;
60 polyscope::SurfaceMesh *psDualMesh;
61 polyscope::SurfaceMesh *psOuterDualMesh;
62 polyscope::SurfaceMesh *psTriangle;
63 polyscope::SurfaceMesh *psDelaunay;
64 polyscope::PointCloud* psCloud;
65 polyscope::PointCloud* psCloudRemaining;
66 polyscope::PointCloud* psCloudCvx;
67 polyscope::PointCloud* psCloudProc;
68 polyscope::PointCloud* psOuterCloudProc;
69 polyscope::PointCloud* psCloudInnerDelaunay;
70 
72 SurfMesh dual_surfmesh;
73 SurfMesh outer_dual_surfmesh;
74 float gridstep = 1.0;
75 int vertex_idx = -1;
76 float Time = 0.0;
77 int nb_cones = 10;
78 bool remove_empty_cells = false;
79 bool local_tangency = false;
80 
81 Size current = 0;
82 Size outer_current = 0;
83 std::vector< Size > order;
84 std::vector< Size > outer_order;
85 std::vector< Point > digital_points;
86 std::vector< Point > outer_digital_points;
87 std::set< Point > useless_points;
88 std::vector< double > intercepts;
89 std::vector< double > outer_intercepts;
90 std::vector<std::vector<SH3::SurfaceMesh::Vertex>> dual_faces;
91 std::vector<RealPoint> dual_positions;
92 std::vector<std::vector<SH3::SurfaceMesh::Vertex>> outer_dual_faces;
93 std::vector<RealPoint> outer_dual_positions;
94 std::vector<SCell> surfels;
95 std::map<Point,Size> surfel2idx;
96 
97 KSpace K;
103 
105 typedef std::vector< Index > Indices;
106 typedef double Scalar;
107 typedef std::vector< Scalar > Scalars;
108 
109 int pickPoint()
110 {
111  if ( order.size() != digital_points.size() )
112  {
113  auto rng = std::default_random_engine {};
114  order.resize( digital_points.size() );
115  for ( Size i = 0; i < order.size(); i++ ) order[ i ] = i;
116  std::shuffle( order.begin(), order.end(), rng);
117  current = 0;
118  }
119  if ( current == order.size() ) current = 0;
120  return static_cast<int>(order[ current++ ]);
121 }
122 int pickOuterPoint()
123 {
124  if ( outer_order.size() != outer_digital_points.size() )
125  {
126  auto rng = std::default_random_engine {};
127  outer_order.resize( outer_digital_points.size() );
128  for ( Size i = 0; i < outer_order.size(); i++ ) outer_order[ i ] = i;
129  std::shuffle( outer_order.begin(), outer_order.end(), rng);
130  outer_current = 0;
131  }
132  if ( outer_current == outer_order.size() ) outer_current = 0;
133  return static_cast<int>(outer_order[ outer_current++ ]);
134 }
135 
136 // ----------------------------------------------------------------------
137 // utilities pointel
138 Point pointelRealPoint2Point( RealPoint p )
139 {
140  RealPoint sp = RealPoint( round( p[ 0 ] / gridstep + 0.5 ),
141  round( p[ 1 ] / gridstep + 0.5 ),
142  round( p[ 2 ] / gridstep + 0.5 ) );
143  return Point( sp[ 0 ], sp[ 1 ], sp[ 2 ] );
144 }
145 RealPoint pointelPoint2RealPoint( Point q )
146 {
147  return RealPoint( gridstep * ( q[ 0 ] - 0.5 ),
148  gridstep * ( q[ 1 ] - 0.5 ),
149  gridstep * ( q[ 2 ] - 0.5 ) );
150 }
151 void embedPointels( const std::vector< Point >& vq, std::vector< RealPoint >& vp )
152 {
153  vp.resize( vq.size() );
154  for ( auto i = 0; i < vp.size(); ++i )
155  vp[ i ] = pointelPoint2RealPoint( vq[ i ] );
156 }
157 void digitizePointels( const std::vector< RealPoint >& vp, std::vector< Point >& vq )
158 {
159  vq.resize( vp.size() );
160  for ( auto i = 0; i < vq.size(); ++i )
161  vq[ i ] = pointelRealPoint2Point( vp[ i ] );
162 }
163 
164 // ----------------------------------------------------------------------
165 // utilities voxel
166 Point voxelRealPoint2Point( RealPoint p )
167 {
168  RealPoint sp = RealPoint( round( p[ 0 ] / gridstep ),
169  round( p[ 1 ] / gridstep ),
170  round( p[ 2 ] / gridstep ) );
171  return Point( sp[ 0 ], sp[ 1 ], sp[ 2 ] );
172 }
173 RealPoint voxelPoint2RealPoint( Point q )
174 {
175  return RealPoint( gridstep * ( q[ 0 ] ),
176  gridstep * ( q[ 1 ] ),
177  gridstep * ( q[ 2 ] ) );
178 }
179 void embedVoxels( const std::vector< Point >& vq, std::vector< RealPoint >& vp )
180 {
181  vp.resize( vq.size() );
182  for ( auto i = 0; i < vp.size(); ++i )
183  vp[ i ] = voxelPoint2RealPoint( vq[ i ] );
184 }
185 void digitizeVoxels( const std::vector< RealPoint >& vp, std::vector< Point >& vq )
186 {
187  vq.resize( vp.size() );
188  for ( auto i = 0; i < vq.size(); ++i )
189  vq[ i ] = voxelRealPoint2Point( vp[ i ] );
190 }
191 
192 // @return the indices of all the points of X different from p that are cotangent to p.
193 Indices tangentCone( const Point& p )
194 {
195  return TC.getCotangentPoints( p );
196 }
197 
198 Scalars distances( const Point& p, const Indices& idx )
199 {
200  Scalars D( idx.size() );
201  for ( Index i = 0; i < idx.size(); i++ )
202  D[ i ] = ( TC.point( i ) - p ).norm();
203  return D;
204 }
205 
207 std::vector< Point >
208 findCorners( const std::unordered_set< Point >& S,
209  const std::vector< Vector >& In,
210  const std::vector< Vector >& Out )
211 {
212  std::vector< Point > C;
213  for ( auto&& p : S )
214  {
215  bool corner = true;
216  for ( auto&& n : In )
217  if ( ! S.count( p+n ) ) { corner = false; break; }
218  if ( ! corner ) continue;
219  for ( auto&& n : Out )
220  if ( S.count( p+n ) ) { corner = false; break; }
221  if ( corner ) C.push_back( p );
222  }
223  return C;
224 }
225 
226 void computeQuadrant( int q,
227  std::vector< Vector >& In,
228  std::vector< Vector >& Out )
229 {
230  In.clear();
231  Out.clear();
232  In.push_back( Vector( q & 0x1 ? 1 : -1, 0, 0 ) );
233  In.push_back( Vector( 0, q & 0x2 ? 1 : -1, 0 ) );
234  In.push_back( Vector( 0, 0, q & 0x4 ? 1 : -1 ) );
235  Out.push_back( Vector( q & 0x1 ? 1 : -1, q & 0x2 ? 1 : -1, q & 0x4 ? 1 : -1 ) );
236  Vector D = ( In[ 1 ] - In[ 0 ] ).crossProduct( In[ 2 ] - In[ 0 ] );
237  if ( D.dot( Out[ 0 ] ) < 0.0 ) std::swap( In[ 1 ], In[ 2 ] );
238  In.push_back( In[ 0 ]+In[ 1 ] );
239  In.push_back( In[ 0 ]+In[ 2 ] );
240  In.push_back( In[ 1 ]+In[ 2 ] );
241 }
242 
243 struct UnorderedPointSetPredicate
244 {
245  typedef DGtal::Z3i::Point Point;
246  const std::unordered_set< Point >* myS;
247  explicit UnorderedPointSetPredicate( const std::unordered_set< Point >& S )
248  : myS( &S ) {}
249  bool operator()( const Point& p ) const
250  { return myS->count( p ) != 0; }
251 };
252 
253 void computePlanes()
254 {
255  // - mode specifies the candidate set, it is one of { ProbingMode::H, ProbingMode::R, ProbingMode::R1 }.
256  using Estimator
257  = DGtal::PlaneProbingTetrahedronEstimator< UnorderedPointSetPredicate,
258  ProbingMode::R >;
259  trace.beginBlock( "Compute planes" );
260  std::vector< RealPoint > positions;
261  std::vector< Point > vertices;
262  std::vector< std::vector<SH3::SurfaceMesh::Vertex> > faces;
263  std::vector< Vector > In;
264  std::vector< Vector > Out;
265  std::unordered_set< Point > S( digital_points.cbegin(), digital_points.cend() );
266  UnorderedPointSetPredicate predS( S );
267  Index i = 0;
268  for ( int q = 0; q < 8; q++ ) // for each quadrant
269  {
270  computeQuadrant( q, In, Out );
271  std::vector< Point > corners = findCorners( S, In, Out );
272  std::cout << "Found " << corners.size() << " in Q" << q << std::endl;
273  std::array<Point, 3> m = { In[ 0 ], In[ 1 ], In[ 2 ] };
274  for ( auto&& p : corners )
275  {
276  Estimator estimator( p, m, predS );
277  // if ( estimator.hexagonState()
278  // != Estimator::Neighborhood::HexagonState::Planar)
279  // continue;
280  std::vector<SH3::SurfaceMesh::Vertex> triangle { i, i+1, i+2 };
281  auto v = estimator.vertices();
282  faces.push_back( triangle );
283  vertices.push_back( v[ 0 ] );
284  vertices.push_back( v[ 1 ] );
285  vertices.push_back( v[ 2 ] );
286  while (estimator.advance().first) {
287  auto state = estimator.hexagonState();
288  if (state == Estimator::Neighborhood::HexagonState::Planar) {
289  auto v = estimator.vertices();
290  if ( S.count( v[ 0 ] ) && S.count( v[ 1 ] ) && S.count( v[ 2 ] ) )
291  {
292  std::vector< Point > X { v[ 0 ], v[ 1 ], v[ 2 ] };
293  auto P = dconv.makePolytope( X );
294  if ( dconv.isFullySubconvex( P, TC.cellCover() ) )
295  // // && TC.arePointsCotangent( v[ 0 ], v[ 1 ], v[ 2 ] ) )
296  {
297  vertices[ i ] = v[ 0 ];
298  vertices[ i+1 ] = v[ 1 ];
299  vertices[ i+2 ] = v[ 2 ];
300  }
301  }
302  }
303  }
304  i += 3;
305  // }
306  }
307  }
308  Time = trace.endBlock();
309  embedPointels( vertices, positions );
310  psTriangle = polyscope::registerSurfaceMesh("Triangle", positions, faces);
311 }
312 
313 void displayMidReconstruction()
314 {
315  std::vector< RealPoint > mid_dual_positions( dual_positions.size() );
316  for ( auto i = 0; i < surfels.size(); i++ )
317  {
318  auto int_vox = K.interiorVoxel( surfels[ i ] );
319  auto ext_vox = K.exteriorVoxel( surfels[ i ] );
320  auto int_p = voxelPoint2RealPoint( int_vox ); //K.uCoords( int_vox ) );
321  auto ext_p = voxelPoint2RealPoint( ext_vox ); //K.uCoords( ext_vox ) );
322  const double s = intercepts[ i ] + 0.001;
323  const RealPoint q = ( 1.0 - s ) * int_p + s * ext_p;
324  const double ds = outer_intercepts[ i ] + 0.001;
325  const RealPoint dq = ( 1.0 - ds ) * ext_p + ds * int_p;
326  mid_dual_positions[ i ] = 0.5*(q+dq);
327  }
328  psDualMesh = polyscope::registerSurfaceMesh("Mid surface", mid_dual_positions, dual_faces);
329 }
330 
331 void displayReconstruction()
332 {
333  for ( auto i = 0; i < surfels.size(); i++ )
334  {
335  auto int_vox = K.interiorVoxel( surfels[ i ] );
336  auto ext_vox = K.exteriorVoxel( surfels[ i ] );
337  auto int_p = voxelPoint2RealPoint( int_vox ); //K.uCoords( int_vox ) );
338  auto ext_p = voxelPoint2RealPoint( ext_vox ); //K.uCoords( ext_vox ) );
339  const double s = intercepts[ i ] + 0.001;
340  const RealPoint q = ( 1.0 - s ) * int_p + s * ext_p;
341  dual_positions[ i ] = q;
342  }
343  psDualMesh = polyscope::registerSurfaceMesh("Reconstruction surface", dual_positions, dual_faces);
344  std::vector< Point > X;
345  for ( Size i = 0; i < current; i++ )
346  X.push_back( digital_points[ order[ i ] ] );
347  std::vector< RealPoint > emb_X;
348  embedVoxels( X, emb_X );
349  psCloudProc = polyscope::registerPointCloud( "Processed points", emb_X );
350  psCloudProc->setPointRadius( gridstep / 300.0 );
351 
352 }
353 
354 void displayRemainingPoints()
355 {
356  std::vector< Point > X;
357  std::vector< RealPoint > emb_X;
358  for ( auto i = current; i < digital_points.size(); i++ )
359  {
360  Point p = digital_points[ order[ i ] ];
361  if ( ! useless_points.count( p ) )
362  X.push_back( p );
363  }
364  embedVoxels( X, emb_X );
365  psCloudProc = polyscope::registerPointCloud( "Remaining points", emb_X );
366  psCloudProc->setPointRadius( gridstep / 200.0 );
367 }
368 
369 void displayOuterReconstruction()
370 {
371  for ( auto i = 0; i < surfels.size(); i++ )
372  {
373  auto int_vox = K.interiorVoxel( surfels[ i ] );
374  auto ext_vox = K.exteriorVoxel( surfels[ i ] );
375  auto int_p = voxelPoint2RealPoint( int_vox ); //K.uCoords( int_vox ) );
376  auto ext_p = voxelPoint2RealPoint( ext_vox ); //K.uCoords( ext_vox ) );
377  const double s = outer_intercepts[ i ] + 0.001;
378  const RealPoint q = ( 1.0 - s ) * ext_p + s * int_p;
379  outer_dual_positions[ i ] = q;
380  }
381  psOuterDualMesh = polyscope::registerSurfaceMesh("Reconstruction outer surface", outer_dual_positions, outer_dual_faces);
382  std::vector< Point > X;
383  for ( Size i = 0; i < outer_current; i++ )
384  X.push_back( outer_digital_points[ outer_order[ i ] ] );
385  std::vector< RealPoint > emb_X;
386  embedVoxels( X, emb_X );
387  psOuterCloudProc = polyscope::registerPointCloud( "Processed outer points", emb_X );
388  psOuterCloudProc->setPointRadius( gridstep / 300.0 );
389 
390 }
391 
392 void updateReconstructionFromCells( const std::vector< Point >& X,
393  const std::vector< Point >& cells )
394 {
395  // Compute plane
396  const Vector N = ( X[ 1 ] - X[ 0 ] ).crossProduct( X[ 2 ] - X[ 0 ] );
397  const auto a = N.dot( X[ 0 ] );
398  // Extract 1-cells which are dual to surfels
399  for ( auto&& kp : cells )
400  {
401  // Look for dimension 1 cells.
402  const Cell c = K.uCell( kp );
403  if ( K.uDim( c ) != 1 ) continue;
404  // Compute dual surfel
405  const Dimension t = *K.uDirs( c );
406  const Cell p0 = K.uIncident( c, t, false );
407  const Cell p1 = K.uIncident( c, t, true );
408  const Point dual_kp = K.uCoords( p0 ) + K.uCoords( p1 ) + Point::diagonal(1);
409  const auto it = surfel2idx.find( dual_kp );
410  if ( it == surfel2idx.cend() ) continue;
411  // Compute and update intercept
412  const Size idx = it->second;
413  const SCell surfel = surfels[ idx ];
414  const Point int_vox = K.interiorVoxel( surfel );
415  const Point ext_vox = K.exteriorVoxel( surfel );
416  const auto int_val = N.dot( int_vox );
417  const auto ext_val = N.dot( ext_vox );
418  // std::cout << " int_val=" << int_val << " a=" << a << " ext_val=" << ext_val;
419  if ( ( int_val <= a && ext_val <= a ) || ( int_val >= a && ext_val >= a ) )
420  {
421  if ( ( int_val < a && ext_val < a ) || ( int_val > a && ext_val > a ) )
422  trace.warning() << "Bad intersection" << std::endl;
423  continue;
424  }
425  const double s = (double)( a - int_val ) / (double) (ext_val - int_val );
426  const double old_s = intercepts[ idx ];
427  // if ( old_s < s ) std::cout << " s=" << old_s << " -> " << s << std::endl;
428  intercepts[ idx ] = std::max( old_s, s );
429  }
430 }
431 
432 void updateOuterReconstructionFromCells( const std::vector< Point >& X,
433  const std::vector< Point >& cells )
434 {
435  // Compute plane
436  const Vector N = ( X[ 1 ] - X[ 0 ] ).crossProduct( X[ 2 ] - X[ 0 ] );
437  const auto a = N.dot( X[ 0 ] );
438  // Extract 1-cells which are dual to surfels
439  for ( auto&& kp : cells )
440  {
441  // Look for dimension 1 cells.
442  const Cell c = K.uCell( kp );
443  if ( K.uDim( c ) != 1 ) continue;
444  // Compute dual surfel
445  const Dimension t = *K.uDirs( c );
446  const Cell p0 = K.uIncident( c, t, false );
447  const Cell p1 = K.uIncident( c, t, true );
448  const Point dual_kp = K.uCoords( p0 ) + K.uCoords( p1 ) + Point::diagonal(1);
449  const auto it = surfel2idx.find( dual_kp );
450  if ( it == surfel2idx.cend() ) continue;
451  // Compute and update intercept
452  const Size idx = it->second;
453  const SCell surfel = surfels[ idx ];
454  const Point int_vox = K.interiorVoxel( surfel );
455  const Point ext_vox = K.exteriorVoxel( surfel );
456  const auto int_val = N.dot( int_vox );
457  const auto ext_val = N.dot( ext_vox );
458  if ( ( int_val <= a && ext_val <= a ) || ( int_val >= a && ext_val >= a ) )
459  {
460  if ( ( int_val < a && ext_val < a ) || ( int_val > a && ext_val > a ) )
461  trace.warning() << "Bad intersection" << std::endl;
462  continue;
463  }
464  const double s = (double)( a - ext_val ) / (double) (int_val - ext_val );
465  outer_intercepts[ idx ] = std::max( outer_intercepts[ idx ], s );
466  }
467 }
468 
469 void updateReconstructionFromCells( Point x, Point y,
470  const std::vector< Point >& cells )
471 {
472  // Compute plane
473  const Vector U = y - x;
474  if ( U == Vector::zero ) return;
475  const double l = U.norm();
476  const RealVector u = U.getNormalized();
477  const RealPoint p( x[ 0 ], x[ 1 ], x[ 2 ] );
478  // Extract 1-cells which are dual to surfels
479  for ( auto&& kp : cells )
480  {
481  // Look for dimension 1 cells.
482  const Cell c = K.uCell( kp );
483  if ( K.uDim( c ) != 1 ) continue;
484  // Compute dual surfel
485  const Dimension r = *K.uDirs( c );
486  const Cell p0 = K.uIncident( c, r, false );
487  const Cell p1 = K.uIncident( c, r, true );
488  const Point dual_kp = K.uCoords( p0 ) + K.uCoords( p1 ) + Point::diagonal(1);
489  const auto it = surfel2idx.find( dual_kp );
490  if ( it == surfel2idx.cend() ) continue;
491  // Compute and update intercept
492  const Size idx = it->second;
493  const SCell surfel = surfels[ idx ];
494  const Point int_vox = K.interiorVoxel( surfel );
495  const Point ext_vox = K.exteriorVoxel( surfel );
496  const Vector V = ext_vox - int_vox;
497  const RealVector v = V.getNormalized();
498  const RealPoint q( int_vox[ 0 ], int_vox[ 1 ], int_vox[ 2 ] );
499  const auto uv = u.dot( v );
500  if ( uv == 0 ) continue;
501  // Solving system to get closest points.
502  const auto c1 = ( q - p ).dot( u );
503  const auto c2 = ( p - q ).dot( v );
504  const double d = 1.0-uv*uv;
505  const double s = ( c1 + uv * c2 ) / d; // on [xy]
506  const double t = ( c2 + uv * c1 ) / d; // on linel
507  if ( ( s < 0.0 ) || ( s > l ) ) continue;
508  // if ( ( t < 0.0 ) || ( t > 1.0 ) ) continue;
509  intercepts[ idx ] = std::max( intercepts[ idx ], std::min( t, 1.0 ) );
510  }
511 }
512 
513 void updateReconstructionFromTangentConeLines( int vertex_idx )
514 {
515  if ( digital_points.empty() ) return;
516  if ( vertex_idx < 0 || vertex_idx >= digital_points.size() ) return;
517  const auto p = digital_points[ vertex_idx ];
518  // trace.beginBlock( "Compute tangent cone" );
519  auto local_X_idx = TC.getCotangentPoints( p );
520  std::vector< Point > local_X;
521  for ( auto idx : local_X_idx )
522  local_X.push_back( TC.point( idx ) );
523  for ( auto&& q : local_X )
524  {
525  std::vector< Point > X { p, q };
526  const auto line_cells = dconv.StarCvxH( X, LS.axis() );
527  updateReconstructionFromCells( p, q, line_cells.toPointRange() );
528  }
529 }
530 
531 void updateReconstructionFromTangentConeTriangles( int vertex_idx )
532 {
533  typedef QuickHull3D::IndexRange IndexRange;
534  if ( digital_points.empty() ) return;
535  if ( vertex_idx < 0 || vertex_idx >= digital_points.size() ) return;
536  const auto p = digital_points[ vertex_idx ];
537  // trace.beginBlock( "Compute tangent cone" );
538  auto local_X_idx = TC.getCotangentPoints( p );
539  local_X_idx.push_back( vertex_idx );
540  std::vector< Point > local_X;
541  for ( auto idx : local_X_idx )
542  local_X.push_back( TC.point( idx ) );
543  QuickHull3D hull;
544  hull.setInput( local_X, false );
545  hull.computeConvexHull();
546  std::vector< Point > positions;
547  hull.getVertexPositions( positions );
548 
549  std::vector< IndexRange > facet_vertices;
550  bool ok = hull.getFacetVertices( facet_vertices );
551  if ( ! ok ) trace.error() << "Bad facet computation" << std::endl;
552  // Update from all cones
553  std::set< std::pair< Point, Point > > edges;
554  for ( auto&& facet : facet_vertices )
555  {
556  const auto nb = facet.size();
557  for ( auto i = 0; i < nb; i++ )
558  edges.insert( std::make_pair( positions[ facet[ i ] ],
559  positions[ facet[ (i+1)%nb ] ] ) );
560  }
561  for ( auto&& e : edges )
562  {
563  if ( e.second < e.first ) continue;
564  std::vector< Point > X { p, e.first, e.second };
565  const auto triangle_cells = dconv.StarCvxH( X, LS.axis() );
566  if ( LS.includes( triangle_cells ) ) // tangent to shape
567  updateReconstructionFromCells( X, triangle_cells.toPointRange() );
568  }
569 
570  // Miss features
571  // // Update from all facets
572  // for ( auto&& facet : facet_vertices )
573  // {
574  // const auto nb = facet.size();
575  // for ( auto i = 0; i < nb; i++ )
576  // for ( auto j = i+1; j < nb; j++ )
577  // for ( auto k = j+1; k < nb; k++ )
578  // {
579  // std::vector< Point > X
580  // { positions[ facet[ i ] ],
581  // positions[ facet[ j ] ],
582  // positions[ facet[ k ] ] };
583  // const auto triangle_cells = dconv.StarCvxH( X, LS.axis() );
584  // if ( LS.includes( triangle_cells ) ) // tangent to shape
585  // updateReconstructionFromCells( X, triangle_cells.toPointRange() );
586  // }
587  // }
588 
589  // Too costly
590  // // Update from all possible triplets
591  // const auto nb = positions.size();
592  // for ( auto i = 0; i < nb; i++ )
593  // for ( auto j = i+1; j < nb; j++ )
594  // for ( auto k = j+1; k < nb; k++ )
595  // {
596  // std::vector< Point > X { positions[ i ], positions[ j ], positions[ k ] };
597  // const auto triangle_cells = dconv.StarCvxH( X, LS.axis() );
598  // if ( LS.includes( triangle_cells ) ) // tangent to shape
599  // updateReconstructionFromCells( X, triangle_cells.toPointRange() );
600  // }
601  // Time = trace.endBlock();
602 }
603 
604 
606 void computeTangentCone( int vertex_idx)
607 {
608  if ( digital_points.empty() ) return;
609  // if ( vertex_idx < 0 || vertex_idx >= digital_points.size() ) return;
610  const auto p = digital_points[ vertex_idx ];
611  // trace.beginBlock( "Compute tangent cone" );
612  auto local_X_idx = TC.getCotangentPoints( p );
613  std::cout << "#cone=" << local_X_idx.size() << std::endl;
614  local_X_idx.push_back( vertex_idx );
615  std::vector< Point > local_X;
616  std::vector< RealPoint > emb_local_X;
617  for ( auto idx : local_X_idx )
618  local_X.push_back( TC.point( idx ) );
619  std::vector< double > values( local_X.size(), 0.0 );
620  values.back() = 1.0;
621  embedVoxels( local_X, emb_local_X );
622  psCloud = polyscope::registerPointCloud( "Tangent cone", emb_local_X );
623  psCloud->setPointRadius( gridstep / 300.0 );
624  psCloud->addScalarQuantity( "Classification", values );
625  QuickHull3D hull;
626  hull.setInput( local_X, false );
627  hull.computeConvexHull();
628  std::vector< Point > positions;
629  std::vector< RealPoint > emb_positions;
630  hull.getVertexPositions( positions );
631  // Time = trace.endBlock();
632  embedVoxels( positions, emb_positions );
633  psCloudCvx = polyscope::registerPointCloud( "Tangent cone vertices", emb_positions );
634  psCloudCvx->setPointRadius( gridstep / 200.0 );
635 
636  updateReconstructionFromTangentConeTriangles( vertex_idx );
637  displayReconstruction();
638 }
639 
641 void updateReconstructionFromLocalTangentDelaunayComplex( int vertex_idx)
642 {
643  if ( digital_points.empty() ) return;
644  // if ( vertex_idx < 0 || vertex_idx >= digital_points.size() ) return;
645  const auto p = digital_points[ vertex_idx ];
646  auto local_X_idx = TC.getCotangentPoints( p );
647  local_X_idx.push_back( vertex_idx );
648  std::vector< Point > local_X;
649  for ( auto idx : local_X_idx )
650  local_X.push_back( TC.point( idx ) );
652  if ( local_tangency )
653  local_LS = DGtal::LatticeSetByIntervals< Space >( local_X.cbegin(), local_X.cend(), 0 ).starOfPoints();
654 
657  bool ok = CvxHelper::computeDelaunayCellComplex( dcomplex, local_X, false );
658  if ( ! ok )
659  trace.error() << "Input set of points is not full dimensional." << std::endl;
660 
661  dcomplex.requireFaceGeometry();
662  // Filter cells
663  std::vector< bool > is_cell_tangent( dcomplex.nbCells(), false );
664  if ( remove_empty_cells )
665  for ( Index c = 0; c < dcomplex.nbCells(); ++c )
666  {
667  auto Y = dcomplex.cellVertexPositions( c );
668  auto P = dconv.makePolytope( Y );
669  if ( P.countUpTo( (Integer)Y.size()+1 ) >= (Integer)Y.size()+1 ) continue;
670  is_cell_tangent[ c ] =
671  local_tangency
672  ? dconv.isFullySubconvex( P, local_LS )
673  : dconv.isFullySubconvex( P, LS );
674  }
675  else
676  for ( Index c = 0; c < dcomplex.nbCells(); ++c )
677  {
678  auto Y = dcomplex.cellVertexPositions( c );
679  is_cell_tangent[ c ] =
680  local_tangency
681  ? dconv.isFullySubconvex( Y, local_LS )
682  : dconv.isFullySubconvex( Y, LS );
683  }
684  // Get faces
685  std::vector< std::vector<SH3::SurfaceMesh::Vertex> > del_faces;
686  std::vector< RealPoint > del_positions;
687  std::set<Index> boundary_or_ext_points;
688  for ( Index f = 0; f < dcomplex.nbFaces(); ++f )
689  {
690  auto f0 = std::make_pair( f, false );
691  auto c0 = dcomplex.faceCell( f0 );
692  auto f1 = dcomplex.opposite( f0 );
693  auto c1 = dcomplex.faceCell( f1 );
694  if ( dcomplex.isInfinite( c0 ) )
695  {
696  std::swap( f0, f1 );
697  std::swap( c0, c1 );
698  }
699  if ( ! is_cell_tangent[ c0 ] )
700  {
701  auto V = dcomplex.cellVertices( c0 );
702  boundary_or_ext_points.insert( V.cbegin(), V.cend() );
703  }
704  if ( ! dcomplex.isInfinite( c1 ) && ! is_cell_tangent[ c1 ] )
705  {
706  auto V = dcomplex.cellVertices( c1 );
707  boundary_or_ext_points.insert( V.cbegin(), V.cend() );
708  }
709  bool bdry =
710  ( is_cell_tangent[ c0 ] && ( dcomplex.isInfinite( c1 )
711  || ( ! is_cell_tangent[ c1 ] ) ) )
712  ||
713  ( ! is_cell_tangent[ c0 ] && ( ! dcomplex.isInfinite( c1 )
714  && ( is_cell_tangent[ c1 ] ) ) );
715  if ( ! bdry ) continue;
716  std::vector< Point > X = dcomplex.faceVertexPositions( f0 );
717  const auto triangle_cells = dconv.StarCvxH( X, LS.axis() );
718  //if ( LS.includes( triangle_cells ) ) // tangent to shape
719  updateReconstructionFromCells( X, triangle_cells.toPointRange() );
720  }
721  useless_points.insert( p );
722  for ( Index v = 0; v < dcomplex.nbVertices(); v++ )
723  {
724  auto q = dcomplex.position( v );
725  if ( ! boundary_or_ext_points.count( v ) )
726  useless_points.insert( q );
727  }
728 }
729 
730 
732 void updateOuterReconstructionFromLocalTangentDelaunayComplex( int vertex_idx)
733 {
734  if ( outer_digital_points.empty() ) return;
735  // if ( vertex_idx < 0 || vertex_idx >= digital_points.size() ) return;
736  const auto p = outer_digital_points[ vertex_idx ];
737  auto local_X_idx = outer_TC.getCotangentPoints( p );
738  local_X_idx.push_back( vertex_idx );
739  std::vector< Point > local_X;
740  for ( auto idx : local_X_idx )
741  local_X.push_back( outer_TC.point( idx ) );
742 
745  bool ok = CvxHelper::computeDelaunayCellComplex( dcomplex, local_X, false );
746  if ( ! ok )
747  trace.error() << "Input set of points is not full dimensional." << std::endl;
748 
749  dcomplex.requireFaceGeometry();
750  // Filter cells
751  std::vector< bool > is_cell_tangent( dcomplex.nbCells(), false );
752  if ( remove_empty_cells )
753  for ( Index c = 0; c < dcomplex.nbCells(); ++c )
754  {
755  auto Y = dcomplex.cellVertexPositions( c );
756  auto P = dconv.makePolytope( Y );
757  if ( P.countUpTo( (Integer)Y.size()+1 ) >= (Integer)Y.size()+1 ) continue;
758  is_cell_tangent[ c ] = dconv.isFullySubconvex( P, outer_LS );
759  }
760  else
761  for ( Index c = 0; c < dcomplex.nbCells(); ++c )
762  {
763  auto Y = dcomplex.cellVertexPositions( c );
764  is_cell_tangent[ c ] = dconv.isFullySubconvex( Y, outer_LS );
765  }
766  // Get faces
767  std::vector< std::vector<SH3::SurfaceMesh::Vertex> > del_faces;
768  for ( Index f = 0; f < dcomplex.nbFaces(); ++f )
769  {
770  auto f0 = std::make_pair( f, false );
771  auto c0 = dcomplex.faceCell( f0 );
772  auto f1 = dcomplex.opposite( f0 );
773  auto c1 = dcomplex.faceCell( f1 );
774  if ( dcomplex.isInfinite( c0 ) )
775  {
776  std::swap( f0, f1 );
777  std::swap( c0, c1 );
778  }
779  bool bdry =
780  ( is_cell_tangent[ c0 ] && ( dcomplex.isInfinite( c1 )
781  || ( ! is_cell_tangent[ c1 ] ) ) )
782  ||
783  ( ! is_cell_tangent[ c0 ] && ( ! dcomplex.isInfinite( c1 )
784  && ( is_cell_tangent[ c1 ] ) ) );
785  if ( ! bdry ) continue;
786  std::vector< Point > X = dcomplex.faceVertexPositions( f0 );
787  const auto triangle_cells = dconv.StarCvxH( X, outer_LS.axis() );
788  //if ( LS.includes( triangle_cells ) ) // tangent to shape
789  updateOuterReconstructionFromCells( X, triangle_cells.toPointRange() );
790  }
791 }
792 
794 void computeLocalTangentDelaunayComplex( int vertex_idx)
795 {
796  if ( digital_points.empty() ) return;
797  // if ( vertex_idx < 0 || vertex_idx >= digital_points.size() ) return;
798  const auto p = digital_points[ vertex_idx ];
799  trace.beginBlock( "Compute tangent cone" );
800  auto local_X_idx = TC.getCotangentPoints( p );
801  std::cout << "#cone=" << local_X_idx.size() << std::endl;
802  local_X_idx.push_back( vertex_idx );
803  std::vector< Point > local_X;
804  for ( auto idx : local_X_idx )
805  local_X.push_back( TC.point( idx ) );
807  if ( local_tangency )
808  local_LS = DGtal::LatticeSetByIntervals< Space >( local_X.cbegin(), local_X.cend(), 0 ).starOfPoints();
809 
812  bool ok = CvxHelper::computeDelaunayCellComplex( dcomplex, local_X, false );
813  if ( ! ok )
814  trace.error() << "Input set of points is not full dimensional." << std::endl;
815  dcomplex.requireFaceGeometry();
816 
817  // Filter cells
818  std::vector< bool > is_cell_tangent( dcomplex.nbCells(), false );
819  if ( remove_empty_cells )
820  for ( Index c = 0; c < dcomplex.nbCells(); ++c )
821  {
822  auto Y = dcomplex.cellVertexPositions( c );
823  auto P = dconv.makePolytope( Y );
824  if ( P.countUpTo( (Integer)Y.size()+1 ) >= (Integer)Y.size()+1 ) continue;
825  is_cell_tangent[ c ] =
826  local_tangency
827  ? dconv.isFullySubconvex( P, local_LS )
828  : dconv.isFullySubconvex( P, LS );
829  }
830  else
831  for ( Index c = 0; c < dcomplex.nbCells(); ++c )
832  {
833  auto Y = dcomplex.cellVertexPositions( c );
834  is_cell_tangent[ c ] =
835  local_tangency
836  ? dconv.isFullySubconvex( Y, local_LS )
837  : dconv.isFullySubconvex( Y, LS );
838  }
839  // Get faces
840  std::vector< std::vector<SH3::SurfaceMesh::Vertex> > del_faces;
841  std::vector< RealPoint > del_positions;
842  std::vector< RealPoint > del_inner_points;
843  std::set<Index> boundary_or_ext_points;
844  for ( Index f = 0; f < dcomplex.nbFaces(); ++f )
845  {
846  auto f0 = std::make_pair( f, false );
847  auto c0 = dcomplex.faceCell( f0 );
848  auto f1 = dcomplex.opposite( f0 );
849  auto c1 = dcomplex.faceCell( f1 );
850  if ( dcomplex.isInfinite( c0 ) )
851  {
852  std::swap( f0, f1 );
853  std::swap( c0, c1 );
854  }
855  if ( ! is_cell_tangent[ c0 ] )
856  {
857  auto V = dcomplex.cellVertices( c0 );
858  boundary_or_ext_points.insert( V.cbegin(), V.cend() );
859  }
860  if ( ! dcomplex.isInfinite( c1 ) && ! is_cell_tangent[ c1 ] )
861  {
862  auto V = dcomplex.cellVertices( c1 );
863  boundary_or_ext_points.insert( V.cbegin(), V.cend() );
864  }
865  bool bdry =
866  ( is_cell_tangent[ c0 ] && ( dcomplex.isInfinite( c1 )
867  || ( ! is_cell_tangent[ c1 ] ) ) )
868  ||
869  ( ! is_cell_tangent[ c0 ] && ( ! dcomplex.isInfinite( c1 )
870  && ( is_cell_tangent[ c1 ] ) ) );
871  if ( bdry ) del_faces.push_back( dcomplex.faceVertices( f0 ) );
872  }
873  for ( Index v = 0; v < dcomplex.nbVertices(); v++ )
874  {
875  auto p = dcomplex.position( v );
876  del_positions.push_back( voxelPoint2RealPoint( p ) );
877  if ( ! boundary_or_ext_points.count( v ) )
878  del_inner_points.push_back( voxelPoint2RealPoint( p ) );
879  }
880  Time = trace.endBlock();
881 
882  psCloudInnerDelaunay = polyscope::registerPointCloud( "Inner Delaunay points",
883  del_inner_points );
884  psCloudInnerDelaunay->setPointRadius( gridstep / 200.0 );
885 
886  psDelaunay = polyscope::registerSurfaceMesh("Delaunay surface",
887  del_positions, del_faces);
888  updateReconstructionFromLocalTangentDelaunayComplex( vertex_idx );
889  displayReconstruction();
890 }
891 
893 void computeGlobalTangentDelaunayComplex()
894 {
895  if ( digital_points.empty() ) return;
896  // if ( vertex_idx < 0 || vertex_idx >= digital_points.size() ) return;
897  const auto p = digital_points[ vertex_idx ];
898  trace.beginBlock( "Compute global tangent Delaunay complex" );
899 
902  bool ok = CvxHelper::computeDelaunayCellComplex( dcomplex, digital_points, false );
903  if ( ! ok )
904  trace.error() << "Input set of points is not full dimensional." << std::endl;
905 
906  // Reorder vertices of faces
907  dcomplex.requireFaceGeometry();
908 
909  // Filter cells
910  std::vector< bool > is_cell_tangent( dcomplex.nbCells(), false );
911  for ( Index c = 0; c < dcomplex.nbCells(); ++c )
912  {
913  auto Y = dcomplex.cellVertexPositions( c );
914  auto P = dconv.makePolytope( Y );
915  if ( P.countUpTo( (Integer)Y.size()+1 ) >= (Integer)Y.size()+1 ) continue;
916  is_cell_tangent[ c ] = dconv.isFullySubconvex( P, LS );
917  }
918 
919  // Get faces
920  std::vector< std::vector<SH3::SurfaceMesh::Vertex> > del_faces;
921  std::vector< RealPoint > del_positions;
922  std::vector< RealPoint > del_inner_points;
923  std::set<Index> useful_points;
924  std::vector< std::pair< Index, bool > > all_bdry_faces;
925  for ( Index f = 0; f < dcomplex.nbFaces(); ++f )
926  {
927  auto f0 = std::make_pair( f, false );
928  auto c0 = dcomplex.faceCell( f0 );
929  auto f1 = dcomplex.opposite( f0 );
930  auto c1 = dcomplex.faceCell( f1 );
931  if ( dcomplex.isInfinite( c0 ) )
932  {
933  std::swap( c0, c1 ); std::swap( f0, f1 );
934  }
935  const bool inf_c1 = dcomplex.isInfinite( c1 );
936  const bool bdry_f0 = is_cell_tangent[ c0 ] &&
937  ( inf_c1 || ! is_cell_tangent[ c1 ] );
938  const bool bdry_f1 = ( ! is_cell_tangent[ c0 ] )
939  && ( ( ! inf_c1 ) && is_cell_tangent[ c1 ] );
940 
941  if ( bdry_f0 )
942  {
943  const auto V = dcomplex.faceVertices( f0 );
944  del_faces.push_back( V );
945  useful_points.insert( V.cbegin(), V.cend() );
946  all_bdry_faces.push_back( f0 );
947  }
948  else if ( bdry_f1 )
949  {
950  const auto V = dcomplex.faceVertices( f1 );
951  del_faces.push_back( V );
952  useful_points.insert( V.cbegin(), V.cend() );
953  all_bdry_faces.push_back( f1 );
954  }
955  else if ( ( ! ( is_cell_tangent[ c0 ] || inf_c1 || is_cell_tangent[ c0 ] )
956  || ( ! is_cell_tangent[ c0 ] && inf_c1 ) )
957  && dconv.isFullySubconvex( dcomplex.faceVertexPositions( f0 ), LS ) )
958  {
959  const auto V = dcomplex.faceVertices( f0 );
960  del_faces.push_back( V );
961  useful_points.insert( V.cbegin(), V.cend() );
962  //del_faces.push_back( dcomplex.faceVertices( f1 ) );
963  all_bdry_faces.push_back( f0 );
964  }
965  }
966  for ( Index v = 0; v < dcomplex.nbVertices(); v++ )
967  {
968  auto p = dcomplex.position( v );
969  del_positions.push_back( voxelPoint2RealPoint( p ) );
970  if ( ! useful_points.count( v ) )
971  {
972  del_inner_points.push_back( voxelPoint2RealPoint( p ) );
973  useless_points.insert( p );
974  }
975  }
976  Time = trace.endBlock();
977 
978  trace.beginBlock( "Update reconstruction" );
979  const auto axis = LS.axis();
980  for ( auto f : all_bdry_faces )
981  {
982  std::vector< Point > X = dcomplex.faceVertexPositions( f );
983  const auto cells = dconv.StarCvxH( X, axis );
984  // std::cout << "(" << f.first << "," << f.second << ") "
985  // << " #X=" << X.size() << " #cells=" << cells.size() << std::endl;
986  updateReconstructionFromCells( X, cells.toPointRange() );
987  }
988  trace.endBlock();
989 
990  psCloudInnerDelaunay = polyscope::registerPointCloud( "Inner Delaunay points",
991  del_inner_points );
992  psCloudInnerDelaunay->setPointRadius( gridstep / 200.0 );
993 
994  psDelaunay = polyscope::registerSurfaceMesh("Delaunay surface",
995  del_positions, del_faces);
996  displayReconstruction();
997 }
998 
999 void myCallback()
1000 {
1001  // Select a vertex with the mouse
1002  if (polyscope::pick::haveSelection()) {
1003  bool goodSelection = false;
1004  auto selection = polyscope::pick::getSelection();
1005  auto selectedSurface = static_cast<polyscope::SurfaceMesh*>(selection.first);
1006  const auto idx = selection.second;
1007 
1008  // Only authorize selection on the input surface and the reconstruction
1009  auto surf = polyscope::getSurfaceMesh("Input surface");
1010  goodSelection = goodSelection || (selectedSurface == surf);
1011  const auto nv = selectedSurface->nVertices();
1012  // Validate that it its a face index
1013  if ( goodSelection )
1014  {
1015  if ( idx < nv )
1016  {
1017  std::ostringstream otext;
1018  otext << "Selected vertex = " << idx;
1019  ImGui::Text( "%s", otext.str().c_str() );
1020  vertex_idx = (Integer)idx;
1021  }
1022  else vertex_idx = -1;
1023  }
1024  }
1025  if (ImGui::Button("Compute global tangent Delaunay complex"))
1026  {
1027  computeGlobalTangentDelaunayComplex();
1028  }
1029  if (ImGui::Button("Compute tangent cone"))
1030  {
1031  computeTangentCone( pickPoint() );
1032  current--;
1033  }
1034  if (ImGui::Button("Compute local Delaunay complex"))
1035  {
1036  computeLocalTangentDelaunayComplex( pickPoint() );
1037  current--;
1038  }
1039  ImGui::Checkbox( "Remove empty cells", &remove_empty_cells );
1040  ImGui::Checkbox( "Check local tangency", &local_tangency );
1041  if (ImGui::Button("Compute planes"))
1042  computePlanes();
1043  if (ImGui::Button("Compute reconstructions from triangles"))
1044  { // todo
1045  trace.beginBlock( "Compute reconstruction" );
1046  for ( int i = 0; i < nb_cones; ++i )
1047  {
1048  trace.progressBar( (double) i, (double) nb_cones );
1049  updateReconstructionFromTangentConeTriangles( pickPoint() );
1050  }
1051  displayReconstruction();
1052  Time = trace.endBlock();
1053  }
1054  if (ImGui::Button("Compute reconstruction from local Delaunay cplx"))
1055  { // todo
1056  trace.beginBlock( "Compute reconstruction" );
1057  for ( int i = 0; i < nb_cones; ++i )
1058  {
1059  trace.progressBar( (double) i, (double) nb_cones );
1060  auto j = pickPoint();
1061  if ( ! useless_points.count( digital_points[ j ] ) )
1062  updateReconstructionFromLocalTangentDelaunayComplex( j );
1063  if ( current == 0 ) break;
1064  }
1065  displayReconstruction();
1066  Time = trace.endBlock();
1067  }
1068  if (ImGui::Button("Compute outer reconstruction from local Delaunay cplx"))
1069  { // todo
1070  trace.beginBlock( "Compute outer reconstruction" );
1071  for ( int i = 0; i < nb_cones; ++i )
1072  {
1073  trace.progressBar( (double) i, (double) nb_cones );
1074  updateOuterReconstructionFromLocalTangentDelaunayComplex( pickOuterPoint() );
1075  }
1076  displayOuterReconstruction();
1077  Time = trace.endBlock();
1078  }
1079  if (ImGui::Button("Compute reconstructions from lines"))
1080  { // todo
1081  trace.beginBlock( "Compute reconstruction" );
1082  for ( int i = 0; i < nb_cones; ++i )
1083  {
1084  trace.progressBar( (double) i, (double) nb_cones );
1085  updateReconstructionFromTangentConeLines( pickPoint() );
1086  }
1087  displayReconstruction();
1088  Time = trace.endBlock();
1089  }
1090  ImGui::SliderInt("Nb cones for reconstruction", &nb_cones, 1, 100);
1091  ImGui::Text( "Computation time = %f ms", Time );
1092  ImGui::Text( "#X = %ld, #P = %ld, #U = %ld",
1093  digital_points.size(), current, useless_points.size() );
1094  if (ImGui::Button("Mid reconstructions"))
1095  displayMidReconstruction();
1096  if (ImGui::Button("Remaining points"))
1097  displayRemainingPoints();
1098 }
1099 
1100 
1101 int main( int argc, char* argv[] )
1102 {
1103  auto params = SH3::defaultParameters() | SHG3::defaultParameters() | SHG3::parametersGeometryEstimation();
1104  params("surfaceComponents", "All");
1105  if ( argc <= 1 ) return 0;
1106  bool input_polynomial = argc > 2;
1107  std::string filename = std::string( argv[ 1 ] );
1108  std::string polynomial = std::string( argv[ 1 ] );
1109  const double h = argc >= 3 ? atof( argv[ 2 ] ) : 1.0;
1110  const double bounds = argc >= 4 ? atof( argv[ 3 ] ) : 10.0;
1111 
1113  if ( input_polynomial )
1114  {
1115  params( "polynomial", polynomial );
1116  params( "gridstep", h );
1117  params( "minAABB", -bounds );
1118  params( "maxAABB", bounds );
1119  params( "offset", 1.0 );
1120  params( "closed", 1 );
1121  auto implicit_shape = SH3::makeImplicitShape3D ( params );
1122  auto digitized_shape = SH3::makeDigitizedImplicitShape3D( implicit_shape, params );
1123  K = SH3::getKSpace( params );
1124  binary_image = SH3::makeBinaryImage(digitized_shape,
1126  params );
1127  }
1128  else
1129  {
1130  //params( "thresholdMin", vm[ "min" ].as<int>() );
1131  // params( "thresholdMax", vm[ "max" ].as<int>() );
1132  binary_image = SH3::makeBinaryImage( filename, params );
1133  K = SH3::getKSpace( binary_image );
1134  }
1135 
1136  // auto binary_image = SH3::makeBinaryImage(filename, params );
1137  // K = SH3::getKSpace( binary_image, params );
1138  auto surface = SH3::makeDigitalSurface( binary_image, K, params );
1139  surfels = SH3::getSurfelRange( surface, params );
1140  // compute map surfel -> idx
1141  for ( Size i = 0; i < surfels.size(); i++ )
1142  surfel2idx[ K.sKCoords( surfels[ i ] ) ] = i;
1143  // compute inner points
1144  std::set< Point > I;
1145  for ( auto s : surfels ) I.insert( K.interiorVoxel( s ) );
1146  digital_points = std::vector<Point>( I.cbegin(), I.cend() );
1147  // compute outer points
1148  std::set< Point > O;
1149  for ( auto s : surfels ) O.insert( K.exteriorVoxel( s ) );
1150  outer_digital_points = std::vector<Point>( O.cbegin(), O.cend() );
1151  // initializa intercepts
1152  intercepts = std::vector< double >( surfels.size(), 0.0 );
1153  outer_intercepts = std::vector< double >( surfels.size(), 0.0 );
1154  auto primalSurface = SH3::makePrimalSurfaceMesh(surface);
1155  SH3::Surfel2Index s2i;
1156  auto dualSurface = SH3::makeDualPolygonalSurface( s2i, surface );
1157  //Need to convert the faces
1158  std::vector<std::vector<SH3::SurfaceMesh::Vertex>> faces;
1159  std::vector<RealPoint> positions;
1160  for(auto face= 0 ; face < primalSurface->nbFaces(); ++face)
1161  faces.push_back(primalSurface->incidentVertices( face ));
1162 
1163  //Recasting to vector of vertices
1164  positions = primalSurface->positions();
1165 
1166  surfmesh = SurfMesh(positions.begin(),
1167  positions.end(),
1168  faces.begin(),
1169  faces.end());
1170  for(auto face= 0 ; face < dualSurface->nbFaces(); ++face)
1171  dual_faces.push_back( dualSurface->verticesAroundFace( face ));
1172 
1173  //Recasting to vector of vertices
1174  for ( auto vtx = 0; vtx < dualSurface->nbVertices(); ++vtx )
1175  dual_positions.push_back( dualSurface->position( vtx ) );
1176 
1177  dual_surfmesh = SurfMesh(dual_positions.begin(),
1178  dual_positions.end(),
1179  dual_faces.begin(),
1180  dual_faces.end());
1181  outer_dual_faces = dual_faces;
1182  outer_dual_positions = dual_positions;
1183  outer_dual_surfmesh = SurfMesh(outer_dual_positions.begin(),
1184  outer_dual_positions.end(),
1185  outer_dual_faces.begin(),
1186  outer_dual_faces.end());
1187  std::cout << surfmesh << std::endl;
1188  // Make digital surface
1189  // digitizePointels( positions, digital_points );
1190  trace.info() << "Inner points has " << digital_points.size() << " points." << std::endl;
1191  trace.info() << "Outer points has " << outer_digital_points.size() << " points." << std::endl;
1194  TC.init( digital_points.cbegin(), digital_points.cend() );
1195  outer_TC = DGtal::TangencyComputer< KSpace >( K );
1196  outer_TC.init( outer_digital_points.cbegin(), outer_digital_points.cend() );
1197  trace.info() << "#cell_cover = " << TC.cellCover().nbCells() << std::endl;
1198  trace.info() << "#outer_cell_cover = " << outer_TC.cellCover().nbCells() << std::endl;
1200  ( digital_points.cbegin(), digital_points.cend(), 0 ).starOfPoints();
1201  trace.info() << "#lattice_cover = " << LS.size() << std::endl;
1203  ( outer_digital_points.cbegin(), outer_digital_points.cend(), 0 ).starOfPoints();
1204  trace.info() << "#outer_lattice_cover = " << outer_LS.size() << std::endl;
1205 
1206  // Initialize polyscope
1207  polyscope::init();
1208 
1209  psMesh = polyscope::registerSurfaceMesh("Input surface", positions, faces);
1210  displayReconstruction();
1211  displayOuterReconstruction();
1212  // psDualMesh = polyscope::registerSurfaceMesh("Input dual surface", dual_positions, dual_faces);
1213 
1214 
1215  // Set the callback function
1216  polyscope::state::userCallback = myCallback;
1217  polyscope::show();
1218  return EXIT_SUCCESS;
1219 }
Aim: Smart pointer based on reference counts.
Definition: CountedPtr.h:80
bool isFullySubconvex(const PointRange &Y, const LatticeSet &StarX) const
LatticeSet StarCvxH(const PointRange &X, Dimension axis=dimension) const
LatticePolytope makePolytope(const PointRange &X, bool make_minkowski_summable=false) const
Aim: This class is a model of CCellularGridSpaceND. It represents the cubical grid as a cell complex,...
Point interiorVoxel(const SCell &c) const
For a given surfel ((n-1)-signed cell), returns its interior voxel (point in Z^d given by the direct ...
Cell uIncident(const Cell &c, Dimension k, bool up) const
Return the forward or backward unsigned cell incident to [c] along axis [k], depending on [up].
DirIterator uDirs(const Cell &p) const
Given an unsigned cell [p], returns an iterator to iterate over each coordinate the cell spans.
const Point & sKCoords(const SCell &c) const
Return its Khalimsky coordinates.
Cell uCell(const PreCell &c) const
From an unsigned cell, returns an unsigned cell lying into this Khalismky space.
const Point & lowerBound() const
Return the lower bound for digital points in this space.
Dimension uDim(const Cell &p) const
Return the dimension of the cell [p].
Point uCoords(const Cell &c) const
Return its digital coordinates.
const Point & upperBound() const
Return the upper bound for digital points in this space.
Point exteriorVoxel(const SCell &c) const
For a given surfel ((n-1)-signed cell), returns its exterior voxel (point in Z^d given by the indirec...
bool includes(const Self &other) const
Aim: A class that locally estimates a normal on a digital set using only a predicate "does a point x ...
auto dot(const PointVector< dim, OtherComponent, OtherStorage > &v) const -> decltype(DGtal::dotProduct(*this, v))
Dot product with a PointVector.
double norm(const NormType type=L_2) const
static Self diagonal(Component val=1)
PointVector< dim, double, std::array< double, dim > > getNormalized() const
static Self zero
Static const for zero PointVector.
Definition: PointVector.h:1595
Aim: This class is used to simplify shape and surface creation. With it, you can create new shapes an...
Aim: This class is used to simplify shape and surface creation. With it, you can create new shapes an...
Definition: Shortcuts.h:105
std::map< Surfel, IdxSurfel > Surfel2Index
Definition: Shortcuts.h:188
Aim: A class that computes tangency to a given digital set. It provides services to compute all the c...
std::ostream & error()
void beginBlock(const std::string &keyword="")
std::ostream & info()
void progressBar(const double currentValue, const double maximalValue)
std::ostream & warning()
double endBlock()
SurfaceMesh< RealPoint, RealVector > SurfMesh
void myCallback()
polyscope::SurfaceMesh * psMesh
SurfMesh surfmesh
CountedPtr< SH3::DigitalSurface > surface
CountedPtr< SH3::BinaryImage > binary_image
SMesh::Index Index
Space::RealPoint RealPoint
Definition: StdDefs.h:170
Space::Point Point
Definition: StdDefs.h:168
Space::Vector Vector
Definition: StdDefs.h:169
DGtal::int32_t Integer
Definition: StdDefs.h:143
DGtal is the top-level namespace which contains all DGtal functions and types.
auto crossProduct(PointVector< 3, LeftEuclideanRing, LeftContainer > const &lhs, PointVector< 3, RightEuclideanRing, RightContainer > const &rhs) -> decltype(DGtal::constructFromArithmeticConversion(lhs, rhs))
Cross product of two 3D Points/Vectors.
DGtal::uint32_t Dimension
Definition: Common.h:136
Trace trace
Definition: Common.h:153
std::pair< typename graph_traits< DGtal::DigitalSurface< TDigitalSurfaceContainer > >::edge_iterator, typename graph_traits< DGtal::DigitalSurface< TDigitalSurfaceContainer > >::edge_iterator > edges(const DGtal::DigitalSurface< TDigitalSurfaceContainer > &digSurf)
std::pair< typename graph_traits< DGtal::DigitalSurface< TDigitalSurfaceContainer > >::vertex_iterator, typename graph_traits< DGtal::DigitalSurface< TDigitalSurfaceContainer > >::vertex_iterator > vertices(const DGtal::DigitalSurface< TDigitalSurfaceContainer > &digSurf)
Aim: represents a d-dimensional complex in a d-dimensional space with the following properties and re...
std::vector< Point > cellVertexPositions(const Cell c) const
Cell faceCell(const Face f) const
const VertexRange & cellVertices(const Cell c) const
Face opposite(const Face f) const
Point position(const Vertex v) const
VertexRange faceVertices(const Face f) const
std::vector< Point > faceVertexPositions(const Face f) const
void requireFaceGeometry()
Forces the computation of face geometry.
bool isInfinite(const Cell c) const
Aim: a geometric kernel to compute the convex hull of digital points with integer-only arithmetic.
Aim: Provides a set of functions to facilitate the computation of convex hulls and polytopes,...
static bool computeDelaunayCellComplex(ConvexCellComplex< Point > &cell_complex, const PointRange &input_points, bool remove_duplicates=true)
Aim: Implements the quickhull algorithm by Barber et al. , a famous arbitrary dimensional convex hull...
Definition: QuickHull.h:140
bool setInput(const std::vector< InputPoint > &input_points, bool remove_duplicates=true)
Definition: QuickHull.h:382
std::vector< Index > IndexRange
Definition: QuickHull.h:149
bool getFacetVertices(std::vector< IndexRange > &facet_vertices) const
Definition: QuickHull.h:666
bool getVertexPositions(std::vector< OutputPoint > &vertex_positions)
Definition: QuickHull.h:612
bool computeConvexHull(Status target=Status::VerticesCompleted)
Definition: QuickHull.h:460
Represents a signed cell in a cellular grid space by its Khalimsky coordinates and a boolean value.
Aim: Represents an embedded mesh as faces and a list of vertices. Vertices may be shared among faces ...
Definition: SurfaceMesh.h:92
int main(int argc, char **argv)
Shortcuts< KSpace > SH3
int max(int a, int b)
K init(Point(0, 0, 0), Point(512, 512, 512), true)
KSpace K
HalfEdgeDataStructure::Size Size
ShortcutsGeometry< Z3i::KSpace > SHG3
TriMesh::Face Face
TriMesh::Vertex Vertex