doublearea.cpp 8.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2013 Alec Jacobson <alecjacobson@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 "doublearea.h"
  9. #include "edge_lengths.h"
  10. #include "sort.h"
  11. #include <cassert>
  12. #include <iostream>
  13. template <typename DerivedV, typename DerivedF, typename DeriveddblA>
  14. IGL_INLINE void igl::doublearea(
  15. const Eigen::PlainObjectBase<DerivedV> & V,
  16. const Eigen::PlainObjectBase<DerivedF> & F,
  17. Eigen::PlainObjectBase<DeriveddblA> & dblA)
  18. {
  19. if (F.cols() == 4) // quads are handled by a specialized function
  20. return doublearea_quad(V,F,dblA);
  21. const int dim = V.cols();
  22. // Only support triangles
  23. assert(F.cols() == 3);
  24. const size_t m = F.rows();
  25. // Compute edge lengths
  26. Eigen::PlainObjectBase<DerivedV> l;
  27. // "Lecture Notes on Geometric Robustness" Shewchuck 09, Section 3.1
  28. // http://www.cs.berkeley.edu/~jrs/meshpapers/robnotes.pdf
  29. // Projected area helper
  30. const auto & proj_doublearea =
  31. [&V,&F](const int x, const int y, const int f)->double
  32. {
  33. auto rx = V(F(f,0),x)-V(F(f,2),x);
  34. auto sx = V(F(f,1),x)-V(F(f,2),x);
  35. auto ry = V(F(f,0),y)-V(F(f,2),y);
  36. auto sy = V(F(f,1),y)-V(F(f,2),y);
  37. return rx*sy - ry*sx;
  38. };
  39. switch(dim)
  40. {
  41. case 3:
  42. {
  43. dblA = Eigen::PlainObjectBase<DeriveddblA>::Zero(m,1);
  44. for(size_t f = 0;f<m;f++)
  45. {
  46. for(int d = 0;d<3;d++)
  47. {
  48. double dblAd = proj_doublearea(d,(d+1)%3,f);
  49. dblA(f) += dblAd*dblAd;
  50. }
  51. }
  52. dblA = dblA.array().sqrt().eval();
  53. break;
  54. }
  55. case 2:
  56. {
  57. dblA.resize(m,1);
  58. for(size_t f = 0;f<m;f++)
  59. {
  60. dblA(f) = proj_doublearea(0,1,f);
  61. }
  62. break;
  63. }
  64. default:
  65. {
  66. edge_lengths(V,F,l);
  67. return doublearea(l,dblA);
  68. }
  69. }
  70. }
  71. template <
  72. typename DerivedA,
  73. typename DerivedB,
  74. typename DerivedC,
  75. typename DerivedD>
  76. IGL_INLINE void igl::doublearea(
  77. const Eigen::PlainObjectBase<DerivedA> & A,
  78. const Eigen::PlainObjectBase<DerivedB> & B,
  79. const Eigen::PlainObjectBase<DerivedC> & C,
  80. Eigen::PlainObjectBase<DerivedD> & D)
  81. {
  82. assert(A.cols() == 2 && "corners should be 2d");
  83. assert(B.cols() == 2 && "corners should be 2d");
  84. assert(C.cols() == 2 && "corners should be 2d");
  85. assert(A.rows() == B.rows() && "corners should have same length");
  86. assert(A.rows() == C.rows() && "corners should have same length");
  87. const auto & R = A-C;
  88. const auto & S = B-C;
  89. D = R.col(0).array()*S.col(1).array() - R.col(1).array()*S.col(0).array();
  90. }
  91. template <
  92. typename DerivedA,
  93. typename DerivedB,
  94. typename DerivedC>
  95. IGL_INLINE typename DerivedA::Scalar igl::doublearea_single(
  96. const Eigen::PlainObjectBase<DerivedA> & A,
  97. const Eigen::PlainObjectBase<DerivedB> & B,
  98. const Eigen::PlainObjectBase<DerivedC> & C)
  99. {
  100. auto r = A-C;
  101. auto s = B-C;
  102. return r(0)*s(1) - r(1)*s(0);
  103. }
  104. template <typename Derivedl, typename DeriveddblA>
  105. IGL_INLINE void igl::doublearea(
  106. const Eigen::PlainObjectBase<Derivedl> & ul,
  107. Eigen::PlainObjectBase<DeriveddblA> & dblA)
  108. {
  109. using namespace Eigen;
  110. using namespace std;
  111. typedef typename Derivedl::Index Index;
  112. // Only support triangles
  113. assert(ul.cols() == 3);
  114. // Number of triangles
  115. const Index m = ul.rows();
  116. Eigen::PlainObjectBase<Derivedl> l;
  117. MatrixXi _;
  118. sort(ul,2,false,l,_);
  119. // semiperimeters
  120. Matrix<typename Derivedl::Scalar,Dynamic,1> s = l.rowwise().sum()*0.5;
  121. assert((size_t)s.rows() == m);
  122. // resize output
  123. dblA.resize(l.rows(),1);
  124. // Minimum number of iterms per openmp thread
  125. #ifndef IGL_OMP_MIN_VALUE
  126. # define IGL_OMP_MIN_VALUE 1000
  127. #endif
  128. #pragma omp parallel for if (m>IGL_OMP_MIN_VALUE)
  129. for(Index i = 0;i<m;i++)
  130. {
  131. //// Heron's formula for area
  132. //const typename Derivedl::Scalar arg =
  133. // s(i)*(s(i)-l(i,0))*(s(i)-l(i,1))*(s(i)-l(i,2));
  134. //assert(arg>=0);
  135. //dblA(i) = 2.0*sqrt(arg);
  136. // Kahan's Heron's formula
  137. const typename Derivedl::Scalar arg =
  138. (l(i,0)+(l(i,1)+l(i,2)))*
  139. (l(i,2)-(l(i,0)-l(i,1)))*
  140. (l(i,2)+(l(i,0)-l(i,1)))*
  141. (l(i,0)+(l(i,1)-l(i,2)));
  142. dblA(i) = 2.0*0.25*sqrt(arg);
  143. assert( l(i,2) - (l(i,0)-l(i,1)) && "FAILED KAHAN'S ASSERTION");
  144. assert(dblA(i) == dblA(i) && "DOUBLEAREA() PRODUCED NaN");
  145. }
  146. }
  147. template <typename DerivedV, typename DerivedF, typename DeriveddblA>
  148. IGL_INLINE void igl::doublearea_quad(
  149. const Eigen::PlainObjectBase<DerivedV> & V,
  150. const Eigen::PlainObjectBase<DerivedF> & F,
  151. Eigen::PlainObjectBase<DeriveddblA> & dblA)
  152. {
  153. assert(V.cols() == 3); // Only supports points in 3D
  154. assert(F.cols() == 4); // Only support quads
  155. const size_t m = F.rows();
  156. // Split the quads into triangles
  157. Eigen::MatrixXi Ft(F.rows()*2,3);
  158. for(size_t i=0; i<m;++i)
  159. {
  160. Ft.row(i*2 ) << F(i,0), F(i,1), F(i,2);
  161. Ft.row(i*2 + 1) << F(i,2), F(i,3), F(i,0);
  162. }
  163. // Compute areas
  164. Eigen::VectorXd doublearea_tri;
  165. igl::doublearea(V,Ft,doublearea_tri);
  166. dblA.resize(F.rows(),1);
  167. for(unsigned i=0; i<F.rows();++i)
  168. dblA(i) = doublearea_tri(i*2) + doublearea_tri(i*2 + 1);
  169. }
  170. #ifdef IGL_STATIC_LIBRARY
  171. // Explicit template specialization
  172. // generated by autoexplicit.sh
  173. template void igl::doublearea<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  174. // generated by autoexplicit.sh
  175. template void igl::doublearea<Eigen::Matrix<float, -1, 3, 1, -1, 3>, Eigen::Matrix<unsigned int, -1, -1, 1, -1, -1>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<float, -1, 3, 1, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<unsigned int, -1, -1, 1, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  176. template void igl::doublearea<Eigen::Matrix<double, -1, 3, 1, -1, 3>, Eigen::Matrix<unsigned int, -1, -1, 1, -1, -1>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 1, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<unsigned int, -1, -1, 1, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  177. template void igl::doublearea<Eigen::Matrix<double, -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<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  178. template void igl::doublearea<Eigen::Matrix<double, -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> >&);
  179. template void igl::doublearea<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::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<int, -1, 1, 0, -1, 1> >&);
  180. template void igl::doublearea<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >(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, 1, 0, -1, 1> >&);
  181. template Eigen::Matrix<double, 2, 1, 0, 2, 1>::Scalar igl::doublearea_single<Eigen::Matrix<double, 2, 1, 0, 2, 1>, Eigen::Matrix<double, 2, 1, 0, 2, 1>, Eigen::Matrix<double, 2, 1, 0, 2, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, 2, 1, 0, 2, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, 2, 1, 0, 2, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, 2, 1, 0, 2, 1> > const&);
  182. template void igl::doublearea<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(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, -1, 0, -1, -1> >&);
  183. template void igl::doublearea<Eigen::Matrix<double, -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<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  184. #endif