Prusa Slicer 2.6.0
Loading...
Searching...
No Matches
igl::PlanarizerShapeUp< DerivedV, DerivedF > Class Template Reference
+ Collaboration diagram for igl::PlanarizerShapeUp< DerivedV, DerivedF >:

Public Member Functions

 PlanarizerShapeUp (const Eigen::PlainObjectBase< DerivedV > &V_, const Eigen::PlainObjectBase< DerivedF > &F_, const int maxIter_, const double &threshold_)
 
void planarize (Eigen::PlainObjectBase< DerivedV > &Vout)
 

Protected Member Functions

void assembleQ ()
 
void assembleP ()
 
void assembleNi ()
 
void assembleSelector (int fi, Eigen::SparseMatrix< typename DerivedV::Scalar > &S)
 

Protected Attributes

long numV
 
long numF
 
const Eigen::PlainObjectBase< DerivedV > & Vin
 
const Eigen::PlainObjectBase< DerivedF > & Fin
 
Eigen::Matrix< typename DerivedV::Scalar, Eigen::Dynamic, 1 > Vv
 
Eigen::Matrix< typename DerivedV::Scalar, Eigen::Dynamic, 1 > weightsSqrt
 
Eigen::Matrix< typename DerivedV::Scalar, Eigen::Dynamic, 1 > P
 
Eigen::SparseMatrix< typename DerivedV::Scalar > Q
 
Eigen::SparseMatrix< typename DerivedV::Scalar > Ni
 
Eigen::SimplicialLDLT< Eigen::SparseMatrix< typename DerivedV::Scalar > > solver
 
int maxIter
 
double threshold
 
const int ni = 4
 

Detailed Description

template<typename DerivedV, typename DerivedF>
class igl::PlanarizerShapeUp< DerivedV, DerivedF >

Constructor & Destructor Documentation

◆ PlanarizerShapeUp()

template<typename DerivedV , typename DerivedF >
igl::PlanarizerShapeUp< DerivedV, DerivedF >::PlanarizerShapeUp ( const Eigen::PlainObjectBase< DerivedV > &  V_,
const Eigen::PlainObjectBase< DerivedF > &  F_,
const int  maxIter_,
const double &  threshold_ 
)
inline
68 :
69numV(V_.rows()),
70numF(F_.rows()),
71Vin(V_),
72Fin(F_),
74maxIter(maxIter_),
75threshold(threshold_)
76{
77 // assemble stacked vertex position vector
78 Vv.setZero(3*numV,1);
79 for (int i =0;i<numV;++i)
80 Vv.segment(3*i,3) = Vin.row(i);
81 // assemble and factorize lhs matrix
82 assembleQ();
83};
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:180
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Index rows() const
Definition PlainObjectBase.h:151
EIGEN_DEVICE_FUNC Derived & setZero(Index size)
Definition CwiseNullaryOp.h:515
Eigen::Matrix< typename DerivedV::Scalar, Eigen::Dynamic, 1 > weightsSqrt
Definition planarize_quad_mesh.cpp:29
long numF
Definition planarize_quad_mesh.cpp:21
const Eigen::PlainObjectBase< DerivedF > & Fin
Definition planarize_quad_mesh.cpp:24
const Eigen::PlainObjectBase< DerivedV > & Vin
Definition planarize_quad_mesh.cpp:23
int maxIter
Definition planarize_quad_mesh.cpp:37
void assembleQ()
Definition planarize_quad_mesh.cpp:86
Eigen::Matrix< typename DerivedV::Scalar, Eigen::Dynamic, 1 > Vv
Definition planarize_quad_mesh.cpp:29
long numV
Definition planarize_quad_mesh.cpp:21
double threshold
Definition planarize_quad_mesh.cpp:38

References igl::PlanarizerShapeUp< DerivedV, DerivedF >::assembleQ(), igl::PlanarizerShapeUp< DerivedV, DerivedF >::numV, Eigen::PlainObjectBase< Derived >::setZero(), igl::PlanarizerShapeUp< DerivedV, DerivedF >::Vin, and igl::PlanarizerShapeUp< DerivedV, DerivedF >::Vv.

+ Here is the call graph for this function:

Member Function Documentation

◆ assembleNi()

template<typename DerivedV , typename DerivedF >
void igl::PlanarizerShapeUp< DerivedV, DerivedF >::assembleNi
inlineprotected
126{
127 std::vector<Eigen::Triplet<typename DerivedV::Scalar>> tripletList;
128 for (int ii = 0; ii< ni; ii++)
129 {
130 for (int jj = 0; jj< ni; jj++)
131 {
132 tripletList.push_back(Eigen::Triplet<typename DerivedV::Scalar>(3*ii+0,3*jj+0,-1./ni));
133 tripletList.push_back(Eigen::Triplet<typename DerivedV::Scalar>(3*ii+1,3*jj+1,-1./ni));
134 tripletList.push_back(Eigen::Triplet<typename DerivedV::Scalar>(3*ii+2,3*jj+2,-1./ni));
135 }
136 tripletList.push_back(Eigen::Triplet<typename DerivedV::Scalar>(3*ii+0,3*ii+0,1.));
137 tripletList.push_back(Eigen::Triplet<typename DerivedV::Scalar>(3*ii+1,3*ii+1,1.));
138 tripletList.push_back(Eigen::Triplet<typename DerivedV::Scalar>(3*ii+2,3*ii+2,1.));
139 }
140 Ni.resize(3*ni,3*ni);
141 Ni.setFromTriplets(tripletList.begin(), tripletList.end());
142}
void setFromTriplets(const InputIterators &begin, const InputIterators &end)
Definition SparseMatrix.h:993
void resize(Index rows, Index cols)
Definition SparseMatrix.h:621
A small structure to hold a non zero as a triplet (i,j,value).
Definition SparseUtil.h:155
const int ni
Definition planarize_quad_mesh.cpp:39
Eigen::SparseMatrix< typename DerivedV::Scalar > Ni
Definition planarize_quad_mesh.cpp:34

◆ assembleP()

template<typename DerivedV , typename DerivedF >
void igl::PlanarizerShapeUp< DerivedV, DerivedF >::assembleP
inlineprotected
167{
168 P.setZero(3*ni*numF);
169 for (int fi = 0; fi< numF; fi++)
170 {
171 // todo: this can be made faster by omitting the selector matrix
173 assembleSelector(fi, Sfi);
175
178 for (int i = 0; i <ni; ++i)
179 CC.col(i) = Vi.segment(3*i, 3);
181
182 // Alec: Doesn't compile
184 // the real() is for compilation purposes
185 Eigen::Matrix<typename DerivedV::Scalar, 3, 1> lambda = es.eigenvalues().real();
186 Eigen::Matrix<typename DerivedV::Scalar, 3, 3> U = es.eigenvectors().real();
187 int min_i;
188 lambda.cwiseAbs().minCoeff(&min_i);
189 U.col(min_i).setZero();
191 for (int i = 0; i <ni; ++i)
192 P.segment(3*ni*fi+3*i, 3) = weightsSqrt[fi]*PP.col(i);
193
194 }
195}
Computes eigenvalues and eigenvectors of general matrices.
Definition EigenSolver.h:65
A versatible sparse matrix representation.
Definition SparseMatrix.h:98
Eigen::Matrix< typename DerivedV::Scalar, Eigen::Dynamic, 1 > P
Definition planarize_quad_mesh.cpp:29
void assembleSelector(int fi, Eigen::SparseMatrix< typename DerivedV::Scalar > &S)
Definition planarize_quad_mesh.cpp:146

References Eigen::EigenSolver< _MatrixType >::eigenvalues(), Eigen::EigenSolver< _MatrixType >::eigenvectors(), and Eigen::PlainObjectBase< Derived >::setZero().

+ Here is the call graph for this function:

◆ assembleQ()

template<typename DerivedV , typename DerivedF >
void igl::PlanarizerShapeUp< DerivedV, DerivedF >::assembleQ
inlineprotected
87{
88 std::vector<Eigen::Triplet<typename DerivedV::Scalar> > tripletList;
89
90 // assemble the Ni matrix
91 assembleNi();
92
93 for (int fi = 0; fi< numF; fi++)
94 {
96 assembleSelector(fi, Sfi);
97
98 // the final matrix per face
100 // put it in the correct block of Q
101 // todo: this can be made faster by omitting the selector matrix
102 for (int k=0; k<Qi.outerSize(); ++k)
104 {
105 typename DerivedV::Scalar val = it.value();
106 int row = it.row();
107 int col = it.col();
108 tripletList.push_back(Eigen::Triplet<typename DerivedV::Scalar>(row+3*ni*fi,col,val));
109 }
110 }
111
112 Q.resize(3*ni*numF,3*numV);
113 Q.setFromTriplets(tripletList.begin(), tripletList.end());
114 // the actual lhs matrix is Q'*Q
115 // prefactor that matrix
118 {
119 std::cerr << "Cholesky failed - PlanarizerShapeUp.cpp" << std::endl;
120 assert(0);
121 }
122}
EIGEN_DEVICE_FUNC RowXpr row(Index i)
This is the const version of row(). *‍/.
Definition BlockMethods.h:859
EIGEN_DEVICE_FUNC ColXpr col(Index i)
This is the const version of col().
Definition BlockMethods.h:838
ComputationInfo info() const
Reports whether previous computation was successful.
Definition SimplicialCholesky.h:107
SimplicialLDLT & compute(const MatrixType &matrix)
Definition SimplicialCholesky.h:461
Definition SparseCompressedBase.h:137
TransposeReturnType transpose()
Definition SparseMatrixBase.h:349
Index outerSize() const
Definition SparseMatrix.h:143
void assembleNi()
Definition planarize_quad_mesh.cpp:125
Eigen::SparseMatrix< typename DerivedV::Scalar > Q
Definition planarize_quad_mesh.cpp:34
Eigen::SimplicialLDLT< Eigen::SparseMatrix< typename DerivedV::Scalar > > solver
Definition planarize_quad_mesh.cpp:35
@ Success
Definition Constants.h:432

References col(), Eigen::SparseMatrix< _Scalar, _Options, _StorageIndex >::outerSize(), row(), and Eigen::Success.

Referenced by igl::PlanarizerShapeUp< DerivedV, DerivedF >::PlanarizerShapeUp().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ assembleSelector()

template<typename DerivedV , typename DerivedF >
void igl::PlanarizerShapeUp< DerivedV, DerivedF >::assembleSelector ( int  fi,
Eigen::SparseMatrix< typename DerivedV::Scalar > &  S 
)
inlineprotected
148{
149
150 std::vector<Eigen::Triplet<typename DerivedV::Scalar>> tripletList;
151 for (int fvi = 0; fvi< ni; fvi++)
152 {
153 int vi = Fin(fi,fvi);
154 tripletList.push_back(Eigen::Triplet<typename DerivedV::Scalar>(3*fvi+0,3*vi+0,1.));
155 tripletList.push_back(Eigen::Triplet<typename DerivedV::Scalar>(3*fvi+1,3*vi+1,1.));
156 tripletList.push_back(Eigen::Triplet<typename DerivedV::Scalar>(3*fvi+2,3*vi+2,1.));
157 }
158
159 S.resize(3*ni,3*numV);
160 S.setFromTriplets(tripletList.begin(), tripletList.end());
161
162}

References Eigen::SparseMatrix< _Scalar, _Options, _StorageIndex >::resize(), and Eigen::SparseMatrix< _Scalar, _Options, _StorageIndex >::setFromTriplets().

+ Here is the call graph for this function:

◆ planarize()

template<typename DerivedV , typename DerivedF >
void igl::PlanarizerShapeUp< DerivedV, DerivedF >::planarize ( Eigen::PlainObjectBase< DerivedV > &  Vout)
inline
200{
202 Vout = Vin;
203
204 for (int iter =0; iter<maxIter; ++iter)
205 {
206 igl::quad_planarity(Vout, Fin, planarity);
207 typename DerivedV::Scalar nonPlanarity = planarity.cwiseAbs().maxCoeff();
208 //std::cerr<<"iter #"<<iter<<": max non-planarity: "<<nonPlanarity<<std::endl;
209 if (nonPlanarity<threshold)
210 break;
211 assembleP();
212 Vv = solver.solve(Q.transpose()*P);
214 {
215 std::cerr << "Linear solve failed - PlanarizerShapeUp.cpp" << std::endl;
216 assert(0);
217 }
218 for (int i =0;i<numV;++i)
219 Vout.row(i) << Vv.segment(3*i,3).transpose();
220 }
221 // set the mean of Vout to the mean of Vin
223 oldMean = Vin.colwise().mean();
224 newMean = Vout.colwise().mean();
225 Vout.rowwise() += (oldMean - newMean);
226
227};
const Solve< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition SparseSolverBase.h:88
void assembleP()
Definition planarize_quad_mesh.cpp:166
IGL_INLINE void quad_planarity(const Eigen::PlainObjectBase< DerivedV > &V, const Eigen::PlainObjectBase< DerivedF > &F, Eigen::PlainObjectBase< DerivedP > &P)
Definition quad_planarity.cpp:12

References igl::quad_planarity(), and Eigen::Success.

Referenced by igl::planarize_quad_mesh().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

Member Data Documentation

◆ Fin

template<typename DerivedV , typename DerivedF >
const Eigen::PlainObjectBase<DerivedF>& igl::PlanarizerShapeUp< DerivedV, DerivedF >::Fin
protected

◆ maxIter

template<typename DerivedV , typename DerivedF >
int igl::PlanarizerShapeUp< DerivedV, DerivedF >::maxIter
protected

◆ Ni

template<typename DerivedV , typename DerivedF >
Eigen::SparseMatrix<typename DerivedV::Scalar > igl::PlanarizerShapeUp< DerivedV, DerivedF >::Ni
protected

◆ ni

template<typename DerivedV , typename DerivedF >
const int igl::PlanarizerShapeUp< DerivedV, DerivedF >::ni = 4
protected

◆ numF

template<typename DerivedV , typename DerivedF >
long igl::PlanarizerShapeUp< DerivedV, DerivedF >::numF
protected

◆ numV

template<typename DerivedV , typename DerivedF >
long igl::PlanarizerShapeUp< DerivedV, DerivedF >::numV
protected

◆ P

template<typename DerivedV , typename DerivedF >
Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 1> igl::PlanarizerShapeUp< DerivedV, DerivedF >::P
protected

◆ Q

template<typename DerivedV , typename DerivedF >
Eigen::SparseMatrix<typename DerivedV::Scalar > igl::PlanarizerShapeUp< DerivedV, DerivedF >::Q
protected

◆ solver

template<typename DerivedV , typename DerivedF >
Eigen::SimplicialLDLT<Eigen::SparseMatrix<typename DerivedV::Scalar > > igl::PlanarizerShapeUp< DerivedV, DerivedF >::solver
protected

◆ threshold

template<typename DerivedV , typename DerivedF >
double igl::PlanarizerShapeUp< DerivedV, DerivedF >::threshold
protected

◆ Vin

template<typename DerivedV , typename DerivedF >
const Eigen::PlainObjectBase<DerivedV>& igl::PlanarizerShapeUp< DerivedV, DerivedF >::Vin
protected

◆ Vv

template<typename DerivedV , typename DerivedF >
Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 1> igl::PlanarizerShapeUp< DerivedV, DerivedF >::Vv
protected

◆ weightsSqrt

template<typename DerivedV , typename DerivedF >
Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 1> igl::PlanarizerShapeUp< DerivedV, DerivedF >::weightsSqrt
protected

The documentation for this class was generated from the following file: