planarize_quad_mesh.cpp 8.4 KB

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