planarize_quad_mesh.cpp 8.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233
  1. #include "planarize_quad_mesh.h"
  2. #include "quad_planarity.h"
  3. #include <Eigen/Sparse>
  4. #include <Eigen/Eigenvalues>
  5. #include <iostream>
  6. namespace igl
  7. {
  8. template <typename DerivedV, typename DerivedF>
  9. class PlanarizerShapeUp
  10. {
  11. protected:
  12. // number of faces, number of vertices
  13. long numV, numF;
  14. // references to the input faces and vertices
  15. const Eigen::PlainObjectBase<DerivedV> &Vin;
  16. const Eigen::PlainObjectBase<DerivedF> &Fin;
  17. // vector consisting of the vertex positions stacked: [x;y;z;x;y;z...]
  18. // vector consisting of a weight per face (currently all set to 1)
  19. // vector consisting of the projected face vertices (might be different for the same vertex belonging to different faces)
  20. Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 1> Vv, weightsSqrt, P;
  21. // Matrices as in the paper
  22. // Q: lhs matrix
  23. // Ni: matrix that subtracts the mean of a face from the 4 vertices of a face
  24. Eigen::SparseMatrix<typename DerivedV::Scalar > Q, Ni;
  25. Eigen::SimplicialLDLT<Eigen::SparseMatrix<typename DerivedV::Scalar > > solver;
  26. int maxIter;
  27. double threshold;
  28. const int ni = 4;
  29. // Matrix assemblers
  30. inline void assembleQ();
  31. inline void assembleP();
  32. inline void assembleNi();
  33. // Selects out of Vv the 4 vertices belonging to face fi
  34. inline void assembleSelector(int fi,
  35. Eigen::SparseMatrix<typename DerivedV::Scalar > &S);
  36. public:
  37. // Init - assemble stacked vector and lhs matrix, factorize
  38. inline PlanarizerShapeUp(const Eigen::PlainObjectBase<DerivedV> &V_,
  39. const Eigen::PlainObjectBase<DerivedF> &F_,
  40. const int maxIter_,
  41. const double &threshold_);
  42. // Planarization - output to Vout
  43. inline void planarize(Eigen::PlainObjectBase<DerivedV> &Vout);
  44. };
  45. }
  46. //Implementation
  47. template <typename DerivedV, typename DerivedF>
  48. inline igl::PlanarizerShapeUp<DerivedV, DerivedF>::PlanarizerShapeUp(const Eigen::PlainObjectBase<DerivedV> &V_,
  49. const Eigen::PlainObjectBase<DerivedF> &F_,
  50. const int maxIter_,
  51. const double &threshold_):
  52. Vin(V_),
  53. Fin(F_),
  54. maxIter(maxIter_),
  55. threshold(threshold_),
  56. numF(F_.rows()),
  57. numV(V_.rows()),
  58. weightsSqrt(Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 1>::Ones(numF,1))
  59. {
  60. // assemble stacked vertex position vector
  61. Vv.setZero(3*numV,1);
  62. for (int i =0;i<numV;++i)
  63. Vv.segment(3*i,3) = Vin.row(i);
  64. // assemble and factorize lhs matrix
  65. assembleQ();
  66. };
  67. template <typename DerivedV, typename DerivedF>
  68. inline void igl::PlanarizerShapeUp<DerivedV, DerivedF>::assembleQ()
  69. {
  70. std::vector<Eigen::Triplet<typename DerivedV::Scalar> > tripletList;
  71. // assemble the Ni matrix
  72. assembleNi();
  73. for (int fi = 0; fi< numF; fi++)
  74. {
  75. Eigen::SparseMatrix<typename DerivedV::Scalar > Sfi;
  76. assembleSelector(fi, Sfi);
  77. // the final matrix per face
  78. Eigen::SparseMatrix<typename DerivedV::Scalar > Qi = weightsSqrt(fi)*Ni*Sfi;
  79. // put it in the correct block of Q
  80. // todo: this can be made faster by omitting the selector matrix
  81. for (int k=0; k<Qi.outerSize(); ++k)
  82. for (typename Eigen::SparseMatrix<typename DerivedV::Scalar >::InnerIterator it(Qi,k); it; ++it)
  83. {
  84. typename DerivedV::Scalar val = it.value();
  85. int row = it.row();
  86. int col = it.col();
  87. tripletList.push_back(Eigen::Triplet<typename DerivedV::Scalar>(row+3*ni*fi,col,val));
  88. }
  89. }
  90. Q.resize(3*ni*numF,3*numV);
  91. Q.setFromTriplets(tripletList.begin(), tripletList.end());
  92. // the actual lhs matrix is Q'*Q
  93. // prefactor that matrix
  94. solver.compute(Q.transpose()*Q);
  95. if(solver.info()!=Eigen::Success)
  96. {
  97. std::cerr << "Cholesky failed - PlanarizerShapeUp.cpp" << std::endl;
  98. assert(0);
  99. }
  100. }
  101. template <typename DerivedV, typename DerivedF>
  102. inline void igl::PlanarizerShapeUp<DerivedV, DerivedF>::assembleNi()
  103. {
  104. std::vector<Eigen::Triplet<typename DerivedV::Scalar>> tripletList;
  105. for (int ii = 0; ii< ni; ii++)
  106. {
  107. for (int jj = 0; jj< ni; jj++)
  108. {
  109. tripletList.push_back(Eigen::Triplet<typename DerivedV::Scalar>(3*ii+0,3*jj+0,-1./ni));
  110. tripletList.push_back(Eigen::Triplet<typename DerivedV::Scalar>(3*ii+1,3*jj+1,-1./ni));
  111. tripletList.push_back(Eigen::Triplet<typename DerivedV::Scalar>(3*ii+2,3*jj+2,-1./ni));
  112. }
  113. tripletList.push_back(Eigen::Triplet<typename DerivedV::Scalar>(3*ii+0,3*ii+0,1.));
  114. tripletList.push_back(Eigen::Triplet<typename DerivedV::Scalar>(3*ii+1,3*ii+1,1.));
  115. tripletList.push_back(Eigen::Triplet<typename DerivedV::Scalar>(3*ii+2,3*ii+2,1.));
  116. }
  117. Ni.resize(3*ni,3*ni);
  118. Ni.setFromTriplets(tripletList.begin(), tripletList.end());
  119. }
  120. //assumes V stacked [x;y;z;x;y;z...];
  121. template <typename DerivedV, typename DerivedF>
  122. inline void igl::PlanarizerShapeUp<DerivedV, DerivedF>::assembleSelector(int fi,
  123. Eigen::SparseMatrix<typename DerivedV::Scalar > &S)
  124. {
  125. std::vector<Eigen::Triplet<typename DerivedV::Scalar>> tripletList;
  126. for (int fvi = 0; fvi< ni; fvi++)
  127. {
  128. int vi = Fin(fi,fvi);
  129. tripletList.push_back(Eigen::Triplet<typename DerivedV::Scalar>(3*fvi+0,3*vi+0,1.));
  130. tripletList.push_back(Eigen::Triplet<typename DerivedV::Scalar>(3*fvi+1,3*vi+1,1.));
  131. tripletList.push_back(Eigen::Triplet<typename DerivedV::Scalar>(3*fvi+2,3*vi+2,1.));
  132. }
  133. S.resize(3*ni,3*numV);
  134. S.setFromTriplets(tripletList.begin(), tripletList.end());
  135. }
  136. //project all faces to their closest planar face
  137. template <typename DerivedV, typename DerivedF>
  138. inline void igl::PlanarizerShapeUp<DerivedV, DerivedF>::assembleP()
  139. {
  140. P.setZero(3*ni*numF);
  141. for (int fi = 0; fi< numF; fi++)
  142. {
  143. // todo: this can be made faster by omitting the selector matrix
  144. Eigen::SparseMatrix<typename DerivedV::Scalar > Sfi;
  145. assembleSelector(fi, Sfi);
  146. Eigen::SparseMatrix<typename DerivedV::Scalar > NSi = Ni*Sfi;
  147. Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 1> Vi = NSi*Vv;
  148. Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, Eigen::Dynamic> CC(3,ni);
  149. for (int i = 0; i <ni; ++i)
  150. CC.col(i) = Vi.segment(3*i, 3);
  151. Eigen::Matrix<typename DerivedV::Scalar, 3, 3> C = CC*CC.transpose();
  152. // Alec: Doesn't compile
  153. Eigen::EigenSolver<Eigen::Matrix<typename DerivedV::Scalar, 3, 3>> es(C);
  154. // the real() is for compilation purposes
  155. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> lambda = es.eigenvalues().real();
  156. Eigen::Matrix<typename DerivedV::Scalar, 3, 3> U = es.eigenvectors().real();
  157. int min_i;
  158. lambda.cwiseAbs().minCoeff(&min_i);
  159. U.col(min_i).setZero();
  160. Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, Eigen::Dynamic> PP = U*U.transpose()*CC;
  161. for (int i = 0; i <ni; ++i)
  162. P.segment(3*ni*fi+3*i, 3) = weightsSqrt[fi]*PP.col(i);
  163. }
  164. }
  165. template <typename DerivedV, typename DerivedF>
  166. inline void igl::PlanarizerShapeUp<DerivedV, DerivedF>::planarize(Eigen::PlainObjectBase<DerivedV> &Vout)
  167. {
  168. Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 1> planarity;
  169. Vout = Vin;
  170. for (int iter =0; iter<maxIter; ++iter)
  171. {
  172. igl::quad_planarity(Vout, Fin, planarity);
  173. typename DerivedV::Scalar nonPlanarity = planarity.cwiseAbs().maxCoeff();
  174. std::cerr<<"iter #"<<iter<<": max non-planarity: "<<nonPlanarity<<std::endl;
  175. if (nonPlanarity<threshold)
  176. break;
  177. assembleP();
  178. Vv = solver.solve(Q.transpose()*P);
  179. if(solver.info()!=Eigen::Success)
  180. {
  181. std::cerr << "Linear solve failed - PlanarizerShapeUp.cpp" << std::endl;
  182. assert(0);
  183. }
  184. for (int i =0;i<numV;++i)
  185. Vout.row(i) << Vv.segment(3*i,3).transpose();
  186. }
  187. // set the mean of Vout to the mean of Vin
  188. Eigen::Matrix<typename DerivedV::Scalar, 1, 3> oldMean, newMean;
  189. oldMean = Vin.colwise().mean();
  190. newMean = Vout.colwise().mean();
  191. Vout.rowwise() += (oldMean - newMean);
  192. };
  193. template <typename DerivedV, typename DerivedF>
  194. IGL_INLINE void igl::planarize_quad_mesh(const Eigen::PlainObjectBase<DerivedV> &Vin,
  195. const Eigen::PlainObjectBase<DerivedF> &Fin,
  196. const int maxIter,
  197. const double &threshold,
  198. Eigen::PlainObjectBase<DerivedV> &Vout)
  199. {
  200. PlanarizerShapeUp<DerivedV, DerivedF> planarizer(Vin, Fin, maxIter, threshold);
  201. planarizer.planarize(Vout);
  202. }