Added initial Unittest for Holefiller and some first cleanupo
This commit is contained in:
@@ -1,13 +1,13 @@
|
||||
/*===========================================================================*\
|
||||
* *
|
||||
* OpenFlipper *
|
||||
* Copyright (c) 2001-2015, RWTH-Aachen University *
|
||||
/* ========================================================================= *
|
||||
* *
|
||||
* OpenMesh *
|
||||
* Copyright (c) 2001-2023, RWTH-Aachen University *
|
||||
* Department of Computer Graphics and Multimedia *
|
||||
* All rights reserved. *
|
||||
* www.openflipper.org *
|
||||
* www.openmesh.org *
|
||||
* *
|
||||
*---------------------------------------------------------------------------*
|
||||
* This file is part of OpenFlipper. *
|
||||
* This file is part of OpenMesh. *
|
||||
*---------------------------------------------------------------------------*
|
||||
* *
|
||||
* Redistribution and use in source and binary forms, with or without *
|
||||
@@ -36,60 +36,29 @@
|
||||
* LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING *
|
||||
* NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS *
|
||||
* SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. *
|
||||
* *
|
||||
\*===========================================================================*/
|
||||
* *
|
||||
* ========================================================================= */
|
||||
|
||||
#pragma once
|
||||
|
||||
#include <vector>
|
||||
#include <float.h>
|
||||
#include <cmath>
|
||||
#include <OpenMesh/Core/Utils/vector_cast.hh>
|
||||
#include "OpenMeshUtils.hh"
|
||||
#include <OpenMesh/Tools/Smoother/JacobiLaplaceSmootherT.hh>
|
||||
#include <OpenMesh/Core/Mesh/PolyConnectivity.hh>
|
||||
|
||||
//=============================================================================
|
||||
|
||||
template< class TheMesh >
|
||||
namespace OpenMesh {
|
||||
namespace HoleFiller {
|
||||
|
||||
template< class MeshT >
|
||||
class HoleFiller
|
||||
{
|
||||
typedef typename MeshT::Point Point;
|
||||
typedef typename MeshT::Scalar Scalar;
|
||||
|
||||
public:
|
||||
typedef TheMesh Mesh;
|
||||
|
||||
import_om_abbreviations( typename Mesh );
|
||||
|
||||
// A weight is a tuple of area and maximum dihedral angle
|
||||
//
|
||||
|
||||
class Weight {
|
||||
public:
|
||||
|
||||
Weight() : angle_( 180 ), area_( FLT_MAX ) {}
|
||||
Weight( Scalar _angle, Scalar _area ) : angle_( _angle ), area_( _area ) {}
|
||||
~Weight() {}
|
||||
|
||||
Scalar angle() const { return angle_; }
|
||||
Scalar area() const { return area_; }
|
||||
|
||||
Weight operator+( const Weight & _other ) const {
|
||||
return Weight( std::max( angle(), _other.angle() ),
|
||||
area() + _other.area() );
|
||||
}
|
||||
|
||||
bool operator<( const Weight & _rhs ) const {
|
||||
return ( angle() < _rhs.angle() ||
|
||||
( angle() == _rhs.angle() && area() < _rhs.area() ) );
|
||||
}
|
||||
|
||||
private:
|
||||
Scalar angle_;
|
||||
Scalar area_;
|
||||
};
|
||||
|
||||
|
||||
// Ctors
|
||||
explicit HoleFiller( Mesh & _mesh );
|
||||
explicit HoleFiller( MeshT & _mesh );
|
||||
~HoleFiller();
|
||||
|
||||
// Identify and fill all holes of the mesh.
|
||||
@@ -97,19 +66,48 @@ public:
|
||||
|
||||
|
||||
// Fill a hole which is identified by one of its boundary edges.
|
||||
void fill_hole( EH _eh, int _stages = 3 );
|
||||
void fill_hole( typename MeshT::EdgeHandle _eh, int _stages = 3 );
|
||||
|
||||
// Fair a filling
|
||||
//void fairing( std::vector< FH >& _faceHandles );
|
||||
void fairing( std::vector< OpenMesh::SmartFaceHandle >& _faceHandles );
|
||||
|
||||
// Remove degenerated faces from the filling
|
||||
void removeDegeneratedFaces( std::vector< FH >& _faceHandles );
|
||||
void removeDegeneratedFaces( std::vector< typename MeshT::FaceHandle >& _faceHandles );
|
||||
|
||||
private:
|
||||
|
||||
// Refine a face
|
||||
bool refine( FH _fh );
|
||||
|
||||
// A weight is a tuple of area and maximum dihedral angle
|
||||
//
|
||||
|
||||
class Weight {
|
||||
public:
|
||||
|
||||
Weight() : angle_( 180 ), area_( FLT_MAX ) {}
|
||||
Weight( Scalar _angle, Scalar _area ) : angle_( _angle ), area_( _area ) {}
|
||||
~Weight() {}
|
||||
|
||||
Scalar angle() const { return angle_; }
|
||||
Scalar area() const { return area_; }
|
||||
|
||||
Weight operator+( const Weight & _other ) const {
|
||||
return Weight( std::max( angle(), _other.angle() ),
|
||||
area() + _other.area() );
|
||||
}
|
||||
|
||||
bool operator<( const Weight & _rhs ) const {
|
||||
return ( angle() < _rhs.angle() ||
|
||||
( angle() == _rhs.angle() && area() < _rhs.area() ) );
|
||||
}
|
||||
|
||||
private:
|
||||
Scalar angle_;
|
||||
Scalar area_;
|
||||
};
|
||||
|
||||
// Refine a face
|
||||
bool refine( typename MeshT::FaceHandle _fh );
|
||||
|
||||
// Relax an edge
|
||||
bool relax_edge( OpenMesh::SmartEdgeHandle _eh );
|
||||
@@ -127,20 +125,20 @@ private:
|
||||
Weight weight( int _i, int _j, int _k );
|
||||
|
||||
// Does edge (_u,_v) already exist?
|
||||
bool exists_edge( OpenMesh::SmartVertexHandle _u, VH _w );
|
||||
bool exists_edge( OpenMesh::SmartVertexHandle _u, typename MeshT::VertexHandle _w );
|
||||
|
||||
// Compute the area of the triangle (_a,_b,_c).
|
||||
Scalar area( VH _a, VH _b, VH _c );
|
||||
Scalar area( typename MeshT::VertexHandle _a, typename MeshT::VertexHandle _b, typename MeshT::VertexHandle _c );
|
||||
|
||||
// Compute the dihedral angle (in degrees) between triangle
|
||||
// (_u,_v,_a) and triangle (_v,_u,_b).
|
||||
Scalar dihedral_angle( VH _u, VH _v, VH _a, VH _b );
|
||||
Scalar dihedral_angle( typename MeshT::VertexHandle _u, typename MeshT::VertexHandle _v, typename MeshT::VertexHandle _a, typename MeshT::VertexHandle _b );
|
||||
|
||||
|
||||
// The mesh, with each vertex we associate a scale factor that is
|
||||
// needed for remeshing
|
||||
|
||||
Mesh & mesh_;
|
||||
MeshT & mesh_;
|
||||
OpenMesh::VPropHandleT< Scalar > scale_;
|
||||
|
||||
/*
|
||||
@@ -157,13 +155,13 @@ private:
|
||||
*/
|
||||
|
||||
|
||||
typedef std::vector< VH > VHVec;
|
||||
typedef typename std::vector< VH >::iterator VHVecIter;
|
||||
typedef typename std::vector< VH >::const_iterator CVHVecIter;
|
||||
typedef std::vector< typename MeshT::VertexHandle > VHVec;
|
||||
typedef typename std::vector< typename MeshT::VertexHandle >::iterator VHVecIter;
|
||||
typedef typename std::vector< typename MeshT::VertexHandle >::const_iterator CVHVecIter;
|
||||
|
||||
typedef std::vector< FH > FHVec;
|
||||
typedef typename std::vector< FH >::iterator FHVecIter;
|
||||
typedef typename std::vector< FH >::const_iterator CFHVecIter;
|
||||
typedef std::vector< typename MeshT::FaceHandle > FHVec;
|
||||
typedef typename std::vector< typename MeshT::FaceHandle >::iterator FHVecIter;
|
||||
typedef typename std::vector< typename MeshT::FaceHandle >::const_iterator CFHVecIter;
|
||||
|
||||
|
||||
// This vector contains all vertices of the hole (in order)
|
||||
@@ -192,9 +190,13 @@ private:
|
||||
std::vector< std::vector< int > > l_;
|
||||
};
|
||||
|
||||
} // namespace HoleFiller
|
||||
} // namespace OpenMesh
|
||||
|
||||
//=============================================================================
|
||||
#ifndef HOLEFILLER_CC
|
||||
#include "HoleFillerT_impl.hh"
|
||||
#endif
|
||||
//=============================================================================
|
||||
|
||||
|
||||
|
||||
@@ -1,13 +1,13 @@
|
||||
/*===========================================================================*\
|
||||
* *
|
||||
* OpenFlipper *
|
||||
* Copyright (c) 2001-2015, RWTH-Aachen University *
|
||||
/* ========================================================================= *
|
||||
* *
|
||||
* OpenMesh *
|
||||
* Copyright (c) 2001-2023, RWTH-Aachen University *
|
||||
* Department of Computer Graphics and Multimedia *
|
||||
* All rights reserved. *
|
||||
* www.openflipper.org *
|
||||
* www.openmesh.org *
|
||||
* *
|
||||
*---------------------------------------------------------------------------*
|
||||
* This file is part of OpenFlipper. *
|
||||
* This file is part of OpenMesh. *
|
||||
*---------------------------------------------------------------------------*
|
||||
* *
|
||||
* Redistribution and use in source and binary forms, with or without *
|
||||
@@ -36,17 +36,24 @@
|
||||
* LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING *
|
||||
* NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS *
|
||||
* SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. *
|
||||
* *
|
||||
\*===========================================================================*/
|
||||
* *
|
||||
* ========================================================================= */
|
||||
|
||||
|
||||
|
||||
//=============================================================================
|
||||
#include "HoleFillerT.hh"
|
||||
#include <OpenMesh/Tools/Smoother/JacobiLaplaceSmootherT.hh>
|
||||
//=============================================================================
|
||||
|
||||
//== NAMESPACES ===============================================================
|
||||
|
||||
|
||||
namespace OpenMesh {
|
||||
namespace HoleFiller {
|
||||
|
||||
template< class MeshT >
|
||||
HoleFiller< MeshT >::HoleFiller( Mesh & _mesh )
|
||||
HoleFiller< MeshT >::HoleFiller(MeshT &_mesh )
|
||||
: mesh_( _mesh )
|
||||
{
|
||||
mesh_.request_vertex_status();
|
||||
@@ -85,7 +92,7 @@ HoleFiller< MeshT >::fill_all_holes( int _stages )
|
||||
{
|
||||
// Collect all boundary edges
|
||||
|
||||
std::vector< EH > bdry_edge;
|
||||
std::vector< typename MeshT::EdgeHandle > bdry_edge;
|
||||
|
||||
for (auto ei : mesh_.edges())
|
||||
if ( ei.is_boundary() )
|
||||
@@ -110,10 +117,10 @@ HoleFiller< MeshT >::fill_all_holes( int _stages )
|
||||
|
||||
std::cerr << "Stage 3 : Smoothing the hole fillings ... ";
|
||||
|
||||
OpenMesh::Smoother::JacobiLaplaceSmootherT< Mesh > smoother( mesh_ );
|
||||
smoother.initialize( OpenMesh::Smoother::SmootherT< Mesh >::
|
||||
OpenMesh::Smoother::JacobiLaplaceSmootherT< MeshT > smoother( mesh_ );
|
||||
smoother.initialize( OpenMesh::Smoother::SmootherT< MeshT >::
|
||||
Tangential_and_Normal,
|
||||
OpenMesh::Smoother::SmootherT< Mesh >::C1 );
|
||||
OpenMesh::Smoother::SmootherT< MeshT >::C1 );
|
||||
|
||||
smoother.smooth( 500 );
|
||||
std::cerr << "ok\n";
|
||||
@@ -129,7 +136,7 @@ HoleFiller< MeshT >::fill_all_holes( int _stages )
|
||||
|
||||
template< class MeshT >
|
||||
void
|
||||
HoleFiller< MeshT >::fill_hole( EH _eh, int _stages )
|
||||
HoleFiller< MeshT >::fill_hole(typename MeshT::EdgeHandle _eh, int _stages )
|
||||
{
|
||||
std::cerr << " Stage 1 : Computing a minimal triangulation ... ";
|
||||
|
||||
@@ -331,7 +338,7 @@ HoleFiller< MeshT >::fairing( std::vector< OpenMesh::SmartFaceHandle >& _faceHan
|
||||
cnt += 1.0f;
|
||||
Point p0 = mesh_.point( vIt );
|
||||
Point p1 = mesh_.point( voh_it.to() );
|
||||
scale += ( p1 - p0 ).norm();
|
||||
scale += norm( p1 - p0 );
|
||||
|
||||
}
|
||||
|
||||
@@ -378,12 +385,12 @@ HoleFiller< MeshT >::fairing( std::vector< OpenMesh::SmartFaceHandle >& _faceHan
|
||||
|
||||
template< class MeshT >
|
||||
bool
|
||||
HoleFiller< MeshT >::refine( FH _fh )
|
||||
HoleFiller< MeshT >::refine(typename MeshT::FaceHandle _fh )
|
||||
{
|
||||
|
||||
// Collect the three edges of the face into e0, e1, e2
|
||||
|
||||
FEI fei = mesh_.fe_iter( _fh );
|
||||
typename MeshT::FEIter fei = mesh_.fe_iter( _fh );
|
||||
OpenMesh::SmartEdgeHandle e0 = *fei; ++fei;
|
||||
OpenMesh::SmartEdgeHandle e1 = *fei; ++fei;
|
||||
OpenMesh::SmartEdgeHandle e2 = *fei; ++fei;
|
||||
@@ -391,11 +398,11 @@ HoleFiller< MeshT >::refine( FH _fh )
|
||||
|
||||
// Collect the vertices, vertex positions and scale factors of the face.
|
||||
|
||||
FVI fvi = mesh_.fv_iter( _fh );
|
||||
typename MeshT::FVIter fvi = mesh_.fv_iter( _fh );
|
||||
|
||||
VH v0 = *fvi; ++fvi;
|
||||
VH v1 = *fvi; ++fvi;
|
||||
VH v2 = *fvi; ++fvi;
|
||||
typename MeshT::VertexHandle v0 = *fvi; ++fvi;
|
||||
typename MeshT::VertexHandle v1 = *fvi; ++fvi;
|
||||
typename MeshT::VertexHandle v2 = *fvi; ++fvi;
|
||||
|
||||
Point p0 = mesh_.point( v0 );
|
||||
Point p1 = mesh_.point( v1 );
|
||||
@@ -410,9 +417,9 @@ HoleFiller< MeshT >::refine( FH _fh )
|
||||
Scalar scale = ( scale0 + scale1 + scale2 ) / 3.0f;
|
||||
Point center = typename MeshT::Scalar(1.0/3.0) * ( p0 + p1 + p2 );
|
||||
|
||||
Scalar d0 = 1.0f * ( p0 - center ).norm();
|
||||
Scalar d1 = 1.0f * ( p1 - center ).norm();
|
||||
Scalar d2 = 1.0f * ( p2 - center ).norm();
|
||||
Scalar d0 = 1.0f * norm( p0 - center );
|
||||
Scalar d1 = 1.0f * norm( p1 - center );
|
||||
Scalar d2 = 1.0f * norm( p2 - center );
|
||||
|
||||
|
||||
//dont split triangles which tend to degenerate
|
||||
@@ -520,15 +527,15 @@ HoleFiller< MeshT >::in_circumsphere( const Point & _x,
|
||||
Point ab = _b - _a;
|
||||
Point ac = _c - _a;
|
||||
|
||||
Scalar a00 = -2.0f * ( ab | _a );
|
||||
Scalar a01 = -2.0f * ( ab | _b );
|
||||
Scalar a02 = -2.0f * ( ab | _c );
|
||||
Scalar b0 = _a.sqrnorm() - _b.sqrnorm();
|
||||
Scalar a00 = -2.0f * ( dot(ab , _a ) );
|
||||
Scalar a01 = -2.0f * ( dot(ab , _b ) );
|
||||
Scalar a02 = -2.0f * ( dot(ab , _c ) );
|
||||
Scalar b0 = norm(_a)*norm(_a) - norm(_b)*norm(_b);
|
||||
|
||||
Scalar a10 = -2.0f * ( ac | _a );
|
||||
Scalar a11 = -2.0f * ( ac | _b );
|
||||
Scalar a12 = -2.0f * ( ac | _c );
|
||||
Scalar b1 = _a.sqrnorm() - _c.sqrnorm();
|
||||
Scalar a10 = -2.0f * ( dot(ac , _a ) );
|
||||
Scalar a11 = -2.0f * ( dot(ac , _b ) );
|
||||
Scalar a12 = -2.0f * ( dot(ac , _c ) );
|
||||
Scalar b1 = norm(_a)*norm(_a) - norm(_c)*norm(_c);
|
||||
|
||||
typename MeshT::Scalar alpha = -(-a11*a02+a01*a12-a12*b0+b1*a02+a11*b0-a01*b1)
|
||||
/ (-a11*a00+a11*a02-a10*a02+a00*a12+a01*a10-a01*a12);
|
||||
@@ -539,7 +546,7 @@ HoleFiller< MeshT >::in_circumsphere( const Point & _x,
|
||||
|
||||
Point center = alpha * _a + beta * _b + gamma * _c;
|
||||
|
||||
return ( _x - center ).sqrnorm() < ( _a - center ).sqrnorm();
|
||||
return norm( _x - center ) * norm( _x - center ) < norm( _a - center ) * norm( _a - center );
|
||||
}
|
||||
|
||||
|
||||
@@ -669,7 +676,7 @@ HoleFiller< MeshT >::weight( int _i, int _j, int _k )
|
||||
|
||||
template< class MeshT >
|
||||
bool
|
||||
HoleFiller< MeshT >::exists_edge( OpenMesh::SmartVertexHandle _u, VH _w )
|
||||
HoleFiller< MeshT >::exists_edge( OpenMesh::SmartVertexHandle _u, typename MeshT::VertexHandle _w )
|
||||
{
|
||||
for ( auto vohi : _u.outgoing_halfedges() )
|
||||
if ( ! vohi.edge().is_boundary() )
|
||||
@@ -690,15 +697,15 @@ HoleFiller< MeshT >::exists_edge( OpenMesh::SmartVertexHandle _u, VH _w )
|
||||
|
||||
template< class MeshT >
|
||||
typename MeshT::Scalar
|
||||
HoleFiller< MeshT >::area( VH _a, VH _b, VH _c )
|
||||
HoleFiller< MeshT >::area( typename MeshT::VertexHandle _a, typename MeshT::VertexHandle _b, typename MeshT::VertexHandle _c )
|
||||
{
|
||||
Point a( mesh_.point( _a ) );
|
||||
Point b( mesh_.point( _b ) );
|
||||
Point c( mesh_.point( _c ) );
|
||||
|
||||
Point n( ( b - a ) % ( c - b ) );
|
||||
Point n( cross(( b - a ) , ( c - b )) );
|
||||
|
||||
return 0.5 * n.norm();
|
||||
return 0.5 * norm(n);
|
||||
}
|
||||
|
||||
|
||||
@@ -715,27 +722,27 @@ HoleFiller< MeshT >::area( VH _a, VH _b, VH _c )
|
||||
|
||||
template< class MeshT >
|
||||
typename MeshT::Scalar
|
||||
HoleFiller< MeshT >::dihedral_angle( VH _u, VH _v, VH _a, VH _b )
|
||||
HoleFiller< MeshT >::dihedral_angle( typename MeshT::VertexHandle _u, typename MeshT::VertexHandle _v, typename MeshT::VertexHandle _a, typename MeshT::VertexHandle _b )
|
||||
{
|
||||
Point u( mesh_.point( _u ) );
|
||||
Point v( mesh_.point( _v ) );
|
||||
Point a( mesh_.point( _a ) );
|
||||
Point b( mesh_.point( _b ) );
|
||||
|
||||
Point n0( cross(( v - u ) , ( a - v )) );
|
||||
Point n1( cross(( u - v ) , ( b - u )) );
|
||||
|
||||
Point n0( ( v - u ) % ( a - v ) );
|
||||
Point n1( ( u - v ) % ( b - u ) );
|
||||
|
||||
n0.normalize();
|
||||
n1.normalize();
|
||||
|
||||
return acos( n0 | n1 ) * 180.0 / M_PI;
|
||||
normalize(n0);
|
||||
normalize(n1);
|
||||
|
||||
return acos( dot(n0,n1) ) * 180.0 / M_PI;
|
||||
}
|
||||
|
||||
|
||||
/// remove degenerated faces
|
||||
template< class MeshT >
|
||||
void
|
||||
HoleFiller< MeshT >::removeDegeneratedFaces( std::vector< FH >& _faceHandles ){
|
||||
HoleFiller< MeshT >::removeDegeneratedFaces( std::vector< typename MeshT::FaceHandle >& _faceHandles ){
|
||||
|
||||
for (int i = _faceHandles.size()-1; i >= 0 ; i--){
|
||||
|
||||
@@ -747,7 +754,7 @@ HoleFiller< MeshT >::removeDegeneratedFaces( std::vector< FH >& _faceHandles ){
|
||||
}
|
||||
|
||||
//get the vertices (works only on triMeshes)
|
||||
FVI fvi = mesh_.fv_iter( _faceHandles[i] );
|
||||
typename MeshT::FaceVertexIterator fvi = mesh_.fv_iter( _faceHandles[i] );
|
||||
Point v0 = mesh_.point( *fvi);
|
||||
++fvi;
|
||||
Point v1 = mesh_.point( *fvi );
|
||||
@@ -762,7 +769,7 @@ HoleFiller< MeshT >::removeDegeneratedFaces( std::vector< FH >& _faceHandles ){
|
||||
|
||||
if (d < FLT_MIN && d > -FLT_MIN) {
|
||||
// degenerated face found
|
||||
FHI hIt = mesh_.fh_iter( _faceHandles[i] );
|
||||
typename MeshT::FaceHalfedgeIterator hIt = mesh_.fh_iter( _faceHandles[i] );
|
||||
|
||||
//try to collapse one of the edges
|
||||
while (hIt.is_valid()){
|
||||
@@ -782,3 +789,8 @@ HoleFiller< MeshT >::removeDegeneratedFaces( std::vector< FH >& _faceHandles ){
|
||||
|
||||
//=============================================================================
|
||||
|
||||
|
||||
//=============================================================================
|
||||
} // namespace HoleFiller
|
||||
} // namespace OpenMesh
|
||||
//=============================================================================
|
||||
|
||||
@@ -185,14 +185,14 @@ compute_weights(LaplaceWeighting _weighting)
|
||||
|
||||
heh2 = Base::mesh_.next_halfedge_handle(heh0);
|
||||
p2 = &Base::mesh_.point(Base::mesh_.to_vertex_handle(heh2));
|
||||
d0 = (*p0 - *p2); d0.normalize();
|
||||
d1 = (*p1 - *p2); d1.normalize();
|
||||
d0 = (*p0 - *p2); normalize(d0);
|
||||
d1 = (*p1 - *p2); normalize(d1);
|
||||
weight += static_cast<typename Mesh::Scalar>(1.0) / tan(acos(std::max(lb, std::min(ub, dot(d0,d1) ))));
|
||||
|
||||
heh2 = Base::mesh_.next_halfedge_handle(heh1);
|
||||
p2 = &Base::mesh_.point(Base::mesh_.to_vertex_handle(heh2));
|
||||
d0 = (*p0 - *p2); d0.normalize();
|
||||
d1 = (*p1 - *p2); d1.normalize();
|
||||
d0 = (*p0 - *p2); normalize(d0);
|
||||
d1 = (*p1 - *p2); normalize(d1);
|
||||
weight += static_cast<typename Mesh::Scalar>(1.0) / tan(acos(std::max(lb, std::min(ub, dot(d0,d1) ))));
|
||||
|
||||
Base::mesh_.property(edge_weights_, *e_it) = weight;
|
||||
|
||||
@@ -11,6 +11,7 @@ set(UNITTEST_SRC
|
||||
unittests_delete_face.cc
|
||||
unittests_eigen3_type.cc
|
||||
unittests_faceless_mesh.cc
|
||||
unittests_holefiller.cc
|
||||
unittests_mc_decimater.cc
|
||||
unittests_mesh_cast.cc
|
||||
unittests_mesh_dual.cc
|
||||
|
||||
BIN
src/Unittests/TestFiles/cube_2holes.off
Normal file
BIN
src/Unittests/TestFiles/cube_2holes.off
Normal file
Binary file not shown.
63
src/Unittests/unittests_holefiller.cc
Normal file
63
src/Unittests/unittests_holefiller.cc
Normal file
@@ -0,0 +1,63 @@
|
||||
|
||||
#include <gtest/gtest.h>
|
||||
#include <Unittests/unittests_common.hh>
|
||||
#include <OpenMesh/Tools/HoleFiller/HoleFillerT.hh>
|
||||
|
||||
namespace {
|
||||
|
||||
class OpenMeshHoleFiller_Triangle : public OpenMeshBase {
|
||||
|
||||
protected:
|
||||
|
||||
// This function is called before each test is run
|
||||
virtual void SetUp() {
|
||||
|
||||
// Do some initial stuff with the member data here...
|
||||
}
|
||||
|
||||
// This function is called after all tests are through
|
||||
virtual void TearDown() {
|
||||
|
||||
// Do some final stuff with the member data here...
|
||||
}
|
||||
|
||||
// Member already defined in OpenMeshBase
|
||||
//Mesh mesh_;
|
||||
};
|
||||
|
||||
/*
|
||||
* ====================================================================
|
||||
* Define tests below
|
||||
* ====================================================================
|
||||
*/
|
||||
|
||||
/*
|
||||
*/
|
||||
TEST_F(OpenMeshHoleFiller_Triangle,Triangle_Hole_Filling) {
|
||||
|
||||
mesh_.clear();
|
||||
|
||||
|
||||
bool ok = OpenMesh::IO::read_mesh(mesh_, "cube_2holes.off");
|
||||
|
||||
ASSERT_TRUE(ok);
|
||||
|
||||
// Check setup
|
||||
EXPECT_EQ(7050u, mesh_.n_vertices() ) << "Wrong number of vertices";
|
||||
EXPECT_EQ(13996u, mesh_.n_faces() ) << "Wrong number of faces";
|
||||
|
||||
|
||||
// Initialize subdivider
|
||||
OpenMesh::HoleFiller::HoleFiller<Mesh> filler(mesh_);
|
||||
|
||||
|
||||
// Execute the algorithm
|
||||
filler.fill_all_holes();
|
||||
|
||||
|
||||
EXPECT_EQ(7526u, mesh_.n_vertices() ) << "Wrong number of vertices after smoothing?";
|
||||
EXPECT_EQ(15048u, mesh_.n_faces() ) << "Wrong number of faces after smoothing?";
|
||||
|
||||
}
|
||||
|
||||
}
|
||||
Reference in New Issue
Block a user