intersect_other.cpp 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2014 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 "intersect_other.h"
  9. #include "CGAL_includes.hpp"
  10. #include "mesh_to_cgal_triangle_list.h"
  11. #include "remesh_intersections.h"
  12. #include "../../remove_unreferenced.h"
  13. #include "../../PlainMatrix.h"
  14. #ifndef IGL_FIRST_HIT_EXCEPTION
  15. #define IGL_FIRST_HIT_EXCEPTION 10
  16. #endif
  17. // Un-exposed helper functions
  18. namespace igl
  19. {
  20. namespace copyleft
  21. {
  22. namespace cgal
  23. {
  24. template <typename Index>
  25. static IGL_INLINE void push_result(
  26. const int f,
  27. const int f_other,
  28. const CGAL::Object & result,
  29. std::map<
  30. Index,
  31. std::vector<std::pair<Index, CGAL::Object> > > &
  32. offending)
  33. //std::map<
  34. // std::pair<typename DerivedF::Index,typename DerivedF::Index>,
  35. // std::vector<typename DerivedF::Index> > & edge2faces)
  36. {
  37. if(offending.count(f) == 0)
  38. {
  39. // first time marking, initialize with new id and empty list
  40. offending[f] = {};
  41. }
  42. offending[f].push_back({f_other,result});
  43. }
  44. template <
  45. typename Kernel,
  46. typename DerivedVA,
  47. typename DerivedFA,
  48. typename DerivedVB,
  49. typename DerivedFB,
  50. typename DerivedIF,
  51. typename DerivedVVAB,
  52. typename DerivedFFAB,
  53. typename DerivedJAB,
  54. typename DerivedIMAB>
  55. static IGL_INLINE bool intersect_other_helper(
  56. const Eigen::MatrixBase<DerivedVA> & VA,
  57. const Eigen::MatrixBase<DerivedFA> & FA,
  58. const Eigen::MatrixBase<DerivedVB> & VB,
  59. const Eigen::MatrixBase<DerivedFB> & FB,
  60. const RemeshSelfIntersectionsParam & params,
  61. Eigen::PlainObjectBase<DerivedIF> & IF,
  62. Eigen::PlainObjectBase<DerivedVVAB> & VVAB,
  63. Eigen::PlainObjectBase<DerivedFFAB> & FFAB,
  64. Eigen::PlainObjectBase<DerivedJAB> & JAB,
  65. Eigen::PlainObjectBase<DerivedIMAB> & IMAB)
  66. {
  67. typedef typename DerivedFA::Index Index;
  68. typedef CGAL::Triangle_3<Kernel> Triangle_3;
  69. //// Axis-align boxes for all-pairs self-intersection detection
  70. typedef std::vector<Triangle_3> Triangles;
  71. typedef typename Triangles::iterator TrianglesIterator;
  72. typedef
  73. CGAL::Box_intersection_d::Box_with_handle_d<double,3,TrianglesIterator>
  74. Box;
  75. typedef
  76. std::map<Index,std::vector<std::pair<Index,CGAL::Object> > >
  77. OffendingMap;
  78. Triangles TA,TB;
  79. // Compute and process self intersections
  80. mesh_to_cgal_triangle_list(VA,FA,TA);
  81. mesh_to_cgal_triangle_list(VB,FB,TB);
  82. // http://www.cgal.org/Manual/latest/doc_html/cgal_manual/Box_intersection_d/Chapter_main.html#Section_63.5
  83. // Create the corresponding vector of bounding boxes
  84. std::vector<Box> A_boxes,B_boxes;
  85. const auto box_up = [](Triangles & T, std::vector<Box> & boxes) -> void
  86. {
  87. boxes.reserve(T.size());
  88. for (
  89. TrianglesIterator tit = T.begin();
  90. tit != T.end();
  91. ++tit)
  92. {
  93. boxes.push_back(Box(tit->bbox(), tit));
  94. }
  95. };
  96. box_up(TA,A_boxes);
  97. box_up(TB,B_boxes);
  98. OffendingMap offendingA,offendingB;
  99. //EdgeMap edge2facesA,edge2facesB;
  100. std::list<int> lIF;
  101. const auto cb = [&](const Box &a, const Box &b) -> void
  102. {
  103. // index in F and T
  104. int fa = a.handle()-TA.begin();
  105. int fb = b.handle()-TB.begin();
  106. const Triangle_3 & A = *a.handle();
  107. const Triangle_3 & B = *b.handle();
  108. if(CGAL::do_intersect(A,B))
  109. {
  110. // There was an intersection
  111. lIF.push_back(fa);
  112. lIF.push_back(fb);
  113. if(params.first_only)
  114. {
  115. throw IGL_FIRST_HIT_EXCEPTION;
  116. }
  117. if(!params.detect_only)
  118. {
  119. CGAL::Object result = CGAL::intersection(A,B);
  120. push_result(fa,fb,result,offendingA);
  121. push_result(fb,fa,result,offendingB);
  122. }
  123. }
  124. };
  125. try{
  126. CGAL::box_intersection_d(
  127. A_boxes.begin(), A_boxes.end(),
  128. B_boxes.begin(), B_boxes.end(),
  129. cb,
  130. std::ptrdiff_t(params.cutoff));
  131. }catch(int e)
  132. {
  133. // Rethrow if not FIRST_HIT_EXCEPTION
  134. if(e != IGL_FIRST_HIT_EXCEPTION)
  135. {
  136. throw e;
  137. }
  138. // Otherwise just fall through
  139. }
  140. // Convert lIF to Eigen matrix
  141. assert(lIF.size()%2 == 0);
  142. IF.resize(lIF.size()/2,2);
  143. {
  144. int i=0;
  145. for(
  146. std::list<int>::const_iterator ifit = lIF.begin();
  147. ifit!=lIF.end();
  148. )
  149. {
  150. IF(i,0) = (*ifit);
  151. ifit++;
  152. IF(i,1) = (*ifit);
  153. ifit++;
  154. i++;
  155. }
  156. }
  157. if(!params.detect_only)
  158. {
  159. // Obsolete, now remesh_intersections expects a single mesh
  160. // remesh_intersections(VA,FA,TA,offendingA,VVA,FFA,JA,IMA);
  161. // remesh_intersections(VB,FB,TB,offendingB,VVB,FFB,JB,IMB);
  162. // Combine mesh and offending maps
  163. PlainMatrix<DerivedVA,Eigen::Dynamic> VAB(VA.rows()+VB.rows(),3);
  164. VAB<<VA,VB;
  165. PlainMatrix<DerivedFA,Eigen::Dynamic> FAB(FA.rows()+FB.rows(),3);
  166. FAB<<FA,(FB.array()+VA.rows());
  167. Triangles TAB;
  168. TAB.reserve(TA.size()+TB.size());
  169. TAB.insert(TAB.end(),TA.begin(),TA.end());
  170. TAB.insert(TAB.end(),TB.begin(),TB.end());
  171. OffendingMap offending;
  172. //offending.reserve(offendingA.size() + offendingB.size());
  173. for (const auto itr : offendingA)
  174. {
  175. // Remap offenders in FB to FAB
  176. auto offenders = itr.second;
  177. for(auto & offender : offenders)
  178. {
  179. offender.first += FA.rows();
  180. }
  181. offending[itr.first] = offenders;
  182. }
  183. for (const auto itr : offendingB)
  184. {
  185. // Store offenders for FB according to place in FAB
  186. offending[FA.rows() + itr.first] = itr.second;
  187. }
  188. remesh_intersections(
  189. VAB,FAB,TAB,offending,
  190. params.stitch_all,params.slow_and_more_precise_rounding,
  191. VVAB,FFAB,JAB,IMAB);
  192. }
  193. return IF.rows() > 0;
  194. }
  195. }
  196. }
  197. }
  198. template <
  199. typename DerivedVA,
  200. typename DerivedFA,
  201. typename DerivedVB,
  202. typename DerivedFB,
  203. typename DerivedIF,
  204. typename DerivedVVAB,
  205. typename DerivedFFAB,
  206. typename DerivedJAB,
  207. typename DerivedIMAB>
  208. IGL_INLINE bool igl::copyleft::cgal::intersect_other(
  209. const Eigen::MatrixBase<DerivedVA> & VA,
  210. const Eigen::MatrixBase<DerivedFA> & FA,
  211. const Eigen::MatrixBase<DerivedVB> & VB,
  212. const Eigen::MatrixBase<DerivedFB> & FB,
  213. const RemeshSelfIntersectionsParam & params,
  214. Eigen::PlainObjectBase<DerivedIF> & IF,
  215. Eigen::PlainObjectBase<DerivedVVAB> & VVAB,
  216. Eigen::PlainObjectBase<DerivedFFAB> & FFAB,
  217. Eigen::PlainObjectBase<DerivedJAB> & JAB,
  218. Eigen::PlainObjectBase<DerivedIMAB> & IMAB)
  219. {
  220. typedef typename DerivedVA::Scalar VAScalar;
  221. typedef typename DerivedVB::Scalar VBScalar;
  222. if(params.detect_only && ! (
  223. std::is_same<VAScalar,CGAL::Epeck::FT>::value ||
  224. std::is_same<VBScalar,CGAL::Epeck::FT>::value
  225. ))
  226. {
  227. return intersect_other_helper<CGAL::Epick>
  228. (VA,FA,VB,FB,params,IF,VVAB,FFAB,JAB,IMAB);
  229. }else
  230. {
  231. return intersect_other_helper<CGAL::Epeck>
  232. (VA,FA,VB,FB,params,IF,VVAB,FFAB,JAB,IMAB);
  233. }
  234. }
  235. template <
  236. typename DerivedVA,
  237. typename DerivedFA,
  238. typename DerivedVB,
  239. typename DerivedFB,
  240. typename DerivedIF>
  241. IGL_INLINE bool igl::copyleft::cgal::intersect_other(
  242. const Eigen::MatrixBase<DerivedVA> & VA,
  243. const Eigen::MatrixBase<DerivedFA> & FA,
  244. const Eigen::MatrixBase<DerivedVB> & VB,
  245. const Eigen::MatrixBase<DerivedFB> & FB,
  246. const bool first_only,
  247. Eigen::PlainObjectBase<DerivedIF> & IF)
  248. {
  249. PlainMatrix<DerivedVA,Eigen::Dynamic> VVAB;
  250. PlainMatrix<DerivedFA,Eigen::Dynamic> FFAB;
  251. Eigen::VectorXi JAB, IMAB;
  252. return intersect_other(
  253. VA,FA,VB,FB,{true,first_only},IF,VVAB,FFAB,JAB,IMAB);
  254. }
  255. #ifdef IGL_STATIC_LIBRARY
  256. // Explicit template instantiation
  257. template bool igl::copyleft::cgal::intersect_other<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::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<CGAL::Epeck::FT, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(
  258. Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&,
  259. Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&,
  260. Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&,
  261. Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, igl::copyleft::cgal::RemeshSelfIntersectionsParam const&,
  262. Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&,
  263. Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Epeck::FT, -1, 3, 0, -1, 3> >&,
  264. Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&,
  265. Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&,
  266. Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  267. template bool igl::copyleft::cgal::intersect_other<Eigen::Matrix<CGAL::Epeck::FT, -1, -1, 1, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<CGAL::Epeck::FT, -1, -1, 1, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<CGAL::Epeck::FT, -1, -1, 1, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(
  268. Eigen::MatrixBase<Eigen::Matrix<CGAL::Epeck::FT, -1, -1, 1, -1, -1> > const&,
  269. Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&,
  270. Eigen::MatrixBase<Eigen::Matrix<CGAL::Epeck::FT, -1, -1, 1, -1, -1> > const&,
  271. Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, igl::copyleft::cgal::RemeshSelfIntersectionsParam const&,
  272. Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&,
  273. Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Epeck::FT, -1, -1, 1, -1, -1> >&,
  274. Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&,
  275. Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&,
  276. Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  277. #endif