planarize_quad_mesh.cpp 8.5 KB

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