polyvector_field_matchings.cpp 9.6 KB

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