vcglib/vcg/complex/algorithms/implicit_smooth.h

330 lines
12 KiB
C
Raw Normal View History

2014-11-05 16:48:34 +01:00
/****************************************************************************
* VCGLib o o *
* Visual and Computer Graphics Library o o *
* _ O _ *
* Copyright(C) 2004 \/)\/ *
* Visual Computing Lab /\/| *
* ISTI - Italian National Research Council | *
* \ *
* All rights reserved. *
* *
* This program is free software; you can redistribute it and/or modify *
* it under the terms of the GNU General Public License as published by *
* the Free Software Foundation; either version 2 of the License, or *
* (at your option) any later version. *
* *
* This program is distributed in the hope that it will be useful, *
* but WITHOUT ANY WARRANTY; without even the implied warranty of *
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
* GNU General Public License (http://www.gnu.org/licenses/gpl.txt) *
* for more details. *
* *
****************************************************************************/
#ifndef __VCG_IMPLICIT_SMOOTHER
#define __VCG_IMPLICIT_SMOOTHER
#include <eigenlib/Eigen/Sparse>
#include <vcg/complex/algorithms/mesh_to_matrix.h>
#include <vcg/complex/algorithms/update/quality.h>
#include <vcg/complex/algorithms/smooth.h>
#define PENALTY 10000
namespace vcg{
template <class MeshType>
class ImplicitSmoother
{
typedef typename MeshType::FaceType FaceType;
typedef typename MeshType::VertexType VertexType;
typedef typename MeshType::CoordType CoordType;
typedef typename MeshType::ScalarType ScalarType;
typedef typename Eigen::Matrix<ScalarType, Eigen::Dynamic, Eigen::Dynamic> MatrixXm;
public:
struct FaceConstraint
{
int numF;
std::vector<ScalarType > BarycentricW;
2014-11-05 16:48:34 +01:00
CoordType TargetPos;
ScalarType facePenalty;
2014-11-05 16:48:34 +01:00
FaceConstraint()
{
numF=-1;
facePenalty=ScalarType(PENALTY);
2014-11-05 16:48:34 +01:00
}
FaceConstraint(int _numF,
const std::vector<ScalarType > &_BarycentricW,
const CoordType &_TargetPos,
const ScalarType fPenalty = ScalarType(PENALTY))
2014-11-05 16:48:34 +01:00
{
numF=_numF;
BarycentricW= std::vector<ScalarType > (_BarycentricW.begin(),_BarycentricW.end());
2014-11-05 16:48:34 +01:00
TargetPos=_TargetPos;
facePenalty=fPenalty;
2014-11-05 16:48:34 +01:00
}
};
struct Parameter
2014-11-05 16:48:34 +01:00
{
//the amount of smoothness, useful only if we set the mass matrix
ScalarType lambda;
//the use of mass matrix to keep the mesh close to its original position
//(weighted per area distributed on vertices)
bool useMassMatrix;
//this bool is used to fix the border vertices of the mesh or not
bool fixBorder;
//this bool is used to set if cotangent weight is used, this flag to false means uniform laplacian
bool useCotWeight;
//the set of fixed vertices
std::vector<int> FixedV;
//the set of faces for barycentric constraints
std::vector<FaceConstraint> ConstrainedF;
//the degree of laplacian
int degree;
2014-11-05 16:48:34 +01:00
Parameter()
2014-11-05 16:48:34 +01:00
{
degree=1;
2014-11-05 16:48:34 +01:00
lambda=0.2;
useMassMatrix=true;
fixBorder=false;
useCotWeight=false;
}
};
private:
static void InitSparse(const std::vector<std::pair<int,int> > &Index,
const std::vector<ScalarType> &Values,
const int m,
const int n,
Eigen::SparseMatrix<ScalarType>& X)
2014-11-05 16:48:34 +01:00
{
assert(Index.size()==Values.size());
std::vector<Eigen::Triplet<ScalarType> > IJV;
IJV.reserve(Index.size());
for(size_t i= 0;i<Index.size();i++)
{
int row=Index[i].first;
int col=Index[i].second;
ScalarType val=Values[i];
assert(row<m);
assert(col<n);
IJV.push_back(Eigen::Triplet<ScalarType>(row,col,val));
}
X.resize(m,n);
X.setFromTriplets(IJV.begin(),IJV.end());
}
static void CollectHardConstraints(MeshType &mesh,const Parameter &SParam,
std::vector<std::pair<int,int> > &IndexC,
std::vector<ScalarType> &WeightC)
2014-11-05 16:48:34 +01:00
{
std::vector<int> To_Fix;
//collect fixed vert
if (SParam.fixBorder)
{
//add penalization constra
for (size_t i=0;i<mesh.vert.size();i++)
2014-11-05 16:48:34 +01:00
{
if (!mesh.vert[i].IsB())continue;
2014-11-05 16:48:34 +01:00
To_Fix.push_back(i);
}
}
//add additional fixed vertices constraint
To_Fix.insert(To_Fix.end(),SParam.FixedV.begin(),SParam.FixedV.end());
//sort and make them unique
std::sort(To_Fix.begin(),To_Fix.end());
typename std::vector<int>::iterator it= std::unique (To_Fix.begin(), To_Fix.end());
To_Fix.resize( std::distance(To_Fix.begin(),it) );
for (size_t i=0;i<To_Fix.size();i++)
{
for (int j=0;j<3;j++)
{
int IndexV=(To_Fix[i]*3)+j;
2014-11-05 16:48:34 +01:00
IndexC.push_back(std::pair<int,int>(IndexV,IndexV));
WeightC.push_back((ScalarType)PENALTY);
}
}
}
static void CollectBarycentricConstraints(MeshType &mesh,
const Parameter &SParam,
std::vector<std::pair<int,int> > &IndexC,
std::vector<ScalarType> &WeightC,
std::vector<int> &IndexRhs,
std::vector<ScalarType> &ValueRhs)
{
ScalarType penalty;
int baseIndex=mesh.vert.size();
for (size_t i=0;i<SParam.ConstrainedF.size();i++)
{
//get the index of the current constraint
int IndexConstraint=baseIndex+i;
//add one hard constraint
int FaceN=SParam.ConstrainedF[i].numF;
assert(FaceN>=0);
assert(FaceN<(int)mesh.face.size());
assert(mesh.face[FaceN].VN()==(int)SParam.ConstrainedF[i].BarycentricW.size());
penalty=SParam.ConstrainedF[i].facePenalty;
//then add all the weights to impose the constraint
for (int j=0;j<mesh.face[FaceN].VN();j++)
{
//get the current weight
ScalarType currW=SParam.ConstrainedF[i].BarycentricW[j];
2014-11-05 16:48:34 +01:00
//get the index of the current vertex
int FaceVert=vcg::tri::Index(mesh,mesh.face[FaceN].V(j));
2014-11-05 16:48:34 +01:00
//then add the constraints componentwise
for (int k=0;k<3;k++)
{
//multiply times 3 per component
int IndexV=(FaceVert*3)+k;
2014-11-05 16:48:34 +01:00
//get the index of the current constraint
int ComponentConstraint=(IndexConstraint*3)+k;
IndexC.push_back(std::pair<int,int>(ComponentConstraint,IndexV));
2014-11-05 16:48:34 +01:00
WeightC.push_back(currW*penalty);
2014-11-05 16:48:34 +01:00
IndexC.push_back(std::pair<int,int>(IndexV,ComponentConstraint));
WeightC.push_back(currW*penalty);
2014-11-05 16:48:34 +01:00
//this to avoid the 1 on diagonal last entry of mass matrix
IndexC.push_back(std::pair<int,int>(ComponentConstraint,ComponentConstraint));
WeightC.push_back(-1);
}
}
for (int j=0;j<3;j++)
{
//get the index of the current constraint
int ComponentConstraint=(IndexConstraint*3)+j;
2014-11-05 16:48:34 +01:00
//get per component value
ScalarType ComponentV=SParam.ConstrainedF[i].TargetPos.V(j);
2014-11-05 16:48:34 +01:00
//add the diagonal value
IndexRhs.push_back(ComponentConstraint);
ValueRhs.push_back(ComponentV*penalty);
}
2014-11-05 16:48:34 +01:00
}
}
2014-11-05 16:48:34 +01:00
public:
2014-11-05 16:48:34 +01:00
static void Compute(MeshType &mesh, Parameter &SParam)
2014-11-05 16:48:34 +01:00
{
//calculate the size of the system
int matr_size=mesh.vert.size()+SParam.ConstrainedF.size();
2014-11-05 16:48:34 +01:00
//the laplacian and the mass matrix
Eigen::SparseMatrix<ScalarType> L,M,B;
2014-11-05 16:48:34 +01:00
//initialize the mass matrix
std::vector<std::pair<int,int> > IndexM;
std::vector<ScalarType> ValuesM;
//add the entries for mass matrix
if (SParam.useMassMatrix) MeshToMatrix<MeshType>::MassMatrixEntry(mesh,IndexM,ValuesM);
//then add entries for lagrange mult due to barycentric constraints
for (size_t i=0;i<SParam.ConstrainedF.size();i++)
{
int baseIndex=(mesh.vert.size()+i)*3;
for (int j=0;j<3;j++)
{
IndexM.push_back(std::pair<int,int>(baseIndex+j,baseIndex+j));
ValuesM.push_back(1);
}
}
//add the hard constraints
CollectHardConstraints(mesh,SParam,IndexM,ValuesM);
//initialize sparse mass matrix
InitSparse(IndexM,ValuesM,matr_size*3,matr_size*3,M);
2014-11-05 16:48:34 +01:00
//initialize the barycentric matrix
std::vector<std::pair<int,int> > IndexB;
std::vector<ScalarType> ValuesB;
2014-11-05 16:48:34 +01:00
std::vector<int> IndexRhs;
std::vector<ScalarType> ValuesRhs;
//then also collect hard constraints
CollectBarycentricConstraints(mesh,SParam,IndexB,ValuesB,IndexRhs,ValuesRhs);
//initialize sparse constraint matrix
InitSparse(IndexB,ValuesB,matr_size*3,matr_size*3,B);
2014-11-05 16:48:34 +01:00
//get the entries for laplacian matrix
std::vector<std::pair<int,int> > IndexL;
std::vector<ScalarType> ValuesL;
MeshToMatrix<MeshType>::GetLaplacianMatrix(mesh,IndexL,ValuesL,SParam.useCotWeight);
2014-11-05 16:48:34 +01:00
//initialize sparse laplacian matrix
InitSparse(IndexL,ValuesL,matr_size*3,matr_size*3,L);
for (int i=0;i<(SParam.degree-1);i++)L=L*L;
2014-11-05 16:48:34 +01:00
//then solve the system
Eigen::SparseMatrix<ScalarType> S = (M + B + SParam.lambda*L);
2014-11-05 16:48:34 +01:00
//SimplicialLDLT
Eigen::SimplicialCholesky<Eigen::SparseMatrix<ScalarType > > solver(S);
2014-11-05 16:48:34 +01:00
assert(solver.info() == Eigen::Success);
MatrixXm V(matr_size*3,1);
2014-11-05 16:48:34 +01:00
//set the first part of the matrix with vertex values
for (size_t i=0;i<mesh.vert.size();i++)
2014-11-05 16:48:34 +01:00
{
int index=i*3;
V(index,0)=mesh.vert[i].P().X();
V(index+1,0)=mesh.vert[i].P().Y();
V(index+2,0)=mesh.vert[i].P().Z();
2014-11-05 16:48:34 +01:00
}
//then set the second part by considering RHS gien by barycentric constraint
for (size_t i=0;i<IndexRhs.size();i++)
{
int index=IndexRhs[i];
ScalarType val=ValuesRhs[i];
V(index,0)=val;
}
//solve the system
2014-11-05 16:48:34 +01:00
V = solver.solve(M*V).eval();
//then copy back values
for (size_t i=0;i<mesh.vert.size();i++)
2014-11-05 16:48:34 +01:00
{
int index=i*3;
mesh.vert[i].P().X()=V(index,0);
mesh.vert[i].P().Y()=V(index+1,0);
mesh.vert[i].P().Z()=V(index+2,0);
2014-11-05 16:48:34 +01:00
}
}
};
}//end namespace vcg
#endif