cross_field_missmatch.cpp 5.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140
  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_missmatch.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. namespace igl {
  19. template <typename DerivedV, typename DerivedF>
  20. class MissMatchCalculator
  21. {
  22. public:
  23. const Eigen::PlainObjectBase<DerivedV> &V;
  24. const Eigen::PlainObjectBase<DerivedF> &F;
  25. const Eigen::PlainObjectBase<DerivedV> &PD1;
  26. const Eigen::PlainObjectBase<DerivedV> &PD2;
  27. #warning "Constructing Eigen::PlainObjectBase directly is deprecated"
  28. Eigen::PlainObjectBase<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. #warning "Constructing Eigen::PlainObjectBase directly is deprecated"
  35. Eigen::PlainObjectBase<DerivedF> TT;
  36. Eigen::PlainObjectBase<DerivedF> TTi;
  37. private:
  38. ///compute the mismatch between 2 faces
  39. inline int MissMatchByCross(const int f0,
  40. const int f1)
  41. {
  42. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> dir0 = PD1.row(f0);
  43. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> dir1 = PD1.row(f1);
  44. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> n0 = N.row(f0);
  45. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> n1 = N.row(f1);
  46. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> dir1Rot = igl::rotation_matrix_from_directions(n1,n0)*dir1;
  47. dir1Rot.normalize();
  48. // TODO: this should be equivalent to the other code below, to check!
  49. // Compute the angle between the two vectors
  50. // double a0 = atan2(dir0.dot(B2.row(f0)),dir0.dot(B1.row(f0)));
  51. // double a1 = atan2(dir1Rot.dot(B2.row(f0)),dir1Rot.dot(B1.row(f0)));
  52. //
  53. // double angle_diff = a1-a0; //VectToAngle(f0,dir1Rot);
  54. double angle_diff = atan2(dir1Rot.dot(PD2.row(f0)),dir1Rot.dot(PD1.row(f0)));
  55. // std::cerr << "Dani: " << dir0(0) << " " << dir0(1) << " " << dir0(2) << " " << dir1Rot(0) << " " << dir1Rot(1) << " " << dir1Rot(2) << " " << angle_diff << std::endl;
  56. double step=M_PI/2.0;
  57. int i=(int)std::floor((angle_diff/step)+0.5);
  58. int k=0;
  59. if (i>=0)
  60. k=i%4;
  61. else
  62. k=(-(3*i))%4;
  63. return k;
  64. }
  65. public:
  66. inline MissMatchCalculator(const Eigen::PlainObjectBase<DerivedV> &_V,
  67. const Eigen::PlainObjectBase<DerivedF> &_F,
  68. const Eigen::PlainObjectBase<DerivedV> &_PD1,
  69. const Eigen::PlainObjectBase<DerivedV> &_PD2
  70. ):
  71. V(_V),
  72. F(_F),
  73. PD1(_PD1),
  74. PD2(_PD2)
  75. {
  76. igl::per_face_normals(V,F,N);
  77. V_border = igl::is_border_vertex(V,F);
  78. igl::vertex_triangle_adjacency(V,F,VF,VFi);
  79. igl::triangle_triangle_adjacency(F,TT,TTi);
  80. }
  81. inline void calculateMissmatch(Eigen::PlainObjectBase<DerivedF> &Handle_MMatch)
  82. {
  83. Handle_MMatch.setConstant(F.rows(),3,-1);
  84. for (size_t i=0;i<F.rows();i++)
  85. {
  86. for (int j=0;j<3;j++)
  87. {
  88. if (((int)i)==TT(i,j) || TT(i,j) == -1)
  89. Handle_MMatch(i,j)=0;
  90. else
  91. Handle_MMatch(i,j) = MissMatchByCross(i,TT(i,j));
  92. }
  93. }
  94. }
  95. };
  96. }
  97. template <typename DerivedV, typename DerivedF>
  98. IGL_INLINE void igl::cross_field_missmatch(const Eigen::PlainObjectBase<DerivedV> &V,
  99. const Eigen::PlainObjectBase<DerivedF> &F,
  100. const Eigen::PlainObjectBase<DerivedV> &PD1,
  101. const Eigen::PlainObjectBase<DerivedV> &PD2,
  102. const bool isCombed,
  103. Eigen::PlainObjectBase<DerivedF> &missmatch)
  104. {
  105. DerivedV PD1_combed;
  106. DerivedV PD2_combed;
  107. if (!isCombed)
  108. igl::comb_cross_field(V,F,PD1,PD2,PD1_combed,PD2_combed);
  109. else
  110. {
  111. PD1_combed = PD1;
  112. PD2_combed = PD2;
  113. }
  114. igl::MissMatchCalculator<DerivedV, DerivedF> sf(V, F, PD1_combed, PD2_combed);
  115. sf.calculateMissmatch(missmatch);
  116. }
  117. #ifdef IGL_STATIC_LIBRARY
  118. // Explicit template specialization
  119. template void igl::cross_field_missmatch<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&, bool, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> >&);
  120. template void igl::cross_field_missmatch<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&, bool, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  121. #endif