mvc.cpp 6.2 KB

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