doublearea.cpp 13 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248
  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 "parallel_for.h"
  11. #include "sort.h"
  12. #include <cassert>
  13. #include <iostream>
  14. template <typename DerivedV, typename DerivedF, typename DeriveddblA>
  15. IGL_INLINE void igl::doublearea(
  16. const Eigen::MatrixBase<DerivedV> & V,
  17. const Eigen::MatrixBase<DerivedF> & F,
  18. Eigen::PlainObjectBase<DeriveddblA> & dblA)
  19. {
  20. // quads are handled by a specialized function
  21. if (F.cols() == 4) return doublearea_quad(V,F,dblA);
  22. const int dim = V.cols();
  23. // Only support triangles
  24. assert(F.cols() == 3);
  25. const size_t m = F.rows();
  26. // Compute edge lengths
  27. Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 3> l;
  28. // Projected area helper
  29. const auto & proj_doublearea =
  30. [&V,&F](const int x, const int y, const int f)->double
  31. {
  32. auto rx = V(F(f,0),x)-V(F(f,2),x);
  33. auto sx = V(F(f,1),x)-V(F(f,2),x);
  34. auto ry = V(F(f,0),y)-V(F(f,2),y);
  35. auto sy = V(F(f,1),y)-V(F(f,2),y);
  36. return rx*sy - ry*sx;
  37. };
  38. switch(dim)
  39. {
  40. case 3:
  41. {
  42. dblA = DeriveddblA::Zero(m,1);
  43. for(size_t f = 0;f<m;f++)
  44. {
  45. for(int d = 0;d<3;d++)
  46. {
  47. double dblAd = proj_doublearea(d,(d+1)%3,f);
  48. dblA(f) += dblAd*dblAd;
  49. }
  50. }
  51. dblA = dblA.array().sqrt().eval();
  52. break;
  53. }
  54. case 2:
  55. {
  56. dblA.resize(m,1);
  57. for(size_t f = 0;f<m;f++)
  58. {
  59. dblA(f) = proj_doublearea(0,1,f);
  60. }
  61. break;
  62. }
  63. default:
  64. {
  65. edge_lengths(V,F,l);
  66. return doublearea(l,0.,dblA);
  67. }
  68. }
  69. }
  70. template <
  71. typename DerivedA,
  72. typename DerivedB,
  73. typename DerivedC,
  74. typename DerivedD>
  75. IGL_INLINE void igl::doublearea(
  76. const Eigen::MatrixBase<DerivedA> & A,
  77. const Eigen::MatrixBase<DerivedB> & B,
  78. const Eigen::MatrixBase<DerivedC> & C,
  79. Eigen::PlainObjectBase<DerivedD> & D)
  80. {
  81. assert((B.cols() == A.cols()) && "dimensions of A and B should match");
  82. assert((C.cols() == A.cols()) && "dimensions of A and C should match");
  83. assert(A.rows() == B.rows() && "corners should have same length");
  84. assert(A.rows() == C.rows() && "corners should have same length");
  85. switch(A.cols())
  86. {
  87. case 2:
  88. {
  89. // For 2d compute signed area
  90. const auto & R = A-C;
  91. const auto & S = B-C;
  92. D = R.col(0).array()*S.col(1).array() - R.col(1).array()*S.col(0).array();
  93. break;
  94. }
  95. default:
  96. {
  97. Eigen::Matrix<typename DerivedD::Scalar,DerivedD::RowsAtCompileTime,3>
  98. uL(A.rows(),3);
  99. uL.col(0) = (B-C).rowwise().norm();
  100. uL.col(1) = (C-A).rowwise().norm();
  101. uL.col(2) = (A-B).rowwise().norm();
  102. doublearea(uL,D);
  103. }
  104. }
  105. }
  106. template <
  107. typename DerivedA,
  108. typename DerivedB,
  109. typename DerivedC>
  110. IGL_INLINE typename DerivedA::Scalar igl::doublearea_single(
  111. const Eigen::MatrixBase<DerivedA> & A,
  112. const Eigen::MatrixBase<DerivedB> & B,
  113. const Eigen::MatrixBase<DerivedC> & C)
  114. {
  115. assert(A.size() == 2 && "Vertices should be 2D");
  116. assert(B.size() == 2 && "Vertices should be 2D");
  117. assert(C.size() == 2 && "Vertices should be 2D");
  118. auto r = A-C;
  119. auto s = B-C;
  120. return r(0)*s(1) - r(1)*s(0);
  121. }
  122. template <typename Derivedl, typename DeriveddblA>
  123. IGL_INLINE void igl::doublearea(
  124. const Eigen::MatrixBase<Derivedl> & ul,
  125. Eigen::PlainObjectBase<DeriveddblA> & dblA)
  126. {
  127. // Default is to leave NaNs and fire asserts in debug mode
  128. return doublearea(ul,0.0/0.0,dblA);
  129. }
  130. template <typename Derivedl, typename DeriveddblA>
  131. IGL_INLINE void igl::doublearea(
  132. const Eigen::MatrixBase<Derivedl> & ul,
  133. const typename Derivedl::Scalar nan_replacement,
  134. Eigen::PlainObjectBase<DeriveddblA> & dblA)
  135. {
  136. using namespace Eigen;
  137. using namespace std;
  138. typedef typename Derivedl::Index Index;
  139. // Only support triangles
  140. assert(ul.cols() == 3);
  141. // Number of triangles
  142. const Index m = ul.rows();
  143. Eigen::Matrix<typename Derivedl::Scalar, Eigen::Dynamic, 3> l;
  144. MatrixXi _;
  145. //
  146. // "Miscalculating Area and Angles of a Needle-like Triangle"
  147. // https://people.eecs.berkeley.edu/~wkahan/Triangle.pdf
  148. igl::sort(ul,2,false,l,_);
  149. // semiperimeters
  150. //Matrix<typename Derivedl::Scalar,Dynamic,1> s = l.rowwise().sum()*0.5;
  151. //assert((Index)s.rows() == m);
  152. // resize output
  153. dblA.resize(l.rows(),1);
  154. parallel_for(
  155. m,
  156. [&l,&dblA,&nan_replacement](const int i)
  157. {
  158. // Kahan's Heron's formula
  159. typedef typename Derivedl::Scalar Scalar;
  160. const Scalar arg =
  161. (l(i,0)+(l(i,1)+l(i,2)))*
  162. (l(i,2)-(l(i,0)-l(i,1)))*
  163. (l(i,2)+(l(i,0)-l(i,1)))*
  164. (l(i,0)+(l(i,1)-l(i,2)));
  165. dblA(i) = 2.0*0.25*sqrt(arg);
  166. // Alec: If the input edge lengths were computed from floating point
  167. // vertex positions then there's no guarantee that they fulfill the
  168. // triangle inequality (in their floating point approximations). For
  169. // nearly degenerate triangles the round-off error during side-length
  170. // computation may be larger than (or rather smaller) than the height of
  171. // the triangle. In "Lecture Notes on Geometric Robustness" Shewchuck 09,
  172. // Section 3.1 http://www.cs.berkeley.edu/~jrs/meshpapers/robnotes.pdf,
  173. // he recommends computing the triangle areas for 2D and 3D using 2D
  174. // signed areas computed with determinants.
  175. assert(
  176. (nan_replacement == nan_replacement ||
  177. (l(i,2) - (l(i,0)-l(i,1)))>=0)
  178. && "Side lengths do not obey the triangle inequality.");
  179. if(dblA(i) != dblA(i))
  180. {
  181. dblA(i) = nan_replacement;
  182. }
  183. assert(dblA(i) == dblA(i) && "DOUBLEAREA() PRODUCED NaN");
  184. },
  185. 1000l);
  186. }
  187. template <typename DerivedV, typename DerivedF, typename DeriveddblA>
  188. IGL_INLINE void igl::doublearea_quad(
  189. const Eigen::MatrixBase<DerivedV> & V,
  190. const Eigen::MatrixBase<DerivedF> & F,
  191. Eigen::PlainObjectBase<DeriveddblA> & dblA)
  192. {
  193. assert(V.cols() == 3); // Only supports points in 3D
  194. assert(F.cols() == 4); // Only support quads
  195. const size_t m = F.rows();
  196. // Split the quads into triangles
  197. Eigen::MatrixXi Ft(F.rows()*2,3);
  198. for(size_t i=0; i<m;++i)
  199. {
  200. Ft.row(i*2 ) << F(i,0), F(i,1), F(i,2);
  201. Ft.row(i*2 + 1) << F(i,2), F(i,3), F(i,0);
  202. }
  203. // Compute areas
  204. Eigen::VectorXd doublearea_tri;
  205. igl::doublearea(V,Ft,doublearea_tri);
  206. dblA.resize(F.rows(),1);
  207. for(unsigned i=0; i<F.rows();++i)
  208. {
  209. dblA(i) = doublearea_tri(i*2) + doublearea_tri(i*2 + 1);
  210. }
  211. }
  212. #ifdef IGL_STATIC_LIBRARY
  213. // Explicit template instantiation
  214. // generated by autoexplicit.sh
  215. template void igl::doublearea<Eigen::Matrix<double, -1, -1, 1, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 1, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  216. template void igl::doublearea<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  217. template void igl::doublearea<Eigen::Matrix<float, -1, 3, 1, -1, 3>, Eigen::Matrix<unsigned int, -1, 3, 1, -1, 3>, Eigen::Matrix<float, -1, 1, 0, -1, 1> >(Eigen::MatrixBase<Eigen::Matrix<float, -1, 3, 1, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<unsigned int, -1, 3, 1, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<float, -1, 1, 0, -1, 1> >&);
  218. template void igl::doublearea<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  219. 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::MatrixBase<Eigen::Matrix<float, -1, 3, 1, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<unsigned int, -1, -1, 1, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  220. 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::MatrixBase<Eigen::Matrix<double, -1, 3, 1, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<unsigned int, -1, -1, 1, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  221. 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::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  222. template void igl::doublearea<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  223. 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::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  224. 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::MatrixBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  225. 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::MatrixBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  226. 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::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  227. template void igl::doublearea<Eigen::Matrix<double, 1, 3, 1, 1, 3>, Eigen::Matrix<double, 1, 3, 1, 1, 3>, Eigen::Matrix<double, 1, 3, 1, 1, 3>, Eigen::Matrix<double, 1, 1, 0, 1, 1> >(Eigen::MatrixBase<Eigen::Matrix<double, 1, 3, 1, 1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<double, 1, 3, 1, 1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<double, 1, 3, 1, 1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, 1, 1, 0, 1, 1> >&);
  228. template void igl::doublearea<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  229. template void igl::doublearea<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  230. template void igl::doublearea<Eigen::Matrix<double, -1, 3, 1, -1, 3>, Eigen::Matrix<int, -1, 3, 1, -1, 3>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, 3, 1, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 3, 1, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  231. template Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar igl::doublearea_single<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&);
  232. 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::MatrixBase<Eigen::Matrix<double, 2, 1, 0, 2, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, 2, 1, 0, 2, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, 2, 1, 0, 2, 1> > const&);
  233. #endif