bfs_orient.cpp 2.9 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798
  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 "bfs_orient.h"
  9. #include "orientable_patches.h"
  10. #include "parallel_for.h"
  11. #include <Eigen/Sparse>
  12. #include <queue>
  13. template <typename DerivedF, typename DerivedFF, typename DerivedC>
  14. IGL_INLINE void igl::bfs_orient(
  15. const Eigen::MatrixBase<DerivedF> & F,
  16. Eigen::PlainObjectBase<DerivedFF> & FF,
  17. Eigen::PlainObjectBase<DerivedC> & C)
  18. {
  19. Eigen::SparseMatrix<typename DerivedF::Scalar> A;
  20. orientable_patches(F,C,A);
  21. // number of faces
  22. const int m = F.rows();
  23. // number of patches
  24. const int num_cc = C.maxCoeff()+1;
  25. Eigen::VectorXi seen = Eigen::VectorXi::Zero(m);
  26. // Edge sets
  27. const int ES[3][2] = {{1,2},{2,0},{0,1}};
  28. if(((void*)&FF) != ((void*)&F))
  29. {
  30. FF = F;
  31. }
  32. // loop over patches
  33. parallel_for(num_cc,[&](const int c)
  34. {
  35. std::queue<typename DerivedF::Scalar> Q;
  36. // find first member of patch c
  37. for(int f = 0;f<FF.rows();f++)
  38. {
  39. if(C(f) == c)
  40. {
  41. Q.push(f);
  42. break;
  43. }
  44. }
  45. assert(!Q.empty());
  46. while(!Q.empty())
  47. {
  48. const typename DerivedF::Scalar f = Q.front();
  49. Q.pop();
  50. if(seen(f) > 0)
  51. {
  52. continue;
  53. }
  54. seen(f)++;
  55. // loop over neighbors of f
  56. for(typename Eigen::SparseMatrix<typename DerivedF::Scalar>::InnerIterator it (A,f); it; ++it)
  57. {
  58. // might be some lingering zeros, and skip self-adjacency
  59. if(it.value() != 0 && it.row() != f)
  60. {
  61. const int n = it.row();
  62. assert(n != f);
  63. // loop over edges of f
  64. for(int efi = 0;efi<3;efi++)
  65. {
  66. // efi'th edge of face f
  67. Eigen::Vector2i ef(FF(f,ES[efi][0]),FF(f,ES[efi][1]));
  68. // loop over edges of n
  69. for(int eni = 0;eni<3;eni++)
  70. {
  71. // eni'th edge of face n
  72. Eigen::Vector2i en(FF(n,ES[eni][0]),FF(n,ES[eni][1]));
  73. // Match (half-edges go same direction)
  74. if(ef(0) == en(0) && ef(1) == en(1))
  75. {
  76. // flip face n
  77. FF.row(n) = FF.row(n).reverse().eval();
  78. }
  79. }
  80. }
  81. // add neighbor to queue
  82. Q.push(n);
  83. }
  84. }
  85. }
  86. },1000);
  87. // make sure flip is OK if &FF = &F
  88. }
  89. #ifdef IGL_STATIC_LIBRARY
  90. // Explicit template instantiation
  91. template void igl::bfs_orient<Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  92. #endif