per_vertex_point_to_plane_quadrics.cpp 5.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151
  1. #include "per_vertex_point_to_plane_quadrics.h"
  2. #include "quadric_binary_plus_operator.h"
  3. #include <Eigen/QR>
  4. #include <cassert>
  5. #include <cmath>
  6. IGL_INLINE void igl::per_vertex_point_to_plane_quadrics(
  7. const Eigen::MatrixXd & V,
  8. const Eigen::MatrixXi & F,
  9. const Eigen::MatrixXi & EMAP,
  10. const Eigen::MatrixXi & EF,
  11. const Eigen::MatrixXi & EI,
  12. std::vector<
  13. std::tuple<Eigen::MatrixXd,Eigen::RowVectorXd,double> > & quadrics)
  14. {
  15. using namespace std;
  16. typedef std::tuple<Eigen::MatrixXd,Eigen::RowVectorXd,double> Quadric;
  17. const int dim = V.cols();
  18. //// Quadrics per face
  19. //std::vector<Quadric> face_quadrics(F.rows());
  20. // Initialize each vertex quadric to zeros
  21. quadrics.resize(
  22. V.rows(),
  23. // gcc <=4.8 can't handle initializer lists correctly
  24. (Quadric)
  25. {Eigen::MatrixXd::Zero(dim,dim),Eigen::RowVectorXd::Zero(dim),0});
  26. Eigen::MatrixXd I = Eigen::MatrixXd::Identity(dim,dim);
  27. // Rather initial with zeros, initial with a small amount of energy pull
  28. // toward original vertex position
  29. const double w = 1e-10;
  30. for(int v = 0;v<V.rows();v++)
  31. {
  32. std::get<0>(quadrics[v]) = w*I;
  33. Eigen::RowVectorXd Vv = V.row(v);
  34. std::get<1>(quadrics[v]) = w*-Vv;
  35. std::get<2>(quadrics[v]) = w*Vv.dot(Vv);
  36. }
  37. // Generic nD qslim from "Simplifying Surfaces with Color and Texture
  38. // using Quadric Error Metric" (follow up to original QSlim)
  39. for(int f = 0;f<F.rows();f++)
  40. {
  41. int infinite_corner = -1;
  42. for(int c = 0;c<3;c++)
  43. {
  44. if(
  45. std::isinf(V(F(f,c),0)) ||
  46. std::isinf(V(F(f,c),1)) ||
  47. std::isinf(V(F(f,c),2)))
  48. {
  49. assert(infinite_corner == -1 && "Should only be one infinite corner");
  50. infinite_corner = c;
  51. }
  52. }
  53. // Inputs:
  54. // p 1 by n row point on the subspace
  55. // S m by n matrix where rows coorespond to orthonormal spanning
  56. // vectors of the subspace to which we're measuring distance (usually
  57. // a plane, m=2)
  58. // weight scalar weight
  59. // Returns quadric triple {A,b,c} so that A-2*b+c measures the quadric
  60. const auto subspace_quadric = [&I](
  61. const Eigen::RowVectorXd & p,
  62. const Eigen::MatrixXd & S,
  63. const double weight)->Quadric
  64. {
  65. // Dimension of subspace
  66. const int m = S.rows();
  67. // Weight face's quadric (v'*A*v + 2*b'*v + c) by area
  68. // e1 and e2 should be perpendicular
  69. Eigen::MatrixXd A = I;
  70. Eigen::RowVectorXd b = -p;
  71. double c = p.dot(p);
  72. for(int i = 0;i<m;i++)
  73. {
  74. Eigen::RowVectorXd ei = S.row(i);
  75. for(int j = 0;j<i;j++) assert(std::abs(S.row(j).dot(ei)) < 1e-10);
  76. A += -ei.transpose()*ei;
  77. b += p.dot(ei)*ei;
  78. c += -pow(p.dot(ei),2);
  79. }
  80. // gcc <=4.8 can't handle initializer lists correctly: needs explicit
  81. // cast
  82. return (Quadric){ weight*A, weight*b, weight*c };
  83. };
  84. if(infinite_corner == -1)
  85. {
  86. // Finite (non-boundary) face
  87. Eigen::RowVectorXd p = V.row(F(f,0));
  88. Eigen::RowVectorXd q = V.row(F(f,1));
  89. Eigen::RowVectorXd r = V.row(F(f,2));
  90. Eigen::RowVectorXd pq = q-p;
  91. Eigen::RowVectorXd pr = r-p;
  92. // Gram Determinant = squared area of parallelogram
  93. double area = sqrt(pq.squaredNorm()*pr.squaredNorm()-pow(pr.dot(pq),2));
  94. Eigen::RowVectorXd e1 = pq.normalized();
  95. Eigen::RowVectorXd e2 = (pr-e1.dot(pr)*e1).normalized();
  96. Eigen::MatrixXd S(2,V.cols());
  97. S<<e1,e2;
  98. Quadric face_quadric = subspace_quadric(p,S,area);
  99. // Throw at each corner
  100. for(int c = 0;c<3;c++)
  101. {
  102. quadrics[F(f,c)] = quadrics[F(f,c)] + face_quadric;
  103. }
  104. }else
  105. {
  106. // cth corner is infinite --> edge opposite cth corner is boundary
  107. // Boundary edge vector
  108. const Eigen::RowVectorXd p = V.row(F(f,(infinite_corner+1)%3));
  109. Eigen::RowVectorXd ev = V.row(F(f,(infinite_corner+2)%3)) - p;
  110. const double length = ev.norm();
  111. ev /= length;
  112. // Face neighbor across boundary edge
  113. int e = EMAP(f+F.rows()*infinite_corner);
  114. int opp = EF(e,0) == f ? 1 : 0;
  115. int n = EF(e,opp);
  116. int nc = EI(e,opp);
  117. assert(
  118. ((F(f,(infinite_corner+1)%3) == F(n,(nc+1)%3) &&
  119. F(f,(infinite_corner+2)%3) == F(n,(nc+2)%3)) ||
  120. (F(f,(infinite_corner+1)%3) == F(n,(nc+2)%3)
  121. && F(f,(infinite_corner+2)%3) == F(n,(nc+1)%3))) &&
  122. "Edge flaps not agreeing on shared edge");
  123. // Edge vector on opposite face
  124. const Eigen::RowVectorXd eu = V.row(F(n,nc)) - p;
  125. assert(!isinf(eu(0)));
  126. // Matrix with vectors spanning plane as columns
  127. Eigen::MatrixXd A(ev.size(),2);
  128. A<<ev.transpose(),eu.transpose();
  129. // Use QR decomposition to find basis for orthogonal space
  130. Eigen::HouseholderQR<Eigen::MatrixXd> qr(A);
  131. const Eigen::MatrixXd Q = qr.householderQ();
  132. const Eigen::MatrixXd N =
  133. Q.topRightCorner(ev.size(),ev.size()-2).transpose();
  134. assert(N.cols() == ev.size());
  135. assert(N.rows() == ev.size()-2);
  136. Eigen::MatrixXd S(N.rows()+1,ev.size());
  137. S<<ev,N;
  138. Quadric boundary_edge_quadric = subspace_quadric(p,S,length);
  139. for(int c = 0;c<3;c++)
  140. {
  141. if(c != infinite_corner)
  142. {
  143. quadrics[F(f,c)] = quadrics[F(f,c)] + boundary_edge_quadric;
  144. }
  145. }
  146. }
  147. }
  148. }