Significant changes to the simplification/optimization framework:

* TriEdgeCollapse is no more multiply derived but it get its "work on two vertex" behavior from a template parameter called VertexPair.
* EdgeCollapse is no more a class devoted to the simplification but it has been renamed as EdgeCollapser and it is a static class templates over a generic <VertexPair> that offer methods to perform Edge Collapses.
* cleaned up the parameter passing method for local optimization classes (and added a baseParameterClass from which every local optimization method must subclass its own parameters)
This commit is contained in:
Paolo Cignoni 2011-05-20 15:12:09 +00:00
parent eec4f43178
commit d7af2e62b3
4 changed files with 253 additions and 414 deletions

View File

@ -20,24 +20,6 @@
* for more details. * * for more details. *
* * * *
****************************************************************************/ ****************************************************************************/
/****************************************************************************
History
$Log: not supported by cvs2svn $
Revision 1.16 2006/10/07 15:04:25 cignoni
removed a useless include
Revision 1.15 2005/10/12 10:36:26 cignoni
Removed unused local type Edge. Now it use the standard simplex edge.
Revision 1.14 2004/12/10 01:04:42 cignoni
better comments
Revision 1.13 2004/11/23 10:34:45 cignoni
passed parameters by reference in many funcs and gcc cleaning
****************************************************************************/
#ifndef __VCG_TETRA_TRI_COLLAPSE #ifndef __VCG_TETRA_TRI_COLLAPSE
#define __VCG_TETRA_TRI_COLLAPSE #define __VCG_TETRA_TRI_COLLAPSE
@ -56,13 +38,13 @@ namespace tri{
See also the corresponding class in the local optimization framework called TriEdgeCollapse See also the corresponding class in the local optimization framework called TriEdgeCollapse
**/ **/
template <class TRI_MESH_TYPE> template <class TRI_MESH_TYPE, class VertexPair>
class EdgeCollapse class EdgeCollapser
{ {
public: public:
/// The tetrahedral mesh type /// The tetrahedral mesh type
typedef TRI_MESH_TYPE TriMeshType; typedef TRI_MESH_TYPE TriMeshType;
/// The tetrahedron type /// The face type
typedef typename TriMeshType::FaceType FaceType; typedef typename TriMeshType::FaceType FaceType;
/// The vertex type /// The vertex type
typedef typename FaceType::VertexType VertexType; typedef typename FaceType::VertexType VertexType;
@ -80,38 +62,30 @@ class EdgeCollapse
///the container of vertex type ///the container of vertex type
typedef typename TriMeshType::VertContainer VertContainer; typedef typename TriMeshType::VertContainer VertContainer;
///half edge type ///half edge type
typedef typename TriMeshType::FaceType::EdgeType EdgeType; //typedef typename TriMeshType::FaceType::EdgeType EdgeType;
/// vector of pos /// vector of pos
typedef typename std::vector<EdgeType> EdgeVec; // typedef typename std::vector<EdgeType> EdgeVec;
///of VFIterator ///of VFIterator
typedef typename vcg::face::VFIterator<FaceType> VFI; typedef typename vcg::face::VFIterator<FaceType> VFI;
/// vector of VFIterator /// vector of VFIterator
typedef typename std::vector<vcg::face::VFIterator<FaceType> > VFIVec; typedef typename std::vector<vcg::face::VFIterator<FaceType> > VFIVec;
private:
struct EdgeSet
/// Default Constructor
EdgeCollapse()
{ {
VFIVec av0,av1,av01;
VFIVec & AV0() { return av0;}
VFIVec & AV1() { return av1;}
VFIVec & AV01(){ return av01;}
}; };
~EdgeCollapse() static void FindSets(VertexPair &p, EdgeSet &es)
{
};
static VFIVec & AV0(){static VFIVec av0; return av0;}
static VFIVec & AV1(){static VFIVec av1; return av1;}
static VFIVec & AV01(){static VFIVec av01; return av01;}
void FindSets(EdgeType &p)
{ {
VertexType * v0 = p.V(0); VertexType * v0 = p.V(0);
VertexType * v1 = p.V(1); VertexType * v1 = p.V(1);
AV0().clear(); // Facce incidenti in v0 es.AV0().clear(); // Facce incidenti in v0
AV1().clear(); // Facce incidenti in v1 es.AV1().clear(); // Facce incidenti in v1
AV01().clear(); // Facce incidenti in v0 e v1 es.AV01().clear(); // Facce incidenti in v0 e v1
VFI x; VFI x;
@ -124,8 +98,8 @@ class EdgeCollapse
zv1 = j; zv1 = j;
break; break;
} }
if(zv1==-1) AV0().push_back( x ); // la faccia x.f non ha il vertice v1 => e' incidente solo in v0 if(zv1==-1) es.AV0().push_back( x ); // la faccia x.f non ha il vertice v1 => e' incidente solo in v0
else AV01().push_back( x ); else es.AV01().push_back( x );
} }
for( x.f = v1->VFp(), x.z = v1->VFi(); x.f!=0; ++x ) for( x.f = v1->VFp(), x.z = v1->VFi(); x.f!=0; ++x )
@ -137,7 +111,7 @@ class EdgeCollapse
zv0 = j; zv0 = j;
break; break;
} }
if(zv0==-1) AV1().push_back( x ); // la faccia x.f non ha il vertice v1 => e' incidente solo in v0 if(zv0==-1) es.AV1().push_back( x ); // la faccia x.f non ha il vertice v1 => e' incidente solo in v0
} }
} }
/* /*
@ -163,7 +137,9 @@ class EdgeCollapse
*/ */
bool LinkConditions(EdgeType pos)
public:
static bool LinkConditions(VertexPair &pos)
{ {
typedef typename vcg::face::VFIterator<FaceType> VFIterator; typedef typename vcg::face::VFIterator<FaceType> VFIterator;
// at the end of the loop each vertex must be counted twice // at the end of the loop each vertex must be counted twice
@ -240,117 +216,25 @@ class EdgeCollapse
return true; return true;
} }
static int Do(TriMeshType &m, VertexPair & c, const Point3<ScalarType> &p)
bool LinkConditionsOld(EdgeType pos){
const int ADJ_1 = TriMeshType::VertexType::NewBitFlag();
const int ADJ_E = TriMeshType::VertexType::NewBitFlag();
//enum {ADJ_1= MeshType::VertexType::USER0,
// ADJ_E= MeshType::VertexType::USER0<<1} ;
// const int ALLADJ = ADJ_1|ADJ_E;
const int NOTALLADJ = ~(ADJ_1 | ADJ_E | TriMeshType::VertexType::VISITED);
const int NOTALLADJ1 = ~(ADJ_E | TriMeshType::VertexType::VISITED);
//EdgePosB<MeshType::face_type::face_base> x;
typename vcg::face::VFIterator<FaceType> x;
// Clear visited and adj flag for all vertices adj to v0;
for(x.f = pos.V(0)->VFp(), x.z = pos.V(0)->VFi(); x.f!=0; ++x ) {
x.f->V1(x.z)->Flags() &= NOTALLADJ;
x.f->V2(x.z)->Flags() &= NOTALLADJ;
}
// Clear visited flag for all vertices adj to v1 and set them adj1 to v1;
for(x.f = pos.V(1)->VFp(), x.z = pos.V(1)->VFi(); x.f!=0; ++x ) {
x.f->V1(x.z)->Flags() &= NOTALLADJ1;
x.f->V2(x.z)->Flags() &= NOTALLADJ1;
}
// Mark vertices adj to v1 as ADJ_1 and adj1 to v1;
for(x.f = pos.V(1)->VFp(), x.z = pos.V(1)->VFi(); x.f!=0; ++x ) {
if(x.f->V1(x.z)==pos.V(0)) x.f->V2(x.z)->Flags() |= ADJ_E | ADJ_1;
else x.f->V2(x.z)->Flags() |= ADJ_1;
if(x.f->V2(x.z)==pos.V(0)) x.f->V1(x.z)->Flags() |= ADJ_E | ADJ_1;
else x.f->V1(x.z)->Flags() |= ADJ_1;
}
// compute the number of:
int adj01=0; // vertices adjacents to both v0 and v1
int adje=0; // vertices adjacents to an egde (usually 2)
for(x.f = pos.V(0)->VFp(), x.z = pos.V(0)->VFi(); x.f!=0; ++x ) {
if(!x.f->V1(x.z)->IsV()) {
x.f->V1(x.z)->SetV();
if(x.f->V1(x.z)->Flags()&ADJ_1) ++adj01;
if(x.f->V1(x.z)->Flags()&ADJ_E) ++adje;
}
if(!x.f->V2(x.z)->IsV()) {
x.f->V2(x.z)->SetV();
if(x.f->V2(x.z)->Flags()&ADJ_1) ++adj01;
if(x.f->V2(x.z)->Flags()&ADJ_E) ++adje;
}
}
//bool val=TopoCheck2();
//if(val != (adj01==adje)) printf("Wrong topo %i %i\n",adj01,adje);
TriMeshType::VertexType::DeleteBitFlag(ADJ_E);
TriMeshType::VertexType::DeleteBitFlag(ADJ_1);
return (adj01==adje);
}
int DoCollapse(TriMeshType &m, EdgeType & c, const Point3<ScalarType> &p)
{ {
FindSets(c); EdgeSet es;
FindSets(c,es);
typename VFIVec::iterator i; typename VFIVec::iterator i;
int n_face_del =0 ; int n_face_del =0 ;
//set Face Face topology for(i=es.AV01().begin();i!=es.AV01().end();++i)
if (TriMeshType::HasFFTopology())
{
//int e0=c.z;
//int e1=c.f->FFi(c.z); //opposite edge
//FaceType *f0=c.f;
//FaceType *f1=f0->FFp(c.z);
//
////take right indexes
//FaceType *f00=f0->FFp((e0+1)%3);
//FaceType *f01=f0->FFp((e0+2)%3);
//int If00=f0->FFi((e0+1)%3);
//int If01=f0->FFi((e0+2)%3);
//
////then attach faces
//f00->FFp(If00)=f01;
//f00->FFi(If00)=If01;
//f01->FFp(If01)=f00;
//f01->FFi(If01)=If00;
////and the ones of face f1
//f00=f1->FFp((e1+1)%3);
//f01=f1->FFp((e1+2)%3);
//If00=f1->FFi((e1+1)%3);
//If01=f1->FFi((e1+2)%3);
//
////and attach faces
//f00->FFp(If00)=f01;
//f00->FFi(If00)=If01;
//f01->FFp(If01)=f00;
//f01->FFi(If01)=If00;
}
for(i=AV01().begin();i!=AV01().end();++i)
{ {
FaceType & f = *((*i).f); FaceType & f = *((*i).f);
assert(f.V((*i).z) == c.V(0)); assert(f.V((*i).z) == c.V(0));
vcg::face::VFDetach(f,((*i).z+1)%3); vcg::face::VFDetach(f,((*i).z+1)%3);
vcg::face::VFDetach(f,((*i).z+2)%3); vcg::face::VFDetach(f,((*i).z+2)%3);
Allocator<TriMeshType>::DeleteFace(m,f); Allocator<TriMeshType>::DeleteFace(m,f);
//n_face_del++; n_face_del++;
} }
//set Vertex Face topology //set Vertex Face topology
for(i=AV0().begin();i!=AV0().end();++i) for(i=es.AV0().begin();i!=es.AV0().end();++i)
{ {
(*i).f->V((*i).z) = c.V(1); // In tutte le facce incidenti in v0, si sostituisce v0 con v1 (*i).f->V((*i).z) = c.V(1); // In tutte le facce incidenti in v0, si sostituisce v0 con v1
(*i).f->VFp((*i).z) = (*i).f->V((*i).z)->VFp(); // e appendo la lista di facce incidenti in v1 a questa faccia (*i).f->VFp((*i).z) = (*i).f->V((*i).z)->VFp(); // e appendo la lista di facce incidenti in v1 a questa faccia
@ -360,7 +244,6 @@ class EdgeCollapse
} }
Allocator<TriMeshType>::DeleteVertex(m,*(c.V(0))); Allocator<TriMeshType>::DeleteVertex(m,*(c.V(0)));
//c.V(0)->SetD();
c.V(1)->P()=p; c.V(1)->P()=p;
return n_face_del; return n_face_del;
} }

View File

@ -98,6 +98,9 @@
#include<vcg/complex/complex.h> #include<vcg/complex/complex.h>
namespace vcg{ namespace vcg{
// Base class for Parameters
// all parameters must be derived from this.
class BaseParameterClass { };
template<class MeshType> template<class MeshType>
class LocalOptimization; class LocalOptimization;
@ -123,33 +126,33 @@ class LocalModification
virtual ModifierType IsOfType() = 0 ; virtual ModifierType IsOfType() = 0 ;
/// return true if the data have not changed since it was created /// return true if the data have not changed since it was created
virtual bool IsUpToDate() = 0 ; virtual bool IsUpToDate() const = 0 ;
/// return true if no constraint disallow this operation to be performed (ex: change of topology in edge collapses) /// return true if no constraint disallow this operation to be performed (ex: change of topology in edge collapses)
virtual bool IsFeasible() = 0; virtual bool IsFeasible(const BaseParameterClass *pp) = 0;
/// Compute the priority to be used in the heap /// Compute the priority to be used in the heap
virtual ScalarType ComputePriority()=0; virtual ScalarType ComputePriority(BaseParameterClass *pp)=0;
/// Return the priority to be used in the heap (implement static priority) /// Return the priority to be used in the heap (implement static priority)
virtual ScalarType Priority() const =0; virtual ScalarType Priority() const =0;
/// Perform the operation and return the variation in the number of simplicies (>0 is refinement, <0 is simplification) /// Perform the operation
virtual void Execute(MeshType &m)=0; virtual void Execute(MeshType &m, BaseParameterClass *pp)=0;
/// perform initialization /// perform initialization
static void Init(MeshType &m, HeapType&); static void Init(MeshType &m, HeapType&, BaseParameterClass *pp);
/// An approximation of the size of the heap with respect of the number of simplex /// An approximation of the size of the heap with respect of the number of simplex
/// of the mesh. When this number is exceeded a clear heap purging is performed. /// of the mesh. When this number is exceeded a clear heap purging is performed.
/// so it is should be reasonably larger than the minimum expected size to avoid too frequent clear heap /// so it is should be reasonably larger than the minimum expected size to avoid too frequent clear heap
/// For example for symmetric edge collapse a 5 is a good guess. /// For example for symmetric edge collapse a 5 is a good guess.
/// while for non symmetric edge collapse a larger number like 9 is a better choice /// while for non symmetric edge collapse a larger number like 9 is a better choice
static float HeapSimplexRatio() {return 6.0f;} ; static float HeapSimplexRatio(BaseParameterClass *) {return 6.0f;} ;
virtual const char *Info(MeshType &) {return 0;} virtual const char *Info(MeshType &) {return 0;}
/// Update the heap as a consequence of this operation /// Update the heap as a consequence of this operation
virtual void UpdateHeap(HeapType&)=0; virtual void UpdateHeap(HeapType&, BaseParameterClass *pp)=0;
}; //end class local modification }; //end class local modification
@ -164,7 +167,7 @@ template<class MeshType>
class LocalOptimization class LocalOptimization
{ {
public: public:
LocalOptimization(MeshType &mm): m(mm){ ClearTermination();e=0.0;HeapSimplexRatio=5;} LocalOptimization(MeshType &mm, BaseParameterClass *_pp): m(mm){ ClearTermination();e=0.0;HeapSimplexRatio=5; pp=_pp;}
struct HeapElem; struct HeapElem;
// scalar type // scalar type
@ -198,6 +201,7 @@ public:
int start; int start;
ScalarType currMetric; ScalarType currMetric;
ScalarType targetMetric; ScalarType targetMetric;
BaseParameterClass *pp;
// The ratio between Heap size and the number of simplices in the current mesh // The ratio between Heap size and the number of simplices in the current mesh
// When this value is exceeded a ClearHeap Start; // When this value is exceeded a ClearHeap Start;
@ -261,7 +265,7 @@ public:
//return (locModPtr->Priority() < h.locModPtr->Priority()); //return (locModPtr->Priority() < h.locModPtr->Priority());
} }
bool IsUpToDate() bool IsUpToDate() const
{ {
return locModPtr->IsUpToDate(); return locModPtr->IsUpToDate();
} }
@ -295,11 +299,11 @@ public:
{ {
//printf("popped out: %s\n",locMod->Info(m)); //printf("popped out: %s\n",locMod->Info(m));
// check if it is feasible // check if it is feasible
if (locMod->IsFeasible()) if (locMod->IsFeasible(this->pp))
{ {
nPerfmormedOps++; nPerfmormedOps++;
locMod->Execute(m); locMod->Execute(m,this->pp);
locMod->UpdateHeap(h); locMod->UpdateHeap(h,this->pp);
} }
} }
//else printf("popped out unfeasible\n"); //else printf("popped out unfeasible\n");
@ -344,9 +348,9 @@ void ClearHeap()
vcg::tri::InitVertexIMark(m); vcg::tri::InitVertexIMark(m);
// The expected size of heap depends on the type of the local modification we are using.. // The expected size of heap depends on the type of the local modification we are using..
HeapSimplexRatio = LocalModificationType::HeapSimplexRatio(); HeapSimplexRatio = LocalModificationType::HeapSimplexRatio(pp);
LocalModificationType::Init(m,h); LocalModificationType::Init(m,h,pp);
std::make_heap(h.begin(),h.end()); std::make_heap(h.begin(),h.end());
if(!h.empty()) currMetric=h.front().pri; if(!h.empty()) currMetric=h.front().pri;
} }
@ -354,7 +358,7 @@ void ClearHeap()
template <class LocalModificationType> void Finalize() template <class LocalModificationType> void Finalize()
{ {
LocalModificationType::Finalize(m,h); LocalModificationType::Finalize(m,h,pp);
} }

View File

@ -20,46 +20,7 @@
* for more details. * * for more details. *
* * * *
****************************************************************************/ ****************************************************************************/
/****************************************************************************
$Log: not supported by cvs2svn $
Revision 1.19 2006/10/15 07:31:21 cignoni
typenames and qualifiers for gcc compliance
Revision 1.18 2006/10/09 20:09:40 cignoni
Changed some access to VertexFaceIterator to reflect the shorter new operators.
Revision 1.17 2005/10/12 10:44:01 cignoni
Now creation of new edge use Ordered() constructor to comply the fact that the basic collapse is simmetric.
Revision 1.16 2005/01/19 10:35:28 cignoni
Better management of symmetric/asymmetric edge collapses
Revision 1.15 2004/12/10 01:03:53 cignoni
better comments and removed logging
Revision 1.14 2004/11/23 10:34:23 cignoni
passed parameters by reference in many funcs and gcc cleaning
Revision 1.13 2004/10/25 16:28:32 ganovelli
pos to edge
Revision 1.12 2004/09/15 11:16:02 ganovelli
changed P() to cP()
Revision 1.11 2004/09/09 13:23:01 ponchio
Header guards typo
Revision 1.10 2004/09/09 13:01:12 ponchio
Linux compatible path in #include
Revision 1.9 2004/09/08 15:13:29 ganovelli
changes for gc
Revision 1.8 2004/09/08 14:33:31 ganovelli
*** empty log message ***
****************************************************************************/
#ifndef __VCG_DECIMATION_TRICOLLAPSE #ifndef __VCG_DECIMATION_TRICOLLAPSE
#define __VCG_DECIMATION_TRICOLLAPSE #define __VCG_DECIMATION_TRICOLLAPSE
@ -81,8 +42,8 @@ namespace tri{
/// This is the base class of all the specialized collapse classes like for example Quadric Edge Collapse. /// This is the base class of all the specialized collapse classes like for example Quadric Edge Collapse.
/// Each derived class /// Each derived class
template<class TriMeshType, class MYTYPE> template<class TriMeshType, class VertexPair, class MYTYPE>
class TriEdgeCollapse: public LocalOptimization<TriMeshType>::LocModType , public EdgeCollapse<TriMeshType> class TriEdgeCollapse: public LocalOptimization<TriMeshType>::LocModType
{ {
public: public:
/// static data to gather statistical information about the reasons of collapse failures /// static data to gather statistical information about the reasons of collapse failures
@ -107,7 +68,7 @@ public:
protected: protected:
typedef typename TriMeshType::FaceType FaceType; typedef typename TriMeshType::FaceType FaceType;
typedef typename TriMeshType::FaceType::VertexType VertexType; typedef typename TriMeshType::FaceType::VertexType VertexType;
typedef typename VertexType::EdgeType EdgeType; // typedef typename VertexType::EdgeType EdgeType;
typedef typename FaceType::VertexType::CoordType CoordType; typedef typename FaceType::VertexType::CoordType CoordType;
typedef typename TriMeshType::VertexType::ScalarType ScalarType; typedef typename TriMeshType::VertexType::ScalarType ScalarType;
typedef typename LocalOptimization<TriMeshType>::HeapElem HeapElem; typedef typename LocalOptimization<TriMeshType>::HeapElem HeapElem;
@ -115,7 +76,7 @@ protected:
TriMeshType *mt; TriMeshType *mt;
///the pair to collapse ///the pair to collapse
EdgeType pos; VertexPair pos;
///mark for up_dating ///mark for up_dating
static int& GlobalMark(){ static int im=0; return im;} static int& GlobalMark(){ static int im=0; return im;}
@ -131,11 +92,11 @@ protected:
inline TriEdgeCollapse() inline TriEdgeCollapse()
{} {}
///Constructor with postype ///Constructor with postype
inline TriEdgeCollapse(const EdgeType &p, int mark) inline TriEdgeCollapse(const VertexPair &p, int mark, BaseParameterClass *pp)
{ {
localMark = mark; localMark = mark;
pos=p; pos=p;
_priority = ComputePriority(); _priority = ComputePriority(pp);
} }
~TriEdgeCollapse() ~TriEdgeCollapse()
@ -146,8 +107,7 @@ private:
public: public:
inline ScalarType ComputePriority(BaseParameterClass *)
inline ScalarType ComputePriority()
{ {
_priority = Distance(pos.V(0)->cP(),pos.V(1)->cP()); _priority = Distance(pos.V(0)->cP(),pos.V(1)->cP());
return _priority; return _priority;
@ -160,20 +120,19 @@ public:
return buf; return buf;
} }
inline void Execute(TriMeshType &m) inline void Execute(TriMeshType &m, BaseParameterClass *)
{ {
CoordType MidPoint=(pos.V(0)->P()+pos.V(1)->P())/2.0; CoordType MidPoint=(pos.V(0)->P()+pos.V(1)->P())/2.0;
DoCollapse(m, pos, MidPoint); EdgeCollapser<TriMeshType,VertexPair>::Do(m, pos, MidPoint);
} }
static bool IsSymmetric() { return true;} static bool IsSymmetric(BaseParameterClass *) { return true;}
// This function is called after an action to re-add in the heap elements whose priority could have been changed. // This function is called after an action to re-add in the heap elements whose priority could have been changed.
// in the plain case we just put again in the heap all the edges around the vertex resulting from the previous collapse: v[1]. // in the plain case we just put again in the heap all the edges around the vertex resulting from the previous collapse: v[1].
// if the collapse is not symmetric you should add also backward edges (because v0->v1 collapse could be different from v1->v0) // if the collapse is not symmetric you should add also backward edges (because v0->v1 collapse could be different from v1->v0)
inline void UpdateHeap(HeapType & h_ret) inline void UpdateHeap(HeapType & h_ret, BaseParameterClass *pp)
{ {
GlobalMark()++; int nn=0; GlobalMark()++; int nn=0;
VertexType *v[2]; VertexType *v[2];
@ -199,20 +158,20 @@ public:
if( !(vfi.V1()->IsV()) && (vfi.V1()->IsRW())) if( !(vfi.V1()->IsV()) && (vfi.V1()->IsRW()))
{ {
vfi.V1()->SetV(); vfi.V1()->SetV();
h_ret.push_back(HeapElem(new MYTYPE(EdgeType( vfi.V(),vfi.V1() ),GlobalMark()))); h_ret.push_back(HeapElem(new MYTYPE(VertexPair( vfi.V(),vfi.V1() ),GlobalMark(),pp)));
std::push_heap(h_ret.begin(),h_ret.end()); std::push_heap(h_ret.begin(),h_ret.end());
if(! this->IsSymmetric()){ if(! this->IsSymmetric(pp)){
h_ret.push_back(HeapElem(new MYTYPE(EdgeType( vfi.V1(),vfi.V()),GlobalMark()))); h_ret.push_back(HeapElem(new MYTYPE(VertexPair( vfi.V1(),vfi.V()),GlobalMark(),pp)));
std::push_heap(h_ret.begin(),h_ret.end()); std::push_heap(h_ret.begin(),h_ret.end());
} }
} }
if( !(vfi.V2()->IsV()) && (vfi.V2()->IsRW())) if( !(vfi.V2()->IsV()) && (vfi.V2()->IsRW()))
{ {
vfi.V2()->SetV(); vfi.V2()->SetV();
h_ret.push_back(HeapElem(new MYTYPE(EdgeType(vfi.F()->V(vfi.I()),vfi.F()->V2(vfi.I())),GlobalMark()))); h_ret.push_back(HeapElem(new MYTYPE(VertexPair(vfi.F()->V(vfi.I()),vfi.F()->V2(vfi.I())),GlobalMark(),pp)));
std::push_heap(h_ret.begin(),h_ret.end()); std::push_heap(h_ret.begin(),h_ret.end());
if(! this->IsSymmetric()){ if(! this->IsSymmetric(pp)){
h_ret.push_back(HeapElem(new MYTYPE(EdgeType (vfi.F()->V1(vfi.I()),vfi.F()->V(vfi.I())),GlobalMark()))); h_ret.push_back(HeapElem(new MYTYPE(VertexPair (vfi.F()->V1(vfi.I()),vfi.F()->V(vfi.I())),GlobalMark(),pp)));
std::push_heap(h_ret.begin(),h_ret.end()); std::push_heap(h_ret.begin(),h_ret.end());
} }
} }
@ -234,27 +193,14 @@ public:
ModifierType IsOfType(){ return TriEdgeCollapseOp;} ModifierType IsOfType(){ return TriEdgeCollapseOp;}
inline bool IsFeasible(){ inline bool IsFeasible(const BaseParameterClass *){
return LinkConditions(pos); return EdgeCollapser<TriMeshType,VertexPair>::LinkConditions(pos);
} }
inline bool IsUpToDate(){ inline bool IsUpToDate() const
// if(pos.V(1)->IsD()) { {
// ++FailStat::OutOfDate(); VertexType *v0=pos.cV(0);
// return false; VertexType *v1=pos.cV(1);
//}
//
// if(pos.V(1)->IsD()) {
// ++FailStat::OutOfDate();
// return false;
//}
VertexType *v0=pos.V(0);
VertexType *v1=pos.V(1);
//if(! (( (!v0->IsD()) && (!v1->IsD())) &&
// localMark>=v0->IMark() &&
// localMark>=v1->IMark()))
if( v0->IsD() || v1->IsD() || if( v0->IsD() || v1->IsD() ||
localMark < v0->IMark() || localMark < v0->IMark() ||
localMark < v1->IMark() ) localMark < v1->IMark() )
@ -269,15 +215,18 @@ public:
return _priority; return _priority;
} }
static void Init(TriMeshType&m,HeapType&h_ret){ static void Init(TriMeshType &m, HeapType &h_ret, BaseParameterClass *pp)
{
vcg::tri::UpdateTopology<TriMeshType>::VertexFace(m);
h_ret.clear(); h_ret.clear();
typename TriMeshType::FaceIterator fi; typename TriMeshType::FaceIterator fi;
for(fi = m.face.begin(); fi != m.face.end();++fi) for(fi = m.face.begin(); fi != m.face.end();++fi)
if(!(*fi).IsD()){ if(!(*fi).IsD()){
for (int j=0;j<3;j++) for (int j=0;j<3;j++)
{ {
EdgeType p=EdgeType::OrderedEdge((*fi).V(j),(*fi).V((j+1)%3)); VertexPair p((*fi).V0(j), (*fi).V1(j));
h_ret.push_back(HeapElem(new MYTYPE(p, IMark(m)))); p.Sort();
h_ret.push_back(HeapElem(new MYTYPE(p, IMark(m),pp)));
//printf("Inserting in heap coll %3i ->%3i %f\n",p.V()-&m.vert[0],p.VFlip()-&m.vert[0],h_ret.back().locModPtr->Priority()); //printf("Inserting in heap coll %3i ->%3i %f\n",p.V()-&m.vert[0],p.VFlip()-&m.vert[0],h_ret.back().locModPtr->Priority());
} }
} }

View File

@ -133,7 +133,7 @@ namespace tri{
}; };
class TriEdgeCollapseQuadricParameter class TriEdgeCollapseQuadricParameter : public BaseParameterClass
{ {
public: public:
double QualityThr; // all double QualityThr; // all
@ -154,20 +154,39 @@ public:
bool QualityQuadric; // During the initialization manage all the edges as border edges adding a set of additional quadrics that are useful mostly for keeping face aspect ratio good. bool QualityQuadric; // During the initialization manage all the edges as border edges adding a set of additional quadrics that are useful mostly for keeping face aspect ratio good.
bool PreserveTopology; bool PreserveTopology;
bool PreserveBoundary; bool PreserveBoundary;
bool MarkComplex;
bool FastPreserveBoundary; bool FastPreserveBoundary;
bool SafeHeapUpdate; bool SafeHeapUpdate;
void SetDefaultParams()
{
UseArea=true;
UseVertexWeight=false;
NormalCheck=false;
NormalThrRad=M_PI/2;
QualityCheck=true;
QualityThr=.1;
BoundaryWeight=.5;
QualityQuadric=false;
OptimalPlacement=true;
ScaleIndependent=true;
QualityWeight=false;
QuadricEpsilon =1e-15;
ScaleFactor=1.0;
PreserveTopology = false;
}
TriEdgeCollapseQuadricParameter() {SetDefaultParams();}
}; };
template<class TriMeshType,class MYTYPE, class HelperType = QInfoStandard<typename TriMeshType::VertexType> > template<class TriMeshType, class VertexPair, class MYTYPE, class HelperType = QInfoStandard<typename TriMeshType::VertexType> >
class TriEdgeCollapseQuadric: public TriEdgeCollapse< TriMeshType, MYTYPE> class TriEdgeCollapseQuadric: public TriEdgeCollapse< TriMeshType, VertexPair, MYTYPE>
{ {
public: public:
typedef typename vcg::tri::TriEdgeCollapse< TriMeshType, MYTYPE > TEC; typedef typename vcg::tri::TriEdgeCollapse< TriMeshType, VertexPair, MYTYPE > TEC;
typedef typename TEC::EdgeType EdgeType; // typedef typename TEC::EdgeType EdgeType;
typedef typename TriEdgeCollapse<TriMeshType, MYTYPE>::HeapType HeapType; typedef typename TriEdgeCollapse<TriMeshType, VertexPair, MYTYPE>::HeapType HeapType;
typedef typename TriEdgeCollapse<TriMeshType, MYTYPE>::HeapElem HeapElem; typedef typename TriEdgeCollapse<TriMeshType, VertexPair, MYTYPE>::HeapElem HeapElem;
typedef typename TriMeshType::CoordType CoordType; typedef typename TriMeshType::CoordType CoordType;
typedef typename TriMeshType::ScalarType ScalarType; typedef typename TriMeshType::ScalarType ScalarType;
typedef math::Quadric< double > QuadricType; typedef math::Quadric< double > QuadricType;
@ -176,10 +195,11 @@ public:
typedef TriEdgeCollapseQuadricParameter QParameter; typedef TriEdgeCollapseQuadricParameter QParameter;
typedef HelperType QH; typedef HelperType QH;
static QParameter & Params(){ // static QParameter & Params(){
static QParameter p; // static QParameter p;
return p; // return p;
} // }
enum Hint { enum Hint {
HNHasFFTopology = 0x0001, // La mesh arriva con la topologia ff gia'fatta HNHasFFTopology = 0x0001, // La mesh arriva con la topologia ff gia'fatta
HNHasVFTopology = 0x0002, // La mesh arriva con la topologia bf gia'fatta HNHasVFTopology = 0x0002, // La mesh arriva con la topologia bf gia'fatta
@ -197,51 +217,51 @@ public:
static std::vector<typename TriMeshType::VertexPointer> _WV; return _WV; static std::vector<typename TriMeshType::VertexPointer> _WV; return _WV;
}; };
inline TriEdgeCollapseQuadric(const EdgeType &p, int i) inline TriEdgeCollapseQuadric(const VertexPair &p, int i, BaseParameterClass *pp)
//:TEC(p,i){}
{ {
this->localMark = i; this->localMark = i;
this->pos=p; this->pos=p;
this->_priority = ComputePriority(); this->_priority = ComputePriority(pp);
} }
inline bool IsFeasible(){ inline bool IsFeasible(BaseParameterClass *_pp){
bool res = ( !Params().PreserveTopology || LinkConditions(this->pos) ); QParameter *pp=(QParameter *)_pp;
if(!res) ++( TriEdgeCollapse< TriMeshType,MYTYPE>::FailStat::LinkConditionEdge() ); bool res = ( !pp->PreserveTopology || EdgeCollapser<TriMeshType, VertexPair>::LinkConditions(this->pos) );
if(!res) ++( TEC::FailStat::LinkConditionEdge() );
return res; return res;
} }
void Execute(TriMeshType &m) void Execute(TriMeshType &m, BaseParameterClass *_pp)
{ CoordType newPos; {
if(Params().OptimalPlacement) newPos= static_cast<MYTYPE*>(this)->ComputeMinimal(); QParameter *pp=(QParameter *)_pp;
CoordType newPos;
if(pp->OptimalPlacement) newPos= static_cast<MYTYPE*>(this)->ComputeMinimal();
else newPos=this->pos.V(1)->P(); else newPos=this->pos.V(1)->P();
//this->pos.V(1)->Qd()+=this->pos.V(0)->Qd();
QH::Qd(this->pos.V(1))+=QH::Qd(this->pos.V(0)); QH::Qd(this->pos.V(1))+=QH::Qd(this->pos.V(0));
//int FaceDel= EdgeCollapser<TriMeshType,VertexPair>::Do(m, this->pos, newPos); // v0 is deleted and v1 take the new position
DoCollapse(m, this->pos, newPos); // v0 is deleted and v1 take the new position
//m.fn-=FaceDel;
//--m.vn;
} }
// Final Clean up after the end of the simplification process // Final Clean up after the end of the simplification process
static void Finalize(TriMeshType &m, HeapType& /*h_ret*/) static void Finalize(TriMeshType &m, HeapType& /*h_ret*/, BaseParameterClass *_pp)
{ {
QParameter *pp=(QParameter *)_pp;
// if the mesh was prepared with precomputed borderflags // if the mesh was prepared with precomputed borderflags
// correctly set them again. // correctly set them again.
if(IsSetHint(HNHasBorderFlag) ) if(IsSetHint(HNHasBorderFlag) )
vcg::tri::UpdateFlags<TriMeshType>::FaceBorderFromVF(m); vcg::tri::UpdateFlags<TriMeshType>::FaceBorderFromVF(m);
// If we had the boundary preservation we should clean up the writable flags // If we had the boundary preservation we should clean up the writable flags
if(Params().FastPreserveBoundary) if(pp->FastPreserveBoundary)
{ {
typename TriMeshType::VertexIterator vi; typename TriMeshType::VertexIterator vi;
for(vi=m.vert.begin();vi!=m.vert.end();++vi) for(vi=m.vert.begin();vi!=m.vert.end();++vi)
if(!(*vi).IsD()) (*vi).SetW(); if(!(*vi).IsD()) (*vi).SetW();
} }
if(Params().PreserveBoundary) if(pp->PreserveBoundary)
{ {
typename std::vector<typename TriMeshType::VertexPointer>::iterator wvi; typename std::vector<typename TriMeshType::VertexPointer>::iterator wvi;
for(wvi=WV().begin();wvi!=WV().end();++wvi) for(wvi=WV().begin();wvi!=WV().end();++wvi)
@ -249,26 +269,21 @@ public:
} }
} }
static void Init(TriMeshType &m,HeapType&h_ret){ static void Init(TriMeshType &m, HeapType &h_ret, BaseParameterClass *_pp)
{
QParameter *pp=(QParameter *)_pp;
typename TriMeshType::VertexIterator vi; typename TriMeshType::VertexIterator vi;
typename TriMeshType::FaceIterator pf; typename TriMeshType::FaceIterator pf;
EdgeType av0,av1,av01; pp->CosineThr=cos(pp->NormalThrRad);
Params().CosineThr=cos(Params().NormalThrRad);
if(!IsSetHint(HNHasVFTopology) ) vcg::tri::UpdateTopology<TriMeshType>::VertexFace(m); if(!IsSetHint(HNHasVFTopology) ) vcg::tri::UpdateTopology<TriMeshType>::VertexFace(m);
if(Params().MarkComplex) {
vcg::tri::UpdateTopology<TriMeshType>::FaceFace(m);
vcg::tri::UpdateFlags<TriMeshType>::FaceBorderFromFF(m);
vcg::tri::UpdateTopology<TriMeshType>::VertexFace(m);
} // e' un po' piu' lenta ma marca i vertici complex
else
if(!IsSetHint(HNHasBorderFlag) ) if(!IsSetHint(HNHasBorderFlag) )
vcg::tri::UpdateFlags<TriMeshType>::FaceBorderFromVF(m); vcg::tri::UpdateFlags<TriMeshType>::FaceBorderFromVF(m);
if(Params().FastPreserveBoundary) if(pp->FastPreserveBoundary)
{ {
for(pf=m.face.begin();pf!=m.face.end();++pf) for(pf=m.face.begin();pf!=m.face.end();++pf)
if( !(*pf).IsD() && (*pf).IsW() ) if( !(*pf).IsD() && (*pf).IsW() )
@ -280,7 +295,7 @@ public:
} }
} }
if(Params().PreserveBoundary) if(pp->PreserveBoundary)
{ {
WV().clear(); WV().clear();
for(pf=m.face.begin();pf!=m.face.end();++pf) for(pf=m.face.begin();pf!=m.face.end();++pf)
@ -293,10 +308,10 @@ public:
} }
} }
InitQuadric(m); InitQuadric(m,pp);
// Initialize the heap with all the possible collapses // Initialize the heap with all the possible collapses
if(IsSymmetric()) if(IsSymmetric(pp))
{ // if the collapse is symmetric (e.g. u->v == v->u) { // if the collapse is symmetric (e.g. u->v == v->u)
for(vi=m.vert.begin();vi!=m.vert.end();++vi) for(vi=m.vert.begin();vi!=m.vert.end();++vi)
if(!(*vi).IsD() && (*vi).IsRW()) if(!(*vi).IsD() && (*vi).IsRW())
@ -311,11 +326,11 @@ public:
assert(x.F()->V(x.I())==&(*vi)); assert(x.F()->V(x.I())==&(*vi));
if((x.V0()<x.V1()) && x.V1()->IsRW() && !x.V1()->IsV()){ if((x.V0()<x.V1()) && x.V1()->IsRW() && !x.V1()->IsV()){
x.V1()->SetV(); x.V1()->SetV();
h_ret.push_back(HeapElem(new MYTYPE(EdgeType(x.V0(),x.V1()),TriEdgeCollapse< TriMeshType,MYTYPE>::GlobalMark() ))); h_ret.push_back(HeapElem(new MYTYPE(VertexPair(x.V0(),x.V1()),TriEdgeCollapse< TriMeshType,VertexPair,MYTYPE>::GlobalMark(),_pp )));
} }
if((x.V0()<x.V2()) && x.V2()->IsRW()&& !x.V2()->IsV()){ if((x.V0()<x.V2()) && x.V2()->IsRW()&& !x.V2()->IsV()){
x.V2()->SetV(); x.V2()->SetV();
h_ret.push_back(HeapElem(new MYTYPE(EdgeType(x.V0(),x.V2()),TriEdgeCollapse< TriMeshType,MYTYPE>::GlobalMark() ))); h_ret.push_back(HeapElem(new MYTYPE(VertexPair(x.V0(),x.V2()),TriEdgeCollapse< TriMeshType,VertexPair,MYTYPE>::GlobalMark(),_pp )));
} }
} }
} }
@ -331,35 +346,19 @@ public:
{ {
assert(x.F()->V(x.I())==&(*vi)); assert(x.F()->V(x.I())==&(*vi));
if(x.V()->IsRW() && x.V1()->IsRW() && !IsMarked(m,x.F()->V1(x.I()))){ if(x.V()->IsRW() && x.V1()->IsRW() && !IsMarked(m,x.F()->V1(x.I()))){
h_ret.push_back( HeapElem( new MYTYPE( EdgeType (x.V(),x.V1()),TriEdgeCollapse< TriMeshType,MYTYPE>::GlobalMark()))); h_ret.push_back( HeapElem( new MYTYPE( VertexPair (x.V(),x.V1()),TriEdgeCollapse< TriMeshType,VertexPair,MYTYPE>::GlobalMark(),_pp)));
} }
if(x.V()->IsRW() && x.V2()->IsRW() && !IsMarked(m,x.F()->V2(x.I()))){ if(x.V()->IsRW() && x.V2()->IsRW() && !IsMarked(m,x.F()->V2(x.I()))){
h_ret.push_back( HeapElem( new MYTYPE( EdgeType (x.V(),x.V2()),TriEdgeCollapse< TriMeshType,MYTYPE>::GlobalMark()))); h_ret.push_back( HeapElem( new MYTYPE( VertexPair (x.V(),x.V2()),TriEdgeCollapse< TriMeshType,VertexPair,MYTYPE>::GlobalMark(),_pp)));
} }
} }
} }
} }
} }
static float HeapSimplexRatio() {return IsSymmetric()?5.0f:9.0f;} static float HeapSimplexRatio(BaseParameterClass *_pp) {return IsSymmetric(_pp)?5.0f:9.0f;}
static bool IsSymmetric() {return Params().OptimalPlacement;} static bool IsSymmetric(BaseParameterClass *_pp) {return ((QParameter *)_pp)->OptimalPlacement;}
static bool IsVertexStable() {return !Params().OptimalPlacement;} static bool IsVertexStable(BaseParameterClass *_pp) {return !((QParameter *)_pp)->OptimalPlacement;}
static void SetDefaultParams(){
Params().UseArea=true;
Params().UseVertexWeight=false;
Params().NormalCheck=false;
Params().NormalThrRad=M_PI/2;
Params().QualityCheck=true;
Params().QualityThr=.1;
Params().BoundaryWeight=.5;
Params().QualityQuadric=false;
Params().OptimalPlacement=true;
Params().ScaleIndependent=true;
Params().QualityWeight=false;
Params().QuadricEpsilon =1e-15;
Params().ScaleFactor=1.0;
Params().PreserveTopology = false;
}
///* ///*
// Funzione principale di valutazione dell'errore del collasso. // Funzione principale di valutazione dell'errore del collasso.
@ -367,7 +366,9 @@ public:
// //
// Da ottimizzare il ciclo sulle normali (deve sparire on e si deve usare per face normals) // Da ottimizzare il ciclo sulle normali (deve sparire on e si deve usare per face normals)
//*/ //*/
ScalarType ComputePriority() { ScalarType ComputePriority(BaseParameterClass *_pp)
{
QParameter *pp=(QParameter *)_pp;
ScalarType error; ScalarType error;
typename vcg::face::VFIterator<FaceType> x; typename vcg::face::VFIterator<FaceType> x;
std::vector<CoordType> on; // original normals std::vector<CoordType> on; // original normals
@ -375,7 +376,7 @@ public:
v[0] = this->pos.V(0); v[0] = this->pos.V(0);
v[1] = this->pos.V(1); v[1] = this->pos.V(1);
if(Params().NormalCheck){ // Compute maximal normal variation if(pp->NormalCheck){ // Compute maximal normal variation
// store the old normals for non-collapsed face in v0 // store the old normals for non-collapsed face in v0
for(x.F() = v[0]->VFp(), x.I() = v[0]->VFi(); x.F()!=0; ++x ) // for all faces in v0 for(x.F() = v[0]->VFp(), x.I() = v[0]->VFi(); x.F()!=0; ++x ) // for all faces in v0
if(x.F()->V(0)!=v[1] && x.F()->V(1)!=v[1] && x.F()->V(2)!=v[1] ) // skip faces with v1 if(x.F()->V(0)!=v[1] && x.F()->V(1)!=v[1] && x.F()->V(2)!=v[1] ) // skip faces with v1
@ -389,7 +390,7 @@ public:
//// Move the two vertexe into new position (storing the old ones) //// Move the two vertexe into new position (storing the old ones)
CoordType OldPos0=v[0]->P(); CoordType OldPos0=v[0]->P();
CoordType OldPos1=v[1]->P(); CoordType OldPos1=v[1]->P();
if(Params().OptimalPlacement) { v[0]->P() = ComputeMinimal(); v[1]->P()=v[0]->P();} if(pp->OptimalPlacement) { v[0]->P() = ComputeMinimal(); v[1]->P()=v[0]->P();}
else v[0]->P() = v[1]->P(); else v[0]->P() = v[1]->P();
//// Rescan faces and compute quality and difference between normals //// Rescan faces and compute quality and difference between normals
@ -402,12 +403,12 @@ public:
for(x.F() = v[0]->VFp(), x.I() = v[0]->VFi(),i=0; x.F()!=0; ++x ) // for all faces in v0 for(x.F() = v[0]->VFp(), x.I() = v[0]->VFi(),i=0; x.F()!=0; ++x ) // for all faces in v0
if(x.F()->V(0)!=v[1] && x.F()->V(1)!=v[1] && x.F()->V(2)!=v[1] ) // skip faces with v1 if(x.F()->V(0)!=v[1] && x.F()->V(1)!=v[1] && x.F()->V(2)!=v[1] ) // skip faces with v1
{ {
if(Params().NormalCheck){ if(pp->NormalCheck){
nn=NormalizedNormal(*x.F()); nn=NormalizedNormal(*x.F());
ndiff=nn.dot(on[i++]); ndiff=nn.dot(on[i++]);
if(ndiff<MinCos) MinCos=ndiff; if(ndiff<MinCos) MinCos=ndiff;
} }
if(Params().QualityCheck){ if(pp->QualityCheck){
qt= QualityFace(*x.F()); qt= QualityFace(*x.F());
if(qt<MinQual) MinQual=qt; if(qt<MinQual) MinQual=qt;
} }
@ -415,12 +416,12 @@ public:
for(x.F() = v[1]->VFp(), x.I() = v[1]->VFi(),i=0; x.F()!=0; ++x ) // for all faces in v1 for(x.F() = v[1]->VFp(), x.I() = v[1]->VFi(),i=0; x.F()!=0; ++x ) // for all faces in v1
if(x.F()->V(0)!=v[0] && x.F()->V(1)!=v[0] && x.F()->V(2)!=v[0] ) // skip faces with v0 if(x.F()->V(0)!=v[0] && x.F()->V(1)!=v[0] && x.F()->V(2)!=v[0] ) // skip faces with v0
{ {
if(Params().NormalCheck){ if(pp->NormalCheck){
nn=NormalizedNormal(*x.F()); nn=NormalizedNormal(*x.F());
ndiff=nn.dot(on[i++]); ndiff=nn.dot(on[i++]);
if(ndiff<MinCos) MinCos=ndiff; if(ndiff<MinCos) MinCos=ndiff;
} }
if(Params().QualityCheck){ if(pp->QualityCheck){
qt= QualityFace(*x.F()); qt= QualityFace(*x.F());
if(qt<MinQual) MinQual=qt; if(qt<MinQual) MinQual=qt;
} }
@ -429,26 +430,26 @@ public:
QuadricType qq=QH::Qd(v[0]); QuadricType qq=QH::Qd(v[0]);
qq+=QH::Qd(v[1]); qq+=QH::Qd(v[1]);
Point3d tpd=Point3d::Construct(v[1]->P()); Point3d tpd=Point3d::Construct(v[1]->P());
double QuadErr = Params().ScaleFactor*qq.Apply(tpd); double QuadErr = pp->ScaleFactor*qq.Apply(tpd);
// All collapses involving triangles with quality larger than <QualityThr> has no penalty; // All collapses involving triangles with quality larger than <QualityThr> has no penalty;
if(MinQual>Params().QualityThr) MinQual=Params().QualityThr; if(MinQual>pp->QualityThr) MinQual=pp->QualityThr;
if(Params().NormalCheck){ if(pp->NormalCheck){
// All collapses where the normal vary less than <NormalThr> (e.g. more than CosineThr) // All collapses where the normal vary less than <NormalThr> (e.g. more than CosineThr)
// have no penalty // have no penalty
if(MinCos>Params().CosineThr) MinCos=Params().CosineThr; if(MinCos>pp->CosineThr) MinCos=pp->CosineThr;
MinCos=(MinCos+1)/2.0; // Now it is in the range 0..1 with 0 very dangerous! MinCos=(MinCos+1)/2.0; // Now it is in the range 0..1 with 0 very dangerous!
} }
if(QuadErr<Params().QuadricEpsilon) QuadErr=Params().QuadricEpsilon; if(QuadErr<pp->QuadricEpsilon) QuadErr=pp->QuadricEpsilon;
if( Params().UseVertexWeight ) QuadErr *= (QH::W(v[1])+QH::W(v[0]))/2; if( pp->UseVertexWeight ) QuadErr *= (QH::W(v[1])+QH::W(v[0]))/2;
if(!Params().QualityCheck && !Params().NormalCheck) error = (ScalarType)(QuadErr); if(!pp->QualityCheck && !pp->NormalCheck) error = (ScalarType)(QuadErr);
if( Params().QualityCheck && !Params().NormalCheck) error = (ScalarType)(QuadErr / MinQual); if( pp->QualityCheck && !pp->NormalCheck) error = (ScalarType)(QuadErr / MinQual);
if(!Params().QualityCheck && Params().NormalCheck) error = (ScalarType)(QuadErr / MinCos); if(!pp->QualityCheck && pp->NormalCheck) error = (ScalarType)(QuadErr / MinCos);
if( Params().QualityCheck && Params().NormalCheck) error = (ScalarType)(QuadErr / (MinQual*MinCos)); if( pp->QualityCheck && pp->NormalCheck) error = (ScalarType)(QuadErr / (MinQual*MinCos));
//Rrestore old position of v0 and v1 //Rrestore old position of v0 and v1
v[0]->P()=OldPos0; v[0]->P()=OldPos0;
@ -460,8 +461,9 @@ public:
// //
//static double MaxError() {return 1e100;} //static double MaxError() {return 1e100;}
// //
inline void UpdateHeap(HeapType & h_ret) inline void UpdateHeap(HeapType & h_ret,BaseParameterClass *_pp)
{ {
QParameter *pp=(QParameter *)_pp;
this->GlobalMark()++; this->GlobalMark()++;
VertexType *v[2]; VertexType *v[2];
v[0]= this->pos.V(0); v[0]= this->pos.V(0);
@ -486,29 +488,29 @@ public:
if( !(vfi.V1()->IsV()) && vfi.V1()->IsRW()) if( !(vfi.V1()->IsV()) && vfi.V1()->IsRW())
{ {
vfi.V1()->SetV(); vfi.V1()->SetV();
h_ret.push_back(HeapElem(new MYTYPE(EdgeType(vfi.V0(),vfi.V1()), this->GlobalMark()))); h_ret.push_back(HeapElem(new MYTYPE(VertexPair(vfi.V0(),vfi.V1()), this->GlobalMark(),_pp)));
std::push_heap(h_ret.begin(),h_ret.end()); std::push_heap(h_ret.begin(),h_ret.end());
if(!IsSymmetric()){ if(!IsSymmetric(pp)){
h_ret.push_back(HeapElem(new MYTYPE(EdgeType(vfi.V1(),vfi.V0()), this->GlobalMark()))); h_ret.push_back(HeapElem(new MYTYPE(VertexPair(vfi.V1(),vfi.V0()), this->GlobalMark(),_pp)));
std::push_heap(h_ret.begin(),h_ret.end()); std::push_heap(h_ret.begin(),h_ret.end());
} }
} }
if( !(vfi.V2()->IsV()) && vfi.V2()->IsRW()) if( !(vfi.V2()->IsV()) && vfi.V2()->IsRW())
{ {
vfi.V2()->SetV(); vfi.V2()->SetV();
h_ret.push_back(HeapElem(new MYTYPE(EdgeType(vfi.V0(),vfi.V2()),this->GlobalMark()))); h_ret.push_back(HeapElem(new MYTYPE(VertexPair(vfi.V0(),vfi.V2()),this->GlobalMark(),_pp)));
std::push_heap(h_ret.begin(),h_ret.end()); std::push_heap(h_ret.begin(),h_ret.end());
if(!IsSymmetric()){ if(!IsSymmetric(pp)){
h_ret.push_back(HeapElem(new MYTYPE(EdgeType(vfi.V2(),vfi.V0()), this->GlobalMark()))); h_ret.push_back( HeapElem(new MYTYPE(VertexPair(vfi.V2(),vfi.V0()), this->GlobalMark(),_pp) ) );
std::push_heap(h_ret.begin(),h_ret.end()); std::push_heap(h_ret.begin(),h_ret.end());
} }
} }
if(Params().SafeHeapUpdate && vfi.V1()->IsRW() && vfi.V2()->IsRW() ) if(pp->SafeHeapUpdate && vfi.V1()->IsRW() && vfi.V2()->IsRW() )
{ {
h_ret.push_back(HeapElem(new MYTYPE(EdgeType(vfi.V1(),vfi.V2()),this->GlobalMark()))); h_ret.push_back(HeapElem(new MYTYPE(VertexPair(vfi.V1(),vfi.V2()),this->GlobalMark(),_pp)));
std::push_heap(h_ret.begin(),h_ret.end()); std::push_heap(h_ret.begin(),h_ret.end());
if(!IsSymmetric()){ if(!IsSymmetric(pp)){
h_ret.push_back(HeapElem(new MYTYPE(EdgeType(vfi.V2(),vfi.V1()), this->GlobalMark()))); h_ret.push_back(HeapElem(new MYTYPE(VertexPair(vfi.V2(),vfi.V1()), this->GlobalMark(),_pp)));
std::push_heap(h_ret.begin(),h_ret.end()); std::push_heap(h_ret.begin(),h_ret.end());
} }
} }
@ -518,8 +520,9 @@ public:
} }
static void InitQuadric(TriMeshType &m) static void InitQuadric(TriMeshType &m,BaseParameterClass *_pp)
{ {
QParameter *pp=(QParameter *)_pp;
typename TriMeshType::FaceIterator pf; typename TriMeshType::FaceIterator pf;
typename TriMeshType::VertexIterator pv; typename TriMeshType::VertexIterator pv;
int j; int j;
@ -540,7 +543,7 @@ static void InitQuadric(TriMeshType &m)
p.SetDirection( ( (*pf).V(1)->cP() - (*pf).V(0)->cP() ) ^ ( (*pf).V(2)->cP() - (*pf).V(0)->cP() )); p.SetDirection( ( (*pf).V(1)->cP() - (*pf).V(0)->cP() ) ^ ( (*pf).V(2)->cP() - (*pf).V(0)->cP() ));
// Se normalizzo non dipende dall'area // Se normalizzo non dipende dall'area
if(!Params().UseArea) if(!pp->UseArea)
p.Normalize(); p.Normalize();
p.SetOffset( p.Direction().dot((*pf).V(0)->cP())); p.SetOffset( p.Direction().dot((*pf).V(0)->cP()));
@ -551,13 +554,13 @@ static void InitQuadric(TriMeshType &m)
for(j=0;j<3;++j) for(j=0;j<3;++j)
if( (*pf).V(j)->IsW() ) if( (*pf).V(j)->IsW() )
{ {
if(Params().QualityWeight) if(pp->QualityWeight)
q*=(*pf).V(j)->Q(); q*=(*pf).V(j)->Q();
QH::Qd((*pf).V(j)) += q; // Sommo la quadrica ai vertici QH::Qd((*pf).V(j)) += q; // Sommo la quadrica ai vertici
} }
for(j=0;j<3;++j) for(j=0;j<3;++j)
if( (*pf).IsB(j) || Params().QualityQuadric ) // Bordo! if( (*pf).IsB(j) || pp->QualityQuadric ) // Bordo!
{ {
Plane3<ScalarType,false> pb; // Piano di bordo Plane3<ScalarType,false> pb; // Piano di bordo
@ -566,8 +569,8 @@ static void InitQuadric(TriMeshType &m)
// poiche' la pesatura in funzione dell'area e'gia fatta in p.Direction() // poiche' la pesatura in funzione dell'area e'gia fatta in p.Direction()
// Senza la normalize il bordo e' pesato in funzione della grandezza della mesh (mesh grandi non decimano sul bordo) // Senza la normalize il bordo e' pesato in funzione della grandezza della mesh (mesh grandi non decimano sul bordo)
pb.SetDirection(p.Direction() ^ ( (*pf).V1(j)->cP() - (*pf).V(j)->cP() ).normalized()); pb.SetDirection(p.Direction() ^ ( (*pf).V1(j)->cP() - (*pf).V(j)->cP() ).normalized());
if( (*pf).IsB(j) ) pb.SetDirection(pb.Direction()* (ScalarType)Params().BoundaryWeight); // amplify border planes if( (*pf).IsB(j) ) pb.SetDirection(pb.Direction()* (ScalarType)pp->BoundaryWeight); // amplify border planes
else pb.SetDirection(pb.Direction()* (ScalarType)(Params().BoundaryWeight/100.0)); // and consider much less quadric for quality else pb.SetDirection(pb.Direction()* (ScalarType)(pp->BoundaryWeight/100.0)); // and consider much less quadric for quality
pb.SetOffset(pb.Direction().dot((*pf).V(j)->cP())); pb.SetOffset(pb.Direction().dot((*pf).V(j)->cP()));
q.ByPlane(pb); q.ByPlane(pb);
@ -576,14 +579,14 @@ static void InitQuadric(TriMeshType &m)
} }
} }
if(Params().ScaleIndependent) if(pp->ScaleIndependent)
{ {
vcg::tri::UpdateBounding<TriMeshType>::Box(m); vcg::tri::UpdateBounding<TriMeshType>::Box(m);
//Make all quadric independent from mesh size //Make all quadric independent from mesh size
Params().ScaleFactor = 1e8*pow(1.0/m.bbox.Diag(),6); // scaling factor pp->ScaleFactor = 1e8*pow(1.0/m.bbox.Diag(),6); // scaling factor
//Params().ScaleFactor *=Params().ScaleFactor ; //pp->ScaleFactor *=pp->ScaleFactor ;
//Params().ScaleFactor *=Params().ScaleFactor ; //pp->ScaleFactor *=pp->ScaleFactor ;
//printf("Scale factor =%f\n",Params().ScaleFactor ); //printf("Scale factor =%f\n",pp->ScaleFactor );
//printf("bb (%5.2f %5.2f %5.2f)-(%5.2f %5.2f %5.2f) Diag %f\n",m.bbox.min[0],m.bbox.min[1],m.bbox.min[2],m.bbox.max[0],m.bbox.max[1],m.bbox.max[2],m.bbox.Diag()); //printf("bb (%5.2f %5.2f %5.2f)-(%5.2f %5.2f %5.2f) Diag %f\n",m.bbox.min[0],m.bbox.min[1],m.bbox.min[2],m.bbox.max[0],m.bbox.max[1],m.bbox.max[2],m.bbox.Diag());
} }
} }
@ -597,7 +600,7 @@ static void InitQuadric(TriMeshType &m)
// //
// //
//static void InitMesh(MESH_TYPE &m){ //static void InitMesh(MESH_TYPE &m){
// Params().CosineThr=cos(Params().NormalThr); // pp->CosineThr=cos(pp->NormalThr);
// InitQuadric(m); // InitQuadric(m);
// //m.Topology(); // //m.Topology();
// //OldInitQuadric(m,UseArea); // //OldInitQuadric(m,UseArea);