Committed temporary version of the cleaned up curvature computation files

This commit is contained in:
Paolo Cignoni 2012-11-12 11:15:21 +00:00
parent 7a8aba311b
commit ed6042e502
2 changed files with 432 additions and 432 deletions

View File

@ -19,55 +19,20 @@
* GNU General Public License (http://www.gnu.org/licenses/gpl.txt) * * GNU General Public License (http://www.gnu.org/licenses/gpl.txt) *
* for more details. * * for more details. *
* * * *
****************************************************************************/
/****************************************************************************
History
$Log: not supported by cvs2svn $
Revision 1.8 2008/05/14 10:03:29 ganovelli
Point3f->Coordtype
Revision 1.7 2008/04/23 16:37:15 onnis
VertexCurvature method added.
Revision 1.6 2008/04/04 10:26:12 cignoni
Cleaned up names, now Kg() gives back Gaussian Curvature (k1*k2), while Kh() gives back Mean Curvature 1/2(k1+k2)
Revision 1.5 2008/03/25 11:00:56 ganovelli
fixed bugs sign of principal direction and mean curvature value
Revision 1.4 2008/03/17 11:29:59 ganovelli
taubin and desbrun estimates added (-> see vcg/simplex/vertex/component.h [component_ocf.h|component_occ.h ]
Revision 1.3 2006/02/27 18:02:57 ponchio
Area -> doublearea/2
added some typename
Revision 1.2 2005/10/25 09:17:41 spinelli
correct IsBorder
Revision 1.1 2005/02/22 16:40:29 ganovelli
created. This version writes the gaussian curvature on the Q() member of
the vertex
****************************************************************************/ ****************************************************************************/
#ifndef VCGLIB_UPDATE_CURVATURE_ #ifndef VCGLIB_UPDATE_CURVATURE_
#define VCGLIB_UPDATE_CURVATURE_ #define VCGLIB_UPDATE_CURVATURE_
#include <vcg/space/index/grid_static_ptr.h> #include <vcg/space/index/grid_static_ptr.h>
#include <vcg/math/base.h>
#include <vcg/math/matrix.h>
#include <vcg/simplex/face/topology.h> #include <vcg/simplex/face/topology.h>
#include <vcg/simplex/face/pos.h> #include <vcg/simplex/face/pos.h>
#include <vcg/simplex/face/jumping_pos.h> #include <vcg/simplex/face/jumping_pos.h>
#include <vcg/container/simple_temporary_data.h>
#include <vcg/complex/algorithms/update/normal.h> #include <vcg/complex/algorithms/update/normal.h>
#include <vcg/complex/algorithms/point_sampling.h> #include <vcg/complex/algorithms/point_sampling.h>
#include <vcg/complex/append.h>
#include <vcg/complex/algorithms/intersection.h> #include <vcg/complex/algorithms/intersection.h>
#include <vcg/complex/algorithms/inertia.h> #include <vcg/complex/algorithms/inertia.h>
#include <vcg/math/matrix33.h> #include <eigenlib/Eigen/Core>
#include <wrap/callback.h> #include <wrap/callback.h>
namespace vcg { namespace vcg {
@ -100,6 +65,7 @@ public:
private: private:
// aux data struct used by PrincipalDirections
struct AdjVertex { struct AdjVertex {
VertexType * vert; VertexType * vert;
float doubleArea; float doubleArea;
@ -122,14 +88,15 @@ public:
URL = "http://dx.doi.org/10.1109/ICCV.1995.466840", URL = "http://dx.doi.org/10.1109/ICCV.1995.466840",
bibsource = "http://www.visionbib.com/bibliography/describe440.html#TT32253", bibsource = "http://www.visionbib.com/bibliography/describe440.html#TT32253",
*/ */
static void PrincipalDirections(MeshType &m) { static void PrincipalDirections(MeshType &m)
{
assert(tri::HasPerFaceVFAdjacency(m) && tri::HasPerVertexVFAdjacency(m)); assert(tri::HasPerFaceVFAdjacency(m) && tri::HasPerVertexVFAdjacency(m));
vcg::tri::UpdateNormal<MeshType>::PerVertexNormalized(m); vcg::tri::UpdateNormal<MeshType>::PerVertexAngleWeighted(m);
vcg::tri::UpdateNormal<MeshType>::NormalizePerVertex(m);
VertexIterator vi; for (VertexIterator vi =m.vert.begin(); vi !=m.vert.end(); ++vi) {
for (vi =m.vert.begin(); vi !=m.vert.end(); ++vi) {
if ( ! (*vi).IsD() && (*vi).VFp() != NULL) { if ( ! (*vi).IsD() && (*vi).VFp() != NULL) {
VertexType * central_vertex = &(*vi); VertexType * central_vertex = &(*vi);
@ -259,25 +226,25 @@ public:
c = best_c; c = best_c;
s = best_s; s = best_s;
vcg::ndim::MatrixMNf minor2x2 (2,2); Eigen::Matrix2f minor2x2;
vcg::ndim::MatrixMNf S (2,2); Eigen::Matrix2f S;
// diagonalize M // diagonalize M
minor2x2[0][0] = QtMQ[1][1]; minor2x2(0,0) = QtMQ[1][1];
minor2x2[0][1] = QtMQ[1][2]; minor2x2(0,1) = QtMQ[1][2];
minor2x2[1][0] = QtMQ[2][1]; minor2x2(1,0) = QtMQ[2][1];
minor2x2[1][1] = QtMQ[2][2]; minor2x2(1,1) = QtMQ[2][2];
S[0][0] = S[1][1] = c; S(0,0) = S(1,1) = c;
S[0][1] = s; S(0,1) = s;
S[1][0] = -1.0f * s; S(1,0) = -1.0f * s;
vcg::ndim::MatrixMNf StMS(S.transpose() * minor2x2 * S); Eigen::Matrix2f StMS = S.transpose() * minor2x2 * S;
// compute curvatures and curvature directions // compute curvatures and curvature directions
float Principal_Curvature1 = (3.0f * StMS[0][0]) - StMS[1][1]; float Principal_Curvature1 = (3.0f * StMS(0,0)) - StMS(1,1);
float Principal_Curvature2 = (3.0f * StMS[1][1]) - StMS[0][0]; float Principal_Curvature2 = (3.0f * StMS(1,1)) - StMS(0,0);
CoordType Principal_Direction1 = T1 * c - T2 * s; CoordType Principal_Direction1 = T1 * c - T2 * s;
CoordType Principal_Direction2 = T1 * s + T2 * c; CoordType Principal_Direction2 = T1 * s + T2 * c;
@ -300,16 +267,17 @@ public:
}; };
/** Curvature meseaure as described in the paper: /** Curvature meseaure as described in the paper:
Robust principal curvatures on Multiple Scales, Yong-Liang Yang, Yu-Kun Lai, Shi-Min Hu Helmut Pottmann Robust principal curvatures on Multiple Scales, Yong-Liang Yang, Yu-Kun Lai, Shi-Min Hu Helmut Pottmann
SGP 2004 SGP 2004
If pointVSfaceInt==true the covariance is computed by montecarlo sampling on the mesh (faster) If pointVSfaceInt==true the covariance is computed by montecarlo sampling on the mesh (faster)
If pointVSfaceInt==false the covariance is computed by (analytic)integration over the surface (slower) If pointVSfaceInt==false the covariance is computed by (analytic)integration over the surface (slower)
*/ */
typedef vcg::GridStaticPtr <FaceType, ScalarType > MeshGridType; typedef vcg::GridStaticPtr <FaceType, ScalarType > MeshGridType;
typedef vcg::GridStaticPtr <VertexType, ScalarType > PointsGridType; typedef vcg::GridStaticPtr <VertexType, ScalarType > PointsGridType;
static void PrincipalDirectionsPCA(MeshType &m, ScalarType r, bool pointVSfaceInt = true,vcg::CallBackPos * cb = NULL) { static void PrincipalDirectionsPCA(MeshType &m, ScalarType r, bool pointVSfaceInt = true,vcg::CallBackPos * cb = NULL)
{
std::vector<VertexType*> closests; std::vector<VertexType*> closests;
std::vector<ScalarType> distances; std::vector<ScalarType> distances;
std::vector<CoordType> points; std::vector<CoordType> points;
@ -318,23 +286,32 @@ public:
MeshType tmpM; MeshType tmpM;
typename std::vector<CoordType>::iterator ii; typename std::vector<CoordType>::iterator ii;
vcg::tri::TrivialSampler<MeshType> vs; vcg::tri::TrivialSampler<MeshType> vs;
tri::UpdateNormal<MeshType>::PerVertexAngleWeighted(m);
tri::UpdateNormal<MeshType>::NormalizePerVertex(m);
MeshGridType mGrid; MeshGridType mGrid;
PointsGridType pGrid; PointsGridType pGrid;
// Fill the grid used // Fill the grid used
if(pointVSfaceInt){ if(pointVSfaceInt)
{
area = Stat<MeshType>::ComputeMeshArea(m); area = Stat<MeshType>::ComputeMeshArea(m);
vcg::tri::SurfaceSampling<MeshType,vcg::tri::TrivialSampler<MeshType> >::Montecarlo(m,vs,1000 * area / (2*M_PI*r*r )); vcg::tri::SurfaceSampling<MeshType,vcg::tri::TrivialSampler<MeshType> >::Montecarlo(m,vs,1000 * area / (2*M_PI*r*r ));
vi = vcg::tri::Allocator<MeshType>::AddVertices(tmpM,m.vert.size()); vi = vcg::tri::Allocator<MeshType>::AddVertices(tmpM,m.vert.size());
for(size_t y = 0; y < m.vert.size(); ++y,++vi) (*vi).P() = m.vert[y].P(); for(size_t y = 0; y < m.vert.size(); ++y,++vi) (*vi).P() = m.vert[y].P();
pGrid.Set(tmpM.vert.begin(),tmpM.vert.end()); pGrid.Set(tmpM.vert.begin(),tmpM.vert.end());
} else{ mGrid.Set(m.face.begin(),m.face.end()); } }
else
{
mGrid.Set(m.face.begin(),m.face.end());
}
int jj = 0; int jj = 0;
for(vi = m.vert.begin(); vi != m.vert.end(); ++vi){ for(vi = m.vert.begin(); vi != m.vert.end(); ++vi)
vcg::Matrix33<ScalarType> A,eigenvectors; {
vcg::Point3<ScalarType> bp,eigenvalues; vcg::Matrix33<ScalarType> A, eigenvectors;
int nrot; vcg::Point3<ScalarType> bp, eigenvalues;
// int nrot;
// sample the neighborhood // sample the neighborhood
if(pointVSfaceInt) if(pointVSfaceInt)
@ -354,7 +331,22 @@ public:
vcg::tri::Inertia<MeshType>::Covariance(tmpM,_bary,A); vcg::tri::Inertia<MeshType>::Covariance(tmpM,_bary,A);
} }
Jacobi(A, eigenvalues , eigenvectors, nrot); // Eigen::Matrix3f AA; AA=A;
// Eigen::SelfAdjointEigenSolver<Eigen::Matrix3f> eig(AA);
// Eigen::Vector3f c_val = eig.eigenvalues();
// Eigen::Matrix3f c_vec = eig.eigenvectors();
// Jacobi(A, eigenvalues , eigenvectors, nrot);
Eigen::Matrix3d AA;
A.ToEigenMatrix(AA);
Eigen::SelfAdjointEigenSolver<Eigen::Matrix3d> eig(AA);
Eigen::Vector3d c_val = eig.eigenvalues();
Eigen::Matrix3d c_vec = eig.eigenvectors(); // eigenvector are stored as columns.
eigenvectors.FromEigenMatrix(c_vec);
eigenvalues.FromEigenVector(c_val);
// EV.transposeInPlace();
// ev.FromEigenVector(c_val);
// get the estimate of curvatures from eigenvalues and eigenvectors // get the estimate of curvatures from eigenvalues and eigenvectors
// find the 2 most tangent eigenvectors (by finding the one closest to the normal) // find the 2 most tangent eigenvectors (by finding the one closest to the normal)
@ -393,15 +385,21 @@ public:
} }
/// \brief Computes the discrete gaussian curvature.
/** For further details, please, refer to: \n /// \brief Computes the discrete mean gaussian curvature.
/**
The algorithm used is the one Desbrun et al. that is based on a discrete analysis of the angles of the faces around a vertex.
- <em> Discrete Differential-Geometry Operators for Triangulated 2-Manifolds Mark Meyer, It requires FaceFace Adjacency;
Mathieu Desbrun, Peter Schroder, Alan H. Barr VisMath '02, Berlin </em>
For further details, please, refer to: \n
<b>Discrete Differential-Geometry Operators for Triangulated 2-Manifolds </b><br>
<i>Mark Meyer, Mathieu Desbrun, Peter Schroder, Alan H. Barr</i><br>
VisMath '02, Berlin
*/ */
static void MeanAndGaussian(MeshType & m) static void MeanAndGaussian(MeshType & m)
{ {
assert(HasFFAdjacency(m)); assert(HasFFAdjacency(m));
float area0, area1, area2, angle0, angle1, angle2; float area0, area1, area2, angle0, angle1, angle2;
FaceIterator fi; FaceIterator fi;
@ -520,7 +518,7 @@ public:
(*vi).Kg() /= (TDAreaPtr)[*vi].A; (*vi).Kg() /= (TDAreaPtr)[*vi].A;
} }
} }
} }
/// \brief Update the mean and the gaussian curvature of a vertex. /// \brief Update the mean and the gaussian curvature of a vertex.
@ -531,12 +529,8 @@ public:
Based on the paper <a href="http://www2.in.tu-clausthal.de/~hormann/papers/Dyn.2001.OTU.pdf"> <em> "Optimizing 3d triangulations using discrete curvature analysis" </em> </a> Based on the paper <a href="http://www2.in.tu-clausthal.de/~hormann/papers/Dyn.2001.OTU.pdf"> <em> "Optimizing 3d triangulations using discrete curvature analysis" </em> </a>
*/ */
static float VertexCurvature(VertexPointer v, bool norm = true) static float ComputeSingleVertexCurvature(VertexPointer v, bool norm = true)
{ {
// VFAdjacency required!
assert(FaceType::HasVFAdjacency());
assert(VertexType::HasVFAdjacency());
VFIteratorType vfi(v); VFIteratorType vfi(v);
float A = 0; float A = 0;
@ -595,10 +589,14 @@ public:
return A; return A;
} }
static void VertexCurvature(MeshType & m){ static void PerVertex(MeshType & m)
{
// VFAdjacency required!
assert(FaceType::HasVFAdjacency());
assert(VertexType::HasVFAdjacency());
for(VertexIterator vi = m.vert.begin(); vi != m.vert.end(); ++vi) for(VertexIterator vi = m.vert.begin(); vi != m.vert.end(); ++vi)
VertexCurvature(&*vi,false); ComputeSingleVertexCurvature(&*vi,false);
} }
@ -613,12 +611,11 @@ public:
} }
*/ */
static void PrincipalDirectionsNormalCycles(MeshType & m){ static void PrincipalDirectionsNormalCycle(MeshType & m){
assert(VertexType::HasVFAdjacency()); assert(VertexType::HasVFAdjacency());
assert(FaceType::HasFFAdjacency()); assert(FaceType::HasFFAdjacency());
assert(FaceType::HasFaceNormal()); assert(FaceType::HasFaceNormal());
typename MeshType::VertexIterator vi; typename MeshType::VertexIterator vi;
for(vi = m.vert.begin(); vi != m.vert.end(); ++vi) for(vi = m.vert.begin(); vi != m.vert.end(); ++vi)
@ -657,28 +654,36 @@ public:
m33[2][0] = m33[0][2]; m33[2][0] = m33[0][2];
m33[2][1] = m33[1][2]; m33[2][1] = m33[1][2];
Eigen::Matrix3d it;
m33.ToEigenMatrix(it);
Eigen::SelfAdjointEigenSolver<Eigen::Matrix3d> eig(it);
Eigen::Vector3d c_val = eig.eigenvalues();
Eigen::Matrix3d c_vec = eig.eigenvectors();
Point3<ScalarType> lambda; Point3<ScalarType> lambda;
Matrix33<ScalarType> vect; Matrix33<ScalarType> vect;
int n_rot; vect.FromEigenMatrix(c_vec);
Jacobi(m33,lambda, vect,n_rot); lambda.FromEigenVector(c_val);
vect.transposeInPlace(); ScalarType bestNormal = 0;
ScalarType normal = std::numeric_limits<ScalarType>::min(); int bestNormalIndex = -1;
int normI = 0;
for(int i = 0; i < 3; ++i) for(int i = 0; i < 3; ++i)
if( fabs((*vi).N().Normalize().dot(vect.GetRow(i))) > normal )
{ {
normal= fabs((*vi).N().Normalize().dot(vect.GetRow(i))); float agreeWithNormal = fabs((*vi).N().Normalize().dot(vect.GetColumn(i)));
normI = i; if( agreeWithNormal > bestNormal )
{
bestNormal= agreeWithNormal;
bestNormalIndex = i;
} }
int maxI = (normI+2)%3; }
int minI = (normI+1)%3; int maxI = (bestNormalIndex+2)%3;
int minI = (bestNormalIndex+1)%3;
if(fabs(lambda[maxI]) < fabs(lambda[minI])) std::swap(maxI,minI); if(fabs(lambda[maxI]) < fabs(lambda[minI])) std::swap(maxI,minI);
(*vi).PD1() = *(Point3<ScalarType>*)(& vect[maxI][0]); (*vi).PD1() = *(Point3<ScalarType>*)(& vect[maxI][0]);
(*vi).PD2() = *(Point3<ScalarType>*)(& vect[minI][0]); (*vi).PD2() = *(Point3<ScalarType>*)(& vect[minI][0]);
(*vi).K1() = lambda[maxI]; (*vi).K1() = lambda[2];
(*vi).K2() = lambda[minI]; (*vi).K2() = lambda[1];
} }
} }
}; };

View File

@ -28,27 +28,21 @@
#define VCGLIB_UPDATE_CURVATURE_FITTING #define VCGLIB_UPDATE_CURVATURE_FITTING
#include <vcg/space/index/grid_static_ptr.h> #include <vcg/space/index/grid_static_ptr.h>
#include <vcg/math/base.h>
#include <vcg/math/matrix.h>
#include <vcg/simplex/face/topology.h> #include <vcg/simplex/face/topology.h>
#include <vcg/simplex/face/pos.h> #include <vcg/simplex/face/pos.h>
#include <vcg/simplex/face/jumping_pos.h> #include <vcg/simplex/face/jumping_pos.h>
#include <vcg/container/simple_temporary_data.h> #include <vcg/container/simple_temporary_data.h>
#include <vcg/complex/algorithms/update/normal.h> #include <vcg/complex/algorithms/update/normal.h>
#include <vcg/complex/algorithms/point_sampling.h> #include <vcg/complex/algorithms/point_sampling.h>
#include <vcg/complex/append.h>
#include <vcg/complex/algorithms/intersection.h> #include <vcg/complex/algorithms/intersection.h>
#include <vcg/complex/algorithms/inertia.h> #include <vcg/complex/algorithms/inertia.h>
#include <vcg/math/matrix33.h> #include <vcg/complex/algorithms/nring.h>
#include <eigenlib/Eigen/Core> #include <eigenlib/Eigen/Core>
#include <eigenlib/Eigen/QR> #include <eigenlib/Eigen/QR>
#include <eigenlib/Eigen/LU> #include <eigenlib/Eigen/LU>
#include <eigenlib/Eigen/SVD> #include <eigenlib/Eigen/SVD>
#include <eigenlib/Eigen/Eigenvalues> #include <eigenlib/Eigen/Eigenvalues>
// GG include
#include <vector>
#include <vcg/complex/algorithms/nring.h>
namespace vcg { namespace vcg {
@ -116,17 +110,17 @@ class Quadric
return 2.0*c()*v + b()*u + e(); return 2.0*c()*v + b()*u + e();
} }
double duv(double u, double v) double duv(double /*u*/, double /*v*/)
{ {
return b(); return b();
} }
double duu(double u, double v) double duu(double /*u*/, double /*v*/)
{ {
return 2.0*a(); return 2.0*a();
} }
double dvv(double u, double v) double dvv(double /*u*/, double /*v*/)
{ {
return 2.0*c(); return 2.0*c();
} }
@ -231,8 +225,9 @@ class Quadric
static void computeCurvature(MeshType & m) static void computeCurvature(MeshType & m)
{ {
Allocator<MeshType>::CompactVertexVector(m);
assert(tri::HasPerVertexVFAdjacency(m) && tri::HasPerFaceVFAdjacency(m) ); if(!HasFVAdjacency(m)) throw vcg::MissingComponentException("FVAdjacency");
vcg::tri::UpdateTopology<MeshType>::VertexFace(m); vcg::tri::UpdateTopology<MeshType>::VertexFace(m);
@ -357,17 +352,17 @@ class Quadric
return 2.0*c()*v + b()*u + e(); return 2.0*c()*v + b()*u + e();
} }
double duv(double u, double v) double duv(double /*u*/, double /*v*/)
{ {
return b(); return b();
} }
double duu(double u, double v) double duu(double /*u*/, double /*v*/)
{ {
return 2.0*a(); return 2.0*a();
} }
double dvv(double u, double v) double dvv(double /*u*/, double /*v*/)
{ {
return 2.0*c(); return 2.0*c();
} }