Prusa Slicer 2.6.0
Loading...
Searching...
No Matches
igl::Frame_field_deformer Class Reference
+ Collaboration diagram for igl::Frame_field_deformer:

Public Member Functions

IGL_INLINE Frame_field_deformer ()
 
IGL_INLINE ~Frame_field_deformer ()
 
IGL_INLINE void init (const Eigen::MatrixXd &_V, const Eigen::MatrixXi &_F, const Eigen::MatrixXd &_D1, const Eigen::MatrixXd &_D2, double _Lambda, double _perturb_rotations, int _fixed=1)
 
IGL_INLINE void optimize (int N, bool reset=false)
 
IGL_INLINE void reset_opt ()
 
IGL_INLINE void precompute_opt ()
 
IGL_INLINE void precompute_ARAP (Eigen::SparseMatrix< double > &Lff, Eigen::MatrixXd &LfcVc)
 
IGL_INLINE void precompute_SMOOTH (Eigen::SparseMatrix< double > &MS, Eigen::MatrixXd &bS)
 
IGL_INLINE void extractBlock (Eigen::SparseMatrix< double > &mat, int r0, int c0, int r, int c, Eigen::SparseMatrix< double > &m1)
 
IGL_INLINE void compute_optimal_rotations ()
 
IGL_INLINE void compute_optimal_positions ()
 
IGL_INLINE void computeXField (std::vector< Eigen::Matrix< double, 3, 2 > > &XF)
 
IGL_INLINE void compute_idealWarp (std::vector< Eigen::Matrix< double, 3, 3 > > &WW)
 

Public Attributes

Eigen::MatrixXd V
 
Eigen::MatrixXi F
 
std::vector< std::vector< int > > VT
 
std::vector< std::vector< int > > VTi
 
Eigen::MatrixXd V_w
 
std::vector< Eigen::Matrix< double, 3, 2 > > FF
 
std::vector< Eigen::Matrix< double, 3, 3 > > WW
 
std::vector< Eigen::Matrix< double, 3, 2 > > XF
 
int fixed
 
double perturb_rotations
 
int nfree
 
int nconst
 
Eigen::MatrixXd C
 
Eigen::SparseMatrix< double > L
 
Eigen::SparseMatrix< double > M
 
Eigen::MatrixXd RHS
 
std::vector< Eigen::Matrix< double, 3, 3 > > RW
 
Eigen::SimplicialCholesky< Eigen::SparseMatrix< double > > solver
 

Private Attributes

double Lambda = 0.1
 

Detailed Description

Constructor & Destructor Documentation

◆ Frame_field_deformer()

IGL_INLINE igl::Frame_field_deformer::Frame_field_deformer ( )
99{}

◆ ~Frame_field_deformer()

IGL_INLINE igl::Frame_field_deformer::~Frame_field_deformer ( )
101{}

Member Function Documentation

◆ compute_idealWarp()

IGL_INLINE void igl::Frame_field_deformer::compute_idealWarp ( std::vector< Eigen::Matrix< double, 3, 3 > > &  WW)
343{
344 using namespace Eigen;
345
346 WW.resize(F.rows());
347 for (int i=0;i<(int)FF.size();i++)
348 {
349 Vector3d v0,v1,v2;
350 v0 = FF[i].col(0);
351 v1 = FF[i].col(1);
352 v2=v0.cross(v1); v2.normalize(); // normal
353
354 Matrix3d A,AI; // compute affine map A that brings:
355 A << v0[0], v1[0], v2[0], // first vector of FF to x unary vector
356 v0[1], v1[1], v2[1], // second vector of FF to xy plane
357 v0[2], v1[2], v2[2]; // triangle normal to z unary vector
358 AI = A.inverse();
359
360 // polar decomposition to discard rotational component (unnecessary but makes it easier)
362 //Matrix<double,3,3> au = svd.matrixU();
363 Matrix<double,3,3> av = svd.matrixV();
364 DiagonalMatrix<double,3> as(svd.singularValues());
365 WW[i] = av*as*av.transpose();
366 }
367}
Represents a diagonal matrix with its storage.
Definition DiagonalMatrix.h:118
Two-sided Jacobi SVD decomposition of a rectangular matrix.
Definition JacobiSVD.h:489
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:180
std::vector< Eigen::Matrix< double, 3, 2 > > FF
Definition frame_field_deformer.cpp:75
std::vector< Eigen::Matrix< double, 3, 3 > > WW
Definition frame_field_deformer.cpp:76
Eigen::MatrixXi F
Definition frame_field_deformer.cpp:68
@ ComputeFullV
Definition Constants.h:387
@ ComputeFullU
Definition Constants.h:383
Definition LDLT.h:16
Eigen::Matrix< double, 3, 3, Eigen::DontAlign > Matrix3d
Definition Point.hpp:71

References Eigen::ComputeFullU, Eigen::ComputeFullV, F, FF, Eigen::SVDBase< Derived >::matrixV(), Eigen::SVDBase< Derived >::singularValues(), and WW.

Referenced by precompute_opt().

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

◆ compute_optimal_positions()

IGL_INLINE void igl::Frame_field_deformer::compute_optimal_positions ( )
275{
276 using namespace Eigen;
277 // compute variable RHS of ARAP-warp part of the system
278 MatrixXd b(nfree,3); // fx3 known term of the system
279 MatrixXd X; // result
280 int t; // triangles incident to edge (i,j)
281 int vi,i1,i2; // index of vertex i wrt tri t0
282
283 for (int i=0;i<nfree;i++)
284 {
285 b.row(i) << 0.0, 0.0, 0.0;
286 for (int k=0;k<(int)VT[i].size();k++) // for all incident triangles
287 {
288 t = VT[i][k]; // incident tri
289 vi = (i==F(t,0))?0:(i==F(t,1))?1:(i==F(t,2))?2:3; // index of i in t
290 assert(vi!=3);
291 i1 = F(t,(vi+1)%3);
292 i2 = F(t,(vi+2)%3);
293 b.row(i)+=(C(t,(vi+2)%3)*RW[t]*(V.row(i1)-V.row(i)).transpose()).transpose();
294 b.row(i)+=(C(t,(vi+1)%3)*RW[t]*(V.row(i2)-V.row(i)).transpose()).transpose();
295 }
296 }
297 b/=2.0;
298 b=-4*b;
299
300 b*=(1-Lambda); // blend
301
302 b+=RHS; // complete known term
303
304 X = solver.solve(b);
305 if (solver.info()!=Eigen::Success) {printf("Solving linear system failed!\n"); return;}
306
307 // copy result to mw.V
308 for (int i=0;i<nfree;i++)
309 V_w.row(i)=X.row(i);
310
311}
ComputationInfo info() const
Reports whether previous computation was successful.
Definition SimplicialCholesky.h:107
const Solve< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition SparseSolverBase.h:88
Eigen::SimplicialCholesky< Eigen::SparseMatrix< double > > solver
Definition frame_field_deformer.cpp:91
int nfree
Definition frame_field_deformer.cpp:84
std::vector< Eigen::Matrix< double, 3, 3 > > RW
Definition frame_field_deformer.cpp:90
std::vector< std::vector< int > > VT
Definition frame_field_deformer.cpp:70
Eigen::MatrixXd RHS
Definition frame_field_deformer.cpp:89
Eigen::MatrixXd V
Definition frame_field_deformer.cpp:67
Eigen::MatrixXd C
Definition frame_field_deformer.cpp:85
double Lambda
Definition frame_field_deformer.cpp:95
Eigen::MatrixXd V_w
Definition frame_field_deformer.cpp:73
@ Success
Definition Constants.h:432
@ X
Definition libslic3r.h:98
constexpr auto size(const C &c) -> decltype(c.size())
Definition span.hpp:183

References C, F, Eigen::SimplicialCholeskyBase< Derived >::info(), Lambda, nfree, RHS, RW, Eigen::SparseSolverBase< Derived >::solve(), solver, Eigen::Success, V, V_w, and VT.

Referenced by optimize().

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

◆ compute_optimal_rotations()

IGL_INLINE void igl::Frame_field_deformer::compute_optimal_rotations ( )
236{
237 using namespace Eigen;
238 Matrix<double,3,3> r,S,P,PP,D;
239
240 for (int i=0;i<F.rows();i++)
241 {
242 // input tri --- could be done once and saved in a matrix
243 P.col(0) = (V.row(F(i,1))-V.row(F(i,0))).transpose();
244 P.col(1) = (V.row(F(i,2))-V.row(F(i,1))).transpose();
245 P.col(2) = (V.row(F(i,0))-V.row(F(i,2))).transpose();
246
247 P = WW[i] * P; // apply ideal warp
248
249 // current tri
250 PP.col(0) = (V_w.row(F(i,1))-V_w.row(F(i,0))).transpose();
251 PP.col(1) = (V_w.row(F(i,2))-V_w.row(F(i,1))).transpose();
252 PP.col(2) = (V_w.row(F(i,0))-V_w.row(F(i,2))).transpose();
253
254 // cotangents
255 D << C(i,2), 0, 0,
256 0, C(i,0), 0,
257 0, 0, C(i,1);
258
259 S = PP*D*P.transpose();
261 Matrix<double,3,3> su = svd.matrixU();
262 Matrix<double,3,3> sv = svd.matrixV();
263 r = su*sv.transpose();
264
265 if (r.determinant()<0) // correct reflections
266 {
267 su(0,2)=-su(0,2); su(1,2)=-su(1,2); su(2,2)=-su(2,2);
268 r = su*sv.transpose();
269 }
270 RW[i] = r*WW[i]; // RW INCORPORATES IDEAL WARP WW!!!
271 }
272}

References C, Eigen::ComputeFullU, Eigen::ComputeFullV, F, Eigen::SVDBase< Derived >::matrixU(), Eigen::SVDBase< Derived >::matrixV(), RW, V, V_w, and WW.

Referenced by optimize().

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

◆ computeXField()

IGL_INLINE void igl::Frame_field_deformer::computeXField ( std::vector< Eigen::Matrix< double, 3, 2 > > &  XF)
314{
315 using namespace Eigen;
316 Matrix<double,3,3> P,PP,DG;
317 XF.resize(F.rows());
318
319 for (int i=0;i<F.rows();i++)
320 {
321 int i0,i1,i2;
322 // indexes of vertices of face i
323 i0 = F(i,0); i1 = F(i,1); i2 = F(i,2);
324
325 // input frame
326 P.col(0) = (V.row(i1)-V.row(i0)).transpose();
327 P.col(1) = (V.row(i2)-V.row(i0)).transpose();
328 P.col(2) = P.col(0).cross(P.col(1));
329
330 // output triangle brought to origin
331 PP.col(0) = (V_w.row(i1)-V_w.row(i0)).transpose();
332 PP.col(1) = (V_w.row(i2)-V_w.row(i0)).transpose();
333 PP.col(2) = PP.col(0).cross(PP.col(1));
334
335 // deformation gradient
336 DG = PP * P.inverse();
337 XF[i] = DG * FF[i];
338 }
339}
std::vector< Eigen::Matrix< double, 3, 2 > > XF
Definition frame_field_deformer.cpp:77

References F, FF, V, V_w, and XF.

Referenced by optimize().

+ Here is the caller graph for this function:

◆ extractBlock()

IGL_INLINE void igl::Frame_field_deformer::extractBlock ( Eigen::SparseMatrix< double > &  mat,
int  r0,
int  c0,
int  r,
int  c,
Eigen::SparseMatrix< double > &  m1 
)
224{
225 std::vector<Eigen::Triplet<double> > tripletList;
226 for (int k=c0; k<c0+c; ++k)
227 for (Eigen::SparseMatrix<double>::InnerIterator it(mat,k); it; ++it)
228 {
229 if (it.row()>=r0 && it.row()<r0+r)
230 tripletList.push_back(Eigen::Triplet<double>(it.row()-r0,it.col()-c0,it.value()));
231 }
232 m1.setFromTriplets(tripletList.begin(), tripletList.end());
233}
Definition SparseCompressedBase.h:137
void setFromTriplets(const InputIterators &begin, const InputIterators &end)
Definition SparseMatrix.h:993
A small structure to hold a non zero as a triplet (i,j,value).
Definition SparseUtil.h:155

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

Referenced by precompute_ARAP(), and precompute_SMOOTH().

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

◆ init()

IGL_INLINE void igl::Frame_field_deformer::init ( const Eigen::MatrixXd &  _V,
const Eigen::MatrixXi &  _F,
const Eigen::MatrixXd &  _D1,
const Eigen::MatrixXd &  _D2,
double  _Lambda,
double  _perturb_rotations,
int  _fixed = 1 
)
110{
111 V = _V;
112 F = _F;
113
114 assert(_D1.rows() == _D2.rows());
115
116 FF.clear();
117 for (unsigned i=0; i < _D1.rows(); ++i)
118 {
120 ff.col(0) = _D1.row(i);
121 ff.col(1) = _D2.row(i);
122 FF.push_back(ff);
123 }
124
125 fixed = _fixed;
126 Lambda = _Lambda;
127 perturb_rotations = _perturb_rotations;
128
129 reset_opt();
131}
IGL_INLINE void reset_opt()
Definition frame_field_deformer.cpp:149
IGL_INLINE void precompute_opt()
Definition frame_field_deformer.cpp:160
double perturb_rotations
Definition frame_field_deformer.cpp:81
int fixed
Definition frame_field_deformer.cpp:79

References F, FF, fixed, Lambda, perturb_rotations, precompute_opt(), reset_opt(), and V.

Referenced by igl::frame_field_deformer().

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

◆ optimize()

IGL_INLINE void igl::Frame_field_deformer::optimize ( int  N,
bool  reset = false 
)
135{
136 //Reset optimization
137 if (reset)
138 reset_opt();
139
140 // Iterative Local/Global optimization
141 for (int i=0; i<N;i++)
142 {
146 }
147}
IGL_INLINE void computeXField(std::vector< Eigen::Matrix< double, 3, 2 > > &XF)
Definition frame_field_deformer.cpp:313
IGL_INLINE void compute_optimal_rotations()
Definition frame_field_deformer.cpp:235
IGL_INLINE void compute_optimal_positions()
Definition frame_field_deformer.cpp:274

References compute_optimal_positions(), compute_optimal_rotations(), computeXField(), reset_opt(), and XF.

Referenced by igl::frame_field_deformer().

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

◆ precompute_ARAP()

IGL_INLINE void igl::Frame_field_deformer::precompute_ARAP ( Eigen::SparseMatrix< double > &  Lff,
Eigen::MatrixXd &  LfcVc 
)
192{
193 using namespace Eigen;
194 fprintf(stdout,"Precomputing ARAP terms\n");
195 SparseMatrix<double> LL = -4*L;
197 extractBlock(LL,0,0,nfree,nfree,Lff);
199 extractBlock(LL,0,nfree,nfree,nconst,Lfc);
200 LfcVc = - Lfc * V_w.block(nfree,0,nconst,3);
201}
A versatible sparse matrix representation.
Definition SparseMatrix.h:98
IGL_INLINE void extractBlock(Eigen::SparseMatrix< double > &mat, int r0, int c0, int r, int c, Eigen::SparseMatrix< double > &m1)
Definition frame_field_deformer.cpp:223
Eigen::SparseMatrix< double > L
Definition frame_field_deformer.cpp:86
int nconst
Definition frame_field_deformer.cpp:84

References extractBlock(), L, nconst, nfree, and V_w.

Referenced by precompute_opt().

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

◆ precompute_opt()

IGL_INLINE void igl::Frame_field_deformer::precompute_opt ( )
161{
162 using namespace Eigen;
163 nfree = V.rows() - fixed; // free vertices (at the beginning ov m.V) - global
164 nconst = V.rows()-nfree; // #constrained vertices
165 igl::vertex_triangle_adjacency(V,F,VT,VTi); // compute vertex to face relationship
166
167 igl::cotmatrix_entries(V,F,C); // cotangent matrix for opt. rotations - global
168
170
171 SparseMatrix<double> MA; // internal matrix for ARAP-warping energy
172 MatrixXd LfcVc; // RHS (partial) for ARAP-warping energy
173 SparseMatrix<double> MS; // internal matrix for smoothing energy
174 MatrixXd bS; // RHS (full) for smoothing energy
175
176 precompute_ARAP(MA,LfcVc); // precompute terms for the ARAP-warp part
177 precompute_SMOOTH(MS,bS); // precompute terms for the smoothing part
178 compute_idealWarp(WW); // computes the ideal warps
179 RW.resize(F.rows()); // init rotation matrices - global
180
181 M = (1-Lambda)*MA + Lambda*MS; // matrix for linear system - global
182
183 RHS = (1-Lambda)*LfcVc + Lambda*bS; // RHS (partial) for linear system - global
184 solver.compute(M); // system pre-conditioning
185 if (solver.info()!=Eigen::Success) {fprintf(stderr,"Decomposition failed in pre-conditioning!\n"); exit(-1);}
186
187 fprintf(stdout,"Preconditioning done.\n");
188
189}
SimplicialCholesky & compute(const MatrixType &matrix)
Definition SimplicialCholesky.h:553
IGL_INLINE void precompute_ARAP(Eigen::SparseMatrix< double > &Lff, Eigen::MatrixXd &LfcVc)
Definition frame_field_deformer.cpp:191
std::vector< std::vector< int > > VTi
Definition frame_field_deformer.cpp:71
IGL_INLINE void precompute_SMOOTH(Eigen::SparseMatrix< double > &MS, Eigen::MatrixXd &bS)
Definition frame_field_deformer.cpp:203
Eigen::SparseMatrix< double > M
Definition frame_field_deformer.cpp:88
IGL_INLINE void compute_idealWarp(std::vector< Eigen::Matrix< double, 3, 3 > > &WW)
Definition frame_field_deformer.cpp:342
IGL_INLINE void vertex_triangle_adjacency(const typename DerivedF::Scalar n, const Eigen::MatrixBase< DerivedF > &F, std::vector< std::vector< VFType > > &VF, std::vector< std::vector< VFiType > > &VFi)
Definition vertex_triangle_adjacency.cpp:12
IGL_INLINE void cotmatrix(const Eigen::MatrixBase< DerivedV > &V, const Eigen::MatrixBase< DerivedF > &F, Eigen::SparseMatrix< Scalar > &L)
Definition cotmatrix.cpp:19
IGL_INLINE void cotmatrix_entries(const Eigen::MatrixBase< DerivedV > &V, const Eigen::MatrixBase< DerivedF > &F, Eigen::PlainObjectBase< DerivedC > &C)
Definition cotmatrix_entries.cpp:20

References C, Eigen::SimplicialCholesky< _MatrixType, _UpLo, _Ordering >::compute(), compute_idealWarp(), igl::cotmatrix(), igl::cotmatrix_entries(), F, fixed, Eigen::SimplicialCholeskyBase< Derived >::info(), L, Lambda, M, nconst, nfree, precompute_ARAP(), precompute_SMOOTH(), RHS, RW, solver, Eigen::Success, V, igl::vertex_triangle_adjacency(), VT, VTi, and WW.

Referenced by init().

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

◆ precompute_SMOOTH()

IGL_INLINE void igl::Frame_field_deformer::precompute_SMOOTH ( Eigen::SparseMatrix< double > &  MS,
Eigen::MatrixXd &  bS 
)
204{
205 using namespace Eigen;
206 fprintf(stdout,"Precomputing SMOOTH terms\n");
207
208 SparseMatrix<double> LL = 4*L*L;
209
210 // top-left
212 extractBlock(LL,0,0,nfree,nfree,MS);
213
214 // top-right
216 extractBlock(LL,0,nfree,nfree,nconst,Mfc);
217
218 MatrixXd MfcVc = Mfc * V_w.block(nfree,0,nconst,3);
219 bS = (LL*V).block(0,0,nfree,3)-MfcVc;
220
221}
EIGEN_DOC_BLOCK_ADDONS_NOT_INNER_PANEL EIGEN_DEVICE_FUNC BlockXpr block(Index startRow, Index startCol, Index blockRows, Index blockCols)
This is the const version of block(Index,Index,Index,Index). *‍/.
Definition BlockMethods.h:64

References block(), extractBlock(), L, nconst, nfree, V, and V_w.

Referenced by precompute_opt().

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

◆ reset_opt()

IGL_INLINE void igl::Frame_field_deformer::reset_opt ( )
150{
151 V_w = V;
152
153 for (unsigned i=0; i<V_w.rows(); ++i)
154 for (unsigned j=0; j<V_w.cols(); ++j)
155 V_w(i,j) += (double(rand())/double(RAND_MAX))*10e-4*perturb_rotations;
156
157}

References perturb_rotations, V, and V_w.

Referenced by init(), and optimize().

+ Here is the caller graph for this function:

Member Data Documentation

◆ C

Eigen::MatrixXd igl::Frame_field_deformer::C

◆ F

Eigen::MatrixXi igl::Frame_field_deformer::F

◆ FF

std::vector< Eigen::Matrix<double,3,2> > igl::Frame_field_deformer::FF

◆ fixed

int igl::Frame_field_deformer::fixed

Referenced by init(), and precompute_opt().

◆ L

Eigen::SparseMatrix<double> igl::Frame_field_deformer::L

◆ Lambda

double igl::Frame_field_deformer::Lambda = 0.1
private

◆ M

Eigen::SparseMatrix<double> igl::Frame_field_deformer::M

Referenced by precompute_opt().

◆ nconst

int igl::Frame_field_deformer::nconst

◆ nfree

int igl::Frame_field_deformer::nfree

◆ perturb_rotations

double igl::Frame_field_deformer::perturb_rotations

Referenced by init(), and reset_opt().

◆ RHS

Eigen::MatrixXd igl::Frame_field_deformer::RHS

◆ RW

std::vector< Eigen::Matrix<double,3,3> > igl::Frame_field_deformer::RW

◆ solver

Eigen::SimplicialCholesky<Eigen::SparseMatrix<double> > igl::Frame_field_deformer::solver

◆ V

◆ V_w

◆ VT

std::vector<std::vector<int> > igl::Frame_field_deformer::VT

◆ VTi

std::vector<std::vector<int> > igl::Frame_field_deformer::VTi

Referenced by precompute_opt().

◆ WW

std::vector< Eigen::Matrix<double,3,3> > igl::Frame_field_deformer::WW

◆ XF

std::vector< Eigen::Matrix<double,3,2> > igl::Frame_field_deformer::XF

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