cross_field_mismatch.cpp 5.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2014 Daniele Panozzo <daniele.panozzo@gmail.com>, Olga Diamanti <olga.diam@gmail.com>
  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 "cross_field_mismatch.h"
  9. #include <cmath>
  10. #include <vector>
  11. #include <deque>
  12. #include <igl/comb_cross_field.h>
  13. #include <igl/per_face_normals.h>
  14. #include <igl/is_border_vertex.h>
  15. #include <igl/vertex_triangle_adjacency.h>
  16. #include <igl/triangle_triangle_adjacency.h>
  17. #include <igl/rotation_matrix_from_directions.h>
  18. #include <igl/PI.h>
  19. namespace igl {
  20. template <typename DerivedV, typename DerivedF, typename DerivedM>
  21. class MismatchCalculator
  22. {
  23. public:
  24. const Eigen::PlainObjectBase<DerivedV> &V;
  25. const Eigen::PlainObjectBase<DerivedF> &F;
  26. const Eigen::PlainObjectBase<DerivedV> &PD1;
  27. const Eigen::PlainObjectBase<DerivedV> &PD2;
  28. DerivedV N;
  29. private:
  30. // internal
  31. std::vector<bool> V_border; // bool
  32. std::vector<std::vector<int> > VF;
  33. std::vector<std::vector<int> > VFi;
  34. DerivedF TT;
  35. DerivedF TTi;
  36. private:
  37. ///compute the mismatch between 2 faces
  38. inline int mismatchByCross(const int f0,
  39. const int f1)
  40. {
  41. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> dir0 = PD1.row(f0);
  42. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> dir1 = PD1.row(f1);
  43. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> n0 = N.row(f0);
  44. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> n1 = N.row(f1);
  45. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> dir1Rot = igl::rotation_matrix_from_directions(n1,n0)*dir1;
  46. dir1Rot.normalize();
  47. double angle_diff = atan2(dir1Rot.dot(PD2.row(f0)),dir1Rot.dot(PD1.row(f0)));
  48. double step=igl::PI/2.0;
  49. int i=(int)std::floor((angle_diff/step)+0.5);
  50. int k=0;
  51. if (i>=0)
  52. k=i%4;
  53. else
  54. k=(-(3*i))%4;
  55. return k;
  56. }
  57. public:
  58. inline MismatchCalculator(const Eigen::PlainObjectBase<DerivedV> &_V,
  59. const Eigen::PlainObjectBase<DerivedF> &_F,
  60. const Eigen::PlainObjectBase<DerivedV> &_PD1,
  61. const Eigen::PlainObjectBase<DerivedV> &_PD2):
  62. V(_V),
  63. F(_F),
  64. PD1(_PD1),
  65. PD2(_PD2)
  66. {
  67. igl::per_face_normals(V,F,N);
  68. V_border = igl::is_border_vertex(V,F);
  69. igl::vertex_triangle_adjacency(V,F,VF,VFi);
  70. igl::triangle_triangle_adjacency(F,TT,TTi);
  71. }
  72. inline void calculateMismatch(Eigen::PlainObjectBase<DerivedM> &Handle_MMatch)
  73. {
  74. Handle_MMatch.setConstant(F.rows(),3,-1);
  75. for (size_t i=0;i<F.rows();i++)
  76. {
  77. for (int j=0;j<3;j++)
  78. {
  79. if (((int)i)==TT(i,j) || TT(i,j) == -1)
  80. Handle_MMatch(i,j)=0;
  81. else
  82. Handle_MMatch(i,j) = mismatchByCross(i, TT(i, j));
  83. }
  84. }
  85. }
  86. };
  87. }
  88. template <typename DerivedV, typename DerivedF, typename DerivedM>
  89. IGL_INLINE void igl::cross_field_mismatch(const Eigen::PlainObjectBase<DerivedV> &V,
  90. const Eigen::PlainObjectBase<DerivedF> &F,
  91. const Eigen::PlainObjectBase<DerivedV> &PD1,
  92. const Eigen::PlainObjectBase<DerivedV> &PD2,
  93. const bool isCombed,
  94. Eigen::PlainObjectBase<DerivedM> &mismatch)
  95. {
  96. DerivedV PD1_combed;
  97. DerivedV PD2_combed;
  98. if (!isCombed)
  99. igl::comb_cross_field(V,F,PD1,PD2,PD1_combed,PD2_combed);
  100. else
  101. {
  102. PD1_combed = PD1;
  103. PD2_combed = PD2;
  104. }
  105. igl::MismatchCalculator<DerivedV, DerivedF, DerivedM> sf(V, F, PD1_combed, PD2_combed);
  106. sf.calculateMismatch(mismatch);
  107. }
  108. #ifdef IGL_STATIC_LIBRARY
  109. // Explicit template instantiation
  110. template void igl::cross_field_mismatch<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const &, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const &, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const &, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const &, const bool, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > &);
  111. template void igl::cross_field_mismatch<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >( Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const &, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const &, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const &, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const &, const bool, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > &);
  112. template void igl::cross_field_mismatch<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 3, 0, -1, 3> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const &, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const &, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const &, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const &, const bool, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > &);
  113. #endif