mvc.cpp 5.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189
  1. #include "mvc.h"
  2. #include <vector>
  3. #include <cassert>
  4. #include <iostream>
  5. // Broken Implementation
  6. IGL_INLINE void igl::mvc(const Eigen::MatrixXd &V, const Eigen::MatrixXd &C, Eigen::MatrixXd &W)
  7. {
  8. // at least three control points
  9. assert(C.rows()>2);
  10. // dimension of points
  11. assert(C.cols() == 3 || C.cols() == 2);
  12. assert(V.cols() == 3 || V.cols() == 2);
  13. // number of polygon points
  14. int num = C.rows();
  15. Eigen::MatrixXd V1,C1;
  16. int i_prev, i_next;
  17. // check if either are 3D but really all z's are 0
  18. bool V_flat = (V.cols() == 3) && (std::sqrt( (V.col(3)).dot(V.col(3)) ) < 1e-10);
  19. bool C_flat = (C.cols() == 3) && (std::sqrt( (C.col(3)).dot(C.col(3)) ) < 1e-10);
  20. // if both are essentially 2D then ignore z-coords
  21. if((C.cols() == 2 || C_flat) && (V.cols() == 2 || V_flat))
  22. {
  23. // ignore z coordinate
  24. V1 = V.block(0,0,V.rows(),2);
  25. C1 = C.block(0,0,C.rows(),2);
  26. }
  27. else
  28. {
  29. // give dummy z coordinate to either mesh or poly
  30. if(V.rows() == 2)
  31. {
  32. V1 = Eigen::MatrixXd(V.rows(),3);
  33. V1.block(0,0,V.rows(),2) = V;
  34. }
  35. else
  36. V1 = V;
  37. if(C.rows() == 2)
  38. {
  39. C1 = Eigen::MatrixXd(C.rows(),3);
  40. C1.block(0,0,C.rows(),2) = C;
  41. }
  42. else
  43. C1 = C;
  44. // check that C is planar
  45. // average normal around poly corners
  46. Eigen::Vector3d n = Eigen::Vector3d::Zero();
  47. // take centroid as point on plane
  48. Eigen::Vector3d p = Eigen::Vector3d::Zero();
  49. for (int i = 0; i<num; ++i)
  50. {
  51. i_prev = (i>0)?(i-1):(num-1);
  52. i_next = (i<num-1)?(i+1):0;
  53. Eigen::Vector3d vnext = (C1.row(i_next) - C1.row(i)).transpose();
  54. Eigen::Vector3d vprev = (C1.row(i_prev) - C1.row(i)).transpose();
  55. n += vnext.cross(vprev);
  56. p += C1.row(i);
  57. }
  58. p/=num;
  59. n/=num;
  60. // normalize n
  61. n /= std::sqrt(n.dot(n));
  62. // check that poly is really coplanar
  63. #ifndef NDEBUG
  64. for (int i = 0; i<num; ++i)
  65. {
  66. double dist_to_plane_C = std::abs((C1.row(i)-p.transpose()).dot(n));
  67. assert(dist_to_plane_C<1e-10);
  68. }
  69. #endif
  70. // check that poly is really coplanar
  71. for (int i = 0; i<V1.rows(); ++i)
  72. {
  73. double dist_to_plane_V = std::abs((V1.row(i)-p.transpose()).dot(n));
  74. if(dist_to_plane_V>1e-10)
  75. std::cerr<<"Distance from V to plane of C is large..."<<std::endl;
  76. }
  77. // change of basis
  78. Eigen::Vector3d b1 = C1.row(1)-C1.row(0);
  79. Eigen::Vector3d b2 = n.cross(b1);
  80. // normalize basis rows
  81. b1 /= std::sqrt(b1.dot(b1));
  82. b2 /= std::sqrt(b2.dot(b2));
  83. n /= std::sqrt(n.dot(n));
  84. //transpose of the basis matrix in the m-file
  85. Eigen::Matrix3d basis = Eigen::Matrix3d::Zero();
  86. basis.col(0) = b1;
  87. basis.col(1) = b2;
  88. basis.col(2) = n;
  89. // change basis of rows vectors by right multiplying with inverse of matrix
  90. // with basis vectors as rows
  91. Eigen::ColPivHouseholderQR<Eigen::Matrix3d> solver = basis.colPivHouseholderQr();
  92. // Throw away coordinates in normal direction
  93. V1 = solver.solve(V1.transpose()).transpose().block(0,0,V1.rows(),2);
  94. C1 = solver.solve(C1.transpose()).transpose().block(0,0,C1.rows(),2);
  95. }
  96. // vectors from V to every C, where CmV(i,j,:) is the vector from domain
  97. // vertex j to handle i
  98. double EPSILON = 1e-10;
  99. Eigen::MatrixXd WW = Eigen::MatrixXd(C1.rows(), V1.rows());
  100. Eigen::MatrixXd dist_C_V (C1.rows(), V1.rows());
  101. std::vector< std::pair<int,int> > on_corner(0);
  102. std::vector< std::pair<int,int> > on_segment(0);
  103. for (int i = 0; i<C1.rows(); ++i)
  104. {
  105. i_prev = (i>0)?(i-1):(num-1);
  106. i_next = (i<num-1)?(i+1):0;
  107. // distance from each corner in C to the next corner so that edge_length(i)
  108. // is the distance from C(i,:) to C(i+1,:) defined cyclically
  109. double edge_length = std::sqrt((C1.row(i) - C1.row(i_next)).dot(C1.row(i) - C1.row(i_next)));
  110. for (int j = 0; j<V1.rows(); ++j)
  111. {
  112. Eigen::VectorXd v = C1.row(i) - V1.row(j);
  113. Eigen::VectorXd vnext = C1.row(i_next) - V1.row(j);
  114. Eigen::VectorXd vprev = C1.row(i_prev) - V1.row(j);
  115. // distance from V to every C, where dist_C_V(i,j) is the distance from domain
  116. // vertex j to handle i
  117. dist_C_V(i,j) = std::sqrt(v.dot(v));
  118. double dist_C_V_next = std::sqrt(vnext.dot(vnext));
  119. double a_prev = std::atan2(vprev[1],vprev[0]) - std::atan2(v[1],v[0]);
  120. double a_next = std::atan2(v[1],v[0]) - std::atan2(vnext[1],vnext[0]);
  121. // mean value coordinates
  122. WW(i,j) = (std::tan(a_prev/2.0) + std::tan(a_next/2.0)) / dist_C_V(i,j);
  123. if (dist_C_V(i,j) < EPSILON)
  124. on_corner.push_back(std::make_pair(j,i));
  125. else
  126. // only in case of no-corner (no need for checking for multiple segments afterwards --
  127. // should only be on one segment (otherwise must be on a corner and we already
  128. // handled that)
  129. // domain vertex j is on the segment from i to i+1 if the distances from vj to
  130. // pi and pi+1 are about
  131. if(abs((dist_C_V(i,j) + dist_C_V_next) / edge_length - 1) < EPSILON)
  132. on_segment.push_back(std::make_pair(j,i));
  133. }
  134. }
  135. // handle degenerate cases
  136. // snap vertices close to corners
  137. for (unsigned i = 0; i<on_corner.size(); ++i)
  138. {
  139. int vi = on_corner[i].first;
  140. int ci = on_corner[i].second;
  141. for (int ii = 0; ii<C.rows(); ++ii)
  142. WW(ii,vi) = (ii==ci)?1:0;
  143. }
  144. // snap vertices close to segments
  145. for (unsigned i = 0; i<on_segment.size(); ++i)
  146. {
  147. int vi = on_segment[i].first;
  148. int ci = on_segment[i].second;
  149. int ci_next = (ci<num-1)?(ci+1):0;
  150. for (int ii = 0; ii<C.rows(); ++ii)
  151. if (ii == ci)
  152. WW(ii,vi) = dist_C_V(ci_next,vi);
  153. else
  154. {
  155. if ( ii == ci_next)
  156. WW(ii,vi) = dist_C_V(ci,vi);
  157. else
  158. WW(ii,vi) = 0;
  159. }
  160. }
  161. // normalize W
  162. for (int i = 0; i<V.rows(); ++i)
  163. WW.col(i) /= WW.col(i).sum();
  164. // we've made W transpose
  165. W = WW.transpose();
  166. }