polyvector_field_matchings.cpp 10 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197
  1. #include <igl/polyvector_field_matchings.h>
  2. #include <igl/edge_topology.h>
  3. #include <igl/per_face_normals.h>
  4. #include <igl/sort_vectors_ccw.h>
  5. #include <igl/rotation_matrix_from_directions.h>
  6. template <typename DerivedS, typename DerivedV, typename DerivedM>
  7. IGL_INLINE void igl::polyvector_field_matching(
  8. const Eigen::PlainObjectBase<DerivedS>& _ua,
  9. const Eigen::PlainObjectBase<DerivedS>& _ub,
  10. const Eigen::PlainObjectBase<DerivedV>& na,
  11. const Eigen::PlainObjectBase<DerivedV>& nb,
  12. const Eigen::PlainObjectBase<DerivedV>& e,
  13. bool match_with_curl,
  14. Eigen::PlainObjectBase<DerivedM>& mab,
  15. Eigen::PlainObjectBase<DerivedM>& mba,
  16. bool is_symmetric)
  17. {
  18. // make sure the matching preserve ccw order of the vectors across the edge
  19. // 1) order vectors in a, ccw e.g. (0,1,2,3)_a not ccw --> (0,3,2,1)_a ccw
  20. // 2) order vectors in b, ccw e.g. (0,1,2,3)_b not ccw --> (0,2,1,3)_b ccw
  21. // 3) the vectors in b that match the ordered vectors in a (in this case (0,3,2,1)_a ) must be a circular shift of the ccw ordered vectors in b - so we need to explicitely check only these matchings to find the best ccw one, there are N of them
  22. int hN = _ua.cols()/3;
  23. int N = hN;
  24. Eigen::Matrix<typename DerivedS::Scalar,1,Eigen::Dynamic> ua (1,N*3);
  25. Eigen::Matrix<typename DerivedS::Scalar,1,Eigen::Dynamic> ub (1,N*3);
  26. if (is_symmetric)
  27. {
  28. ua <<_ua, -_ua;
  29. ub <<_ub, -_ub;
  30. N = 2*hN;
  31. }
  32. else
  33. {
  34. ua =_ua;
  35. ub =_ub;
  36. }
  37. Eigen::Matrix<typename DerivedM::Scalar,Eigen::Dynamic,1> order_a, order_b;
  38. igl::sort_vectors_ccw(ua, na, order_a);
  39. igl::sort_vectors_ccw(ub, nb, order_b);
  40. //checking all possible circshifts of order_b as matches for order_a
  41. Eigen::Matrix<typename DerivedM::Scalar,Eigen::Dynamic,Eigen::Dynamic> all_matches(N,N);
  42. Eigen::Matrix<typename DerivedS::Scalar,1,Eigen::Dynamic> all_scores(1,N);
  43. for (int s =0; s< N; ++s)
  44. {
  45. all_matches.row(s) = order_b;
  46. typename DerivedS::Scalar current_score=0;
  47. for (int i=0; i< N; ++i)
  48. {
  49. if (match_with_curl)
  50. current_score += fabs(ua.segment(order_a[i]*3, 3).dot(e) - ub.segment(order_b[i]*3, 3).dot(e));
  51. else
  52. {
  53. Eigen::Matrix<typename DerivedS::Scalar,3,1> na_ = na.transpose();
  54. Eigen::Matrix<typename DerivedS::Scalar,3,1> nb_ = nb.transpose();
  55. Eigen::Matrix<typename DerivedS::Scalar,1,3> uaRot = igl::rotation_matrix_from_directions(na_, nb_)*ua.segment(order_a[i]*3, 3).transpose();
  56. current_score += (uaRot-ub.segment(order_b[i]*3, 3)).norm();
  57. }
  58. }
  59. all_scores[s] = current_score;
  60. // do a circshift on order_b to check the next preserving matching
  61. int temp = order_b[0];
  62. for (int i =0; i< N-1; ++i)
  63. order_b[i] = order_b[i+1];
  64. order_b(N-1) = temp;
  65. }
  66. Eigen::Matrix<typename DerivedM::Scalar,1,Eigen::Dynamic> best_matching_for_sorted_a;
  67. int best_i;
  68. all_scores.minCoeff(&best_i);
  69. best_matching_for_sorted_a = all_matches.row(best_i);
  70. // best_matching_for_sorted_a is the matching for the sorted vectors in a
  71. // get the matching for the initial (unsorted) vectors
  72. mab.resize(1,N);
  73. for (int i=0; i< N; ++i)
  74. mab[order_a[i]] = best_matching_for_sorted_a[i];
  75. //mab contains the best matching a->b, get the opposite too
  76. mba.resize(1, N);
  77. for (int i=0; i< N; ++i)
  78. mba[mab[i]] = i;
  79. if (is_symmetric)
  80. {
  81. mab = mab.head(hN);
  82. mba = mba.head(hN);
  83. }
  84. }
  85. template <typename DerivedS, typename DerivedV, typename DerivedF, typename DerivedE, typename DerivedM, typename DerivedC>
  86. IGL_INLINE typename DerivedC::Scalar igl::polyvector_field_matchings(
  87. const Eigen::PlainObjectBase<DerivedS>& sol3D,
  88. const Eigen::PlainObjectBase<DerivedV>&V,
  89. const Eigen::PlainObjectBase<DerivedF>&F,
  90. const Eigen::PlainObjectBase<DerivedE>&E,
  91. const Eigen::PlainObjectBase<DerivedV>& FN,
  92. const Eigen::MatrixXi &E2F,
  93. bool match_with_curl,
  94. Eigen::PlainObjectBase<DerivedM>& match_ab,
  95. Eigen::PlainObjectBase<DerivedM>& match_ba,
  96. Eigen::PlainObjectBase<DerivedC>& curl)
  97. {
  98. int numEdges = E.rows();
  99. int half_degree = sol3D.cols()/3;
  100. Eigen::VectorXi isBorderEdge;
  101. isBorderEdge.setZero(numEdges,1);
  102. for(unsigned i=0; i<numEdges; ++i)
  103. {
  104. if ((E2F(i,0) == -1) || ((E2F(i,1) == -1)))
  105. isBorderEdge[i] = 1;
  106. }
  107. curl.setZero(numEdges,1);
  108. match_ab.setZero(numEdges, half_degree);
  109. match_ba.setZero(numEdges, half_degree);
  110. typename DerivedC::Scalar meanCurl = 0;
  111. for (int ei=0; ei<numEdges; ++ei)
  112. {
  113. if (isBorderEdge[ei])
  114. continue;
  115. // the two faces of the flap
  116. int a = E2F(ei,0);
  117. int b = E2F(ei,1);
  118. Eigen::Matrix<typename DerivedV::Scalar, 1, Eigen::Dynamic> ce = (V.row(E(ei,1))-V.row(E(ei,0))).normalized().template cast<typename DerivedV::Scalar>();
  119. Eigen::Matrix<typename DerivedM::Scalar, 1, Eigen::Dynamic> mab, mba;
  120. igl::polyvector_field_matching(sol3D.row(a).eval(),
  121. sol3D.row(b).eval(),
  122. FN.row(a).eval(),
  123. FN.row(b).eval(),
  124. ce,
  125. match_with_curl,
  126. mab,
  127. mba,
  128. true);
  129. match_ab.row(ei) = mab;
  130. match_ba.row(ei) = mba;
  131. Eigen::Matrix<typename DerivedS::Scalar, 1, Eigen::Dynamic> matched;
  132. matched.resize(1, 3*half_degree);
  133. for (int i = 0; i<half_degree; ++i)
  134. {
  135. int sign = (mab[i]<half_degree)?1:-1;
  136. matched.segment(i*3, 3) = sign*sol3D.row(b).segment((mab[i]%half_degree)*3, 3);
  137. }
  138. typename DerivedC::Scalar avgCurl = 0;
  139. for (int i = 0; i<half_degree; ++i)
  140. avgCurl += fabs(sol3D.row(a).segment(i*3, 3).dot(ce) - matched.segment(i*3, 3).dot(ce));
  141. avgCurl = avgCurl/half_degree;
  142. curl[ ei ] = avgCurl;
  143. meanCurl+= avgCurl;
  144. }
  145. meanCurl /= 1.*(numEdges - isBorderEdge.sum());
  146. return meanCurl;
  147. }
  148. template <typename DerivedS, typename DerivedV, typename DerivedF, typename DerivedM, typename DerivedC>
  149. IGL_INLINE typename DerivedC::Scalar igl::polyvector_field_matchings(
  150. const Eigen::PlainObjectBase<DerivedS>& sol3D,
  151. const Eigen::PlainObjectBase<DerivedV>&V,
  152. const Eigen::PlainObjectBase<DerivedF>&F,
  153. bool match_with_curl,
  154. Eigen::PlainObjectBase<DerivedM>& match_ab,
  155. Eigen::PlainObjectBase<DerivedM>& match_ba,
  156. Eigen::PlainObjectBase<DerivedC>& curl)
  157. {
  158. Eigen::MatrixXi E, E2F, F2E;
  159. igl::edge_topology(V,F,E,F2E,E2F);
  160. Eigen::PlainObjectBase<DerivedV> FN;
  161. igl::per_face_normals(V,F,FN);
  162. return igl::polyvector_field_matchings(sol3D, V, F, E, FN, E2F, match_with_curl, match_ab, match_ba, curl);
  163. }
  164. #ifdef IGL_STATIC_LIBRARY
  165. // Explicit template specialization
  166. template Eigen::Matrix<float, -1, 1, 0, -1, 1>::Scalar igl::polyvector_field_matchings<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<float, -1, 1, 0, -1, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, bool, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<float, -1, 1, 0, -1, 1> >&);
  167. template Eigen::Matrix<double, -1, 1, 0, -1, 1>::Scalar igl::polyvector_field_matchings<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, bool, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  168. template void igl::polyvector_field_matching<Eigen::Matrix<double, 1, -1, 1, 1, -1>, Eigen::Matrix<double, 1, -1, 1, 1, -1>, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, 1, -1, 1, 1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, 1, -1, 1, 1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, 1, -1, 1, 1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, 1, -1, 1, 1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, 1, -1, 1, 1, -1> > const&, bool, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&, bool);
  169. #endif