mesh_boolean.cpp 17 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2015 Alec Jacobson <[email protected]>
  4. // Qingnan Zhou <[email protected]>
  5. //
  6. // This Source Code Form is subject to the terms of the Mozilla Public License
  7. // v. 2.0. If a copy of the MPL was not distributed with this file, You can
  8. // obtain one at http://mozilla.org/MPL/2.0/.
  9. //
  10. #include "mesh_boolean.h"
  11. #include "assign_scalar.h"
  12. #include "extract_cells.h"
  13. #include "mesh_boolean_type_to_funcs.h"
  14. #include "propagate_winding_numbers.h"
  15. #include "relabel_small_immersed_cells.h"
  16. #include "remesh_self_intersections.h"
  17. #include "string_to_mesh_boolean_type.h"
  18. #include "../../combine.h"
  19. #include "../../cumsum.h"
  20. #include "../../extract_manifold_patches.h"
  21. #include "../../get_seconds.h"
  22. #include "../../parallel_for.h"
  23. #include "../../remove_unreferenced.h"
  24. #include "../../resolve_duplicated_faces.h"
  25. #include "../../unique_edge_map.h"
  26. #include "../../unique_simplices.h"
  27. #include "../../C_STR.h"
  28. #include <CGAL/Exact_predicates_exact_constructions_kernel.h>
  29. #include <algorithm>
  30. //#define MESH_BOOLEAN_TIMING
  31. //#define DOUBLE_CHECK_EXACT_OUTPUT
  32. //#define SMALL_CELL_REMOVAL
  33. template <
  34. typename DerivedVA,
  35. typename DerivedFA,
  36. typename DerivedVB,
  37. typename DerivedFB,
  38. typename DerivedVC,
  39. typename DerivedFC,
  40. typename DerivedJ>
  41. IGL_INLINE bool igl::copyleft::cgal::mesh_boolean(
  42. const Eigen::MatrixBase<DerivedVA > & VA,
  43. const Eigen::MatrixBase<DerivedFA > & FA,
  44. const Eigen::MatrixBase<DerivedVB > & VB,
  45. const Eigen::MatrixBase<DerivedFB > & FB,
  46. const MeshBooleanType & type,
  47. Eigen::PlainObjectBase<DerivedVC > & VC,
  48. Eigen::PlainObjectBase<DerivedFC > & FC,
  49. Eigen::PlainObjectBase<DerivedJ > & J)
  50. {
  51. std::function<int(const int, const int)> keep;
  52. std::function<int(const Eigen::Matrix<int,1,Eigen::Dynamic>) > wind_num_op;
  53. mesh_boolean_type_to_funcs(type,wind_num_op,keep);
  54. return mesh_boolean(VA,FA,VB,FB,wind_num_op,keep,VC,FC,J);
  55. }
  56. template <
  57. typename DerivedVA,
  58. typename DerivedFA,
  59. typename DerivedVB,
  60. typename DerivedFB,
  61. typename DerivedVC,
  62. typename DerivedFC,
  63. typename DerivedJ>
  64. IGL_INLINE bool igl::copyleft::cgal::mesh_boolean(
  65. const Eigen::MatrixBase<DerivedVA > & VA,
  66. const Eigen::MatrixBase<DerivedFA > & FA,
  67. const Eigen::MatrixBase<DerivedVB > & VB,
  68. const Eigen::MatrixBase<DerivedFB > & FB,
  69. const std::string & type_str,
  70. Eigen::PlainObjectBase<DerivedVC > & VC,
  71. Eigen::PlainObjectBase<DerivedFC > & FC,
  72. Eigen::PlainObjectBase<DerivedJ > & J)
  73. {
  74. return mesh_boolean(
  75. VA,FA,VB,FB,string_to_mesh_boolean_type(type_str),VC,FC,J);
  76. }
  77. template <
  78. typename DerivedVA,
  79. typename DerivedFA,
  80. typename DerivedVB,
  81. typename DerivedFB,
  82. typename DerivedVC,
  83. typename DerivedFC,
  84. typename DerivedJ>
  85. IGL_INLINE bool igl::copyleft::cgal::mesh_boolean(
  86. const Eigen::MatrixBase<DerivedVA> & VA,
  87. const Eigen::MatrixBase<DerivedFA> & FA,
  88. const Eigen::MatrixBase<DerivedVB> & VB,
  89. const Eigen::MatrixBase<DerivedFB> & FB,
  90. const std::function<int(const Eigen::Matrix<int,1,Eigen::Dynamic>) >& wind_num_op,
  91. const std::function<int(const int, const int)> & keep,
  92. Eigen::PlainObjectBase<DerivedVC > & VC,
  93. Eigen::PlainObjectBase<DerivedFC > & FC,
  94. Eigen::PlainObjectBase<DerivedJ > & J)
  95. {
  96. // Generate combined mesh (VA,FA,VB,FB) -> (V,F)
  97. Eigen::Matrix<size_t,2,1> sizes(FA.rows(),FB.rows());
  98. // TODO: This is a precision template **bug** that results in failure to
  99. // compile. If DerivedVA::Scalar is double and DerivedVB::Scalar is
  100. // CGAL::Epeck::FT then the following assignment will not compile. This
  101. // implies that VA must have the trumping precision (and a valid assignment
  102. // operator from VB's type).
  103. Eigen::Matrix<typename DerivedVA::Scalar,Eigen::Dynamic,3> VV(VA.rows() + VB.rows(), 3);
  104. Eigen::Matrix<typename DerivedFA::Scalar,Eigen::Dynamic,3> FF(FA.rows() + FB.rows(), 3);
  105. // Can't use comma initializer
  106. for(int a = 0;a<VA.rows();a++)
  107. {
  108. for(int d = 0;d<3;d++) VV(a,d) = VA(a,d);
  109. }
  110. for(int b = 0;b<VB.rows();b++)
  111. {
  112. for(int d = 0;d<3;d++) VV(VA.rows()+b,d) = VB(b,d);
  113. }
  114. FF.block(0, 0, FA.rows(), 3) = FA;
  115. // Eigen struggles to assign nothing to nothing and will assert if FB is empty
  116. if(FB.rows() > 0)
  117. {
  118. FF.block(FA.rows(), 0, FB.rows(), 3) = FB.array() + VA.rows();
  119. }
  120. return mesh_boolean(VV,FF,sizes,wind_num_op,keep,VC,FC,J);
  121. }
  122. template <
  123. typename DerivedV,
  124. typename DerivedF,
  125. typename DerivedVC,
  126. typename DerivedFC,
  127. typename DerivedJ>
  128. IGL_INLINE bool igl::copyleft::cgal::mesh_boolean(
  129. const std::vector<DerivedV > & Vlist,
  130. const std::vector<DerivedF > & Flist,
  131. const std::function<int(const Eigen::Matrix<int,1,Eigen::Dynamic>) >& wind_num_op,
  132. const std::function<int(const int, const int)> & keep,
  133. Eigen::PlainObjectBase<DerivedVC > & VC,
  134. Eigen::PlainObjectBase<DerivedFC > & FC,
  135. Eigen::PlainObjectBase<DerivedJ > & J)
  136. {
  137. DerivedV VV;
  138. DerivedF FF;
  139. Eigen::Matrix<size_t,Eigen::Dynamic,1> Vsizes,Fsizes;
  140. igl::combine(Vlist,Flist,VV,FF,Vsizes,Fsizes);
  141. return mesh_boolean(VV,FF,Fsizes,wind_num_op,keep,VC,FC,J);
  142. }
  143. template <
  144. typename DerivedV,
  145. typename DerivedF,
  146. typename DerivedVC,
  147. typename DerivedFC,
  148. typename DerivedJ>
  149. IGL_INLINE bool igl::copyleft::cgal::mesh_boolean(
  150. const std::vector<DerivedV > & Vlist,
  151. const std::vector<DerivedF > & Flist,
  152. const MeshBooleanType & type,
  153. Eigen::PlainObjectBase<DerivedVC > & VC,
  154. Eigen::PlainObjectBase<DerivedFC > & FC,
  155. Eigen::PlainObjectBase<DerivedJ > & J)
  156. {
  157. DerivedV VV;
  158. DerivedF FF;
  159. Eigen::Matrix<size_t,Eigen::Dynamic,1> Vsizes,Fsizes;
  160. igl::combine(Vlist,Flist,VV,FF,Vsizes,Fsizes);
  161. std::function<int(const int, const int)> keep;
  162. std::function<int(const Eigen::Matrix<int,1,Eigen::Dynamic>) > wind_num_op;
  163. mesh_boolean_type_to_funcs(type,wind_num_op,keep);
  164. return mesh_boolean(VV,FF,Fsizes,wind_num_op,keep,VC,FC,J);
  165. }
  166. template <
  167. typename DerivedVV,
  168. typename DerivedFF,
  169. typename Derivedsizes,
  170. typename DerivedVC,
  171. typename DerivedFC,
  172. typename DerivedJ>
  173. IGL_INLINE bool igl::copyleft::cgal::mesh_boolean(
  174. const Eigen::MatrixBase<DerivedVV > & VV,
  175. const Eigen::MatrixBase<DerivedFF > & FF,
  176. const Eigen::MatrixBase<Derivedsizes> & sizes,
  177. const std::function<int(const Eigen::Matrix<int,1,Eigen::Dynamic>) >& wind_num_op,
  178. const std::function<int(const int, const int)> & keep,
  179. Eigen::PlainObjectBase<DerivedVC > & VC,
  180. Eigen::PlainObjectBase<DerivedFC > & FC,
  181. Eigen::PlainObjectBase<DerivedJ > & J)
  182. {
  183. #ifdef MESH_BOOLEAN_TIMING
  184. const auto & tictoc = []() -> double
  185. {
  186. static double t_start = igl::get_seconds();
  187. double diff = igl::get_seconds()-t_start;
  188. t_start += diff;
  189. return diff;
  190. };
  191. const auto log_time = [&](const std::string& label) -> void {
  192. printf("%50s: %0.5lf\n",
  193. C_STR("mesh_boolean." << label),tictoc());
  194. };
  195. tictoc();
  196. #endif
  197. typedef CGAL::Epeck Kernel;
  198. typedef Kernel::FT ExactScalar;
  199. typedef Eigen::Matrix<typename DerivedJ::Scalar,Eigen::Dynamic,1> VectorXJ;
  200. typedef Eigen::Matrix<
  201. ExactScalar,
  202. Eigen::Dynamic,
  203. Eigen::Dynamic,
  204. DerivedVC::IsRowMajor> MatrixXES;
  205. MatrixXES V;
  206. DerivedFC F;
  207. VectorXJ CJ;
  208. {
  209. Eigen::VectorXi I;
  210. igl::copyleft::cgal::RemeshSelfIntersectionsParam params;
  211. params.stitch_all = true;
  212. MatrixXES Vr;
  213. DerivedFC Fr;
  214. Eigen::MatrixXi IF;
  215. igl::copyleft::cgal::remesh_self_intersections(
  216. VV, FF, params, Vr, Fr, IF, CJ, I);
  217. assert(I.size() == Vr.rows());
  218. // Merge coinciding vertices into non-manifold vertices.
  219. std::for_each(Fr.data(), Fr.data()+Fr.size(),
  220. [&I](typename DerivedFC::Scalar& a) { a=I[a]; });
  221. // Remove unreferenced vertices.
  222. Eigen::VectorXi UIM;
  223. igl::remove_unreferenced(Vr, Fr, V, F, UIM);
  224. }
  225. #ifdef MESH_BOOLEAN_TIMING
  226. log_time("resolve_self_intersection");
  227. #endif
  228. Eigen::MatrixXi E, uE;
  229. Eigen::VectorXi EMAP, uEC, uEE;
  230. igl::unique_edge_map(F, E, uE, EMAP, uEC, uEE);
  231. // Compute patches (F,EMAP,uEC,uEE) --> (P)
  232. Eigen::VectorXi P;
  233. const size_t num_patches = igl::extract_manifold_patches(F,EMAP,uEC,uEE,P);
  234. #ifdef MESH_BOOLEAN_TIMING
  235. log_time("patch_extraction");
  236. #endif
  237. // Compute cells (V,F,P,E,uE,EMAP) -> (per_patch_cells)
  238. Eigen::MatrixXi per_patch_cells;
  239. const size_t num_cells =
  240. extract_cells( V, F, P, uE, EMAP, uEC, uEE, per_patch_cells);
  241. #ifdef MESH_BOOLEAN_TIMING
  242. log_time("cell_extraction");
  243. #endif
  244. // Compute winding numbers on each side of each facet.
  245. const size_t num_faces = F.rows();
  246. // W(f,:) --> [w1out,w1in,w2out,w2in, ... wnout,wnint] winding numbers above
  247. // and below each face w.r.t. each input mesh, so that W(f,2*i) is the
  248. // winding number above face f w.r.t. input i, and W(f,2*i+1) is the winding
  249. // number below face f w.r.t. input i.
  250. Eigen::MatrixXi W;
  251. // labels(f) = i means that face f comes from mesh i
  252. Eigen::VectorXi labels(num_faces);
  253. // cumulative sizes
  254. Derivedsizes cumsizes;
  255. igl::cumsum(sizes,1,cumsizes);
  256. const size_t num_inputs = sizes.size();
  257. std::transform(
  258. CJ.data(),
  259. CJ.data()+CJ.size(),
  260. labels.data(),
  261. // Determine which input mesh birth face i comes from
  262. [&num_inputs,&cumsizes](int i)->int
  263. {
  264. for(int k = 0;k<num_inputs;k++)
  265. {
  266. if(i<cumsizes(k)) return k;
  267. }
  268. assert(false && "Birth parent index out of range");
  269. return -1;
  270. });
  271. bool valid = true;
  272. if (num_faces > 0)
  273. {
  274. valid = valid && propagate_winding_numbers(
  275. V,F,uE,uEC,uEE,num_patches,P,num_cells,per_patch_cells,labels,W);
  276. } else
  277. {
  278. W.resize(0, 2*num_inputs);
  279. }
  280. assert((size_t)W.rows() == num_faces);
  281. // If W doesn't have enough columns, pad with zeros
  282. if (W.cols() <= 2*num_inputs)
  283. {
  284. const int old_ncols = W.cols();
  285. W.conservativeResize(num_faces,2*num_inputs);
  286. W.rightCols(2*num_inputs-old_ncols).setConstant(0);
  287. }
  288. assert((size_t)W.cols() == 2*num_inputs);
  289. #ifdef MESH_BOOLEAN_TIMING
  290. log_time("propagate_input_winding_number");
  291. #endif
  292. // Compute resulting winding number.
  293. Eigen::MatrixXi Wr(num_faces, 2);
  294. for (size_t i=0; i<num_faces; i++)
  295. {
  296. // Winding number vectors above and below
  297. Eigen::RowVectorXi w_out(1,num_inputs), w_in(1,num_inputs);
  298. for(size_t k =0;k<num_inputs;k++)
  299. {
  300. w_out(k) = W(i,2*k+0);
  301. w_in(k) = W(i,2*k+1);
  302. }
  303. Wr(i,0) = wind_num_op(w_out);
  304. Wr(i,1) = wind_num_op(w_in);
  305. }
  306. #ifdef MESH_BOOLEAN_TIMING
  307. log_time("compute_output_winding_number");
  308. #endif
  309. #ifdef SMALL_CELL_REMOVAL
  310. igl::copyleft::cgal::relabel_small_immersed_cells(
  311. V, F, num_patches, P, num_cells, per_patch_cells, 1e-3, Wr);
  312. #endif
  313. // Extract boundary separating inside from outside.
  314. auto index_to_signed_index = [&](size_t i, bool ori) -> int
  315. {
  316. return (i+1)*(ori?1:-1);
  317. };
  318. //auto signed_index_to_index = [&](int i) -> size_t {
  319. // return abs(i) - 1;
  320. //};
  321. std::vector<int> selected;
  322. for(size_t i=0; i<num_faces; i++)
  323. {
  324. auto should_keep = keep(Wr(i,0), Wr(i,1));
  325. if (should_keep > 0)
  326. {
  327. selected.push_back(index_to_signed_index(i, true));
  328. } else if (should_keep < 0)
  329. {
  330. selected.push_back(index_to_signed_index(i, false));
  331. }
  332. }
  333. const size_t num_selected = selected.size();
  334. DerivedFC kept_faces(num_selected, 3);
  335. DerivedJ kept_face_indices(num_selected, 1);
  336. for (size_t i=0; i<num_selected; i++)
  337. {
  338. size_t idx = abs(selected[i]) - 1;
  339. if (selected[i] > 0)
  340. {
  341. kept_faces.row(i) = F.row(idx);
  342. } else
  343. {
  344. kept_faces.row(i) = F.row(idx).reverse();
  345. }
  346. kept_face_indices(i, 0) = CJ[idx];
  347. }
  348. #ifdef MESH_BOOLEAN_TIMING
  349. log_time("extract_output");
  350. #endif
  351. // Finally, remove duplicated faces and unreferenced vertices.
  352. {
  353. DerivedFC G;
  354. DerivedJ JJ;
  355. igl::resolve_duplicated_faces(kept_faces, G, JJ);
  356. J = kept_face_indices(JJ);
  357. #ifdef DOUBLE_CHECK_EXACT_OUTPUT
  358. {
  359. // Sanity check on exact output.
  360. igl::copyleft::cgal::RemeshSelfIntersectionsParam params;
  361. params.detect_only = true;
  362. params.first_only = true;
  363. MatrixXES dummy_VV;
  364. DerivedFC dummy_FF, dummy_IF;
  365. Eigen::VectorXi dummy_J, dummy_IM;
  366. igl::copyleft::cgal::SelfIntersectMesh<
  367. Kernel,
  368. MatrixXES, DerivedFC,
  369. MatrixXES, DerivedFC,
  370. DerivedFC,
  371. Eigen::VectorXi,
  372. Eigen::VectorXi
  373. > checker(V, G, params,
  374. dummy_VV, dummy_FF, dummy_IF, dummy_J, dummy_IM);
  375. if (checker.count != 0)
  376. {
  377. throw "Self-intersection not fully resolved.";
  378. }
  379. }
  380. #endif
  381. // `assign` _after_ removing unreferenced vertices
  382. {
  383. Eigen::VectorXi I,J;
  384. igl::remove_unreferenced(V.rows(),G,I,J);
  385. FC = G;
  386. std::for_each(FC.data(),FC.data()+FC.size(),
  387. [&I](typename DerivedFC::Scalar & a){a=I(a);});
  388. const bool slow_and_more_precise = false;
  389. VC.resize(J.size(),3);
  390. igl::parallel_for(J.size(),[&](Eigen::Index i)
  391. {
  392. for(Eigen::Index c = 0;c<3;c++)
  393. {
  394. assign_scalar(V(J(i),c),slow_and_more_precise,VC(i,c));
  395. }
  396. },1000);
  397. }
  398. }
  399. #ifdef MESH_BOOLEAN_TIMING
  400. log_time("clean_up");
  401. #endif
  402. return valid;
  403. }
  404. template <
  405. typename DerivedVA,
  406. typename DerivedFA,
  407. typename DerivedVB,
  408. typename DerivedFB,
  409. typename DerivedVC,
  410. typename DerivedFC>
  411. IGL_INLINE bool igl::copyleft::cgal::mesh_boolean(
  412. const Eigen::MatrixBase<DerivedVA > & VA,
  413. const Eigen::MatrixBase<DerivedFA > & FA,
  414. const Eigen::MatrixBase<DerivedVB > & VB,
  415. const Eigen::MatrixBase<DerivedFB > & FB,
  416. const MeshBooleanType & type,
  417. Eigen::PlainObjectBase<DerivedVC > & VC,
  418. Eigen::PlainObjectBase<DerivedFC > & FC)
  419. {
  420. Eigen::Matrix<typename DerivedFC::Scalar, Eigen::Dynamic,1> J;
  421. return igl::copyleft::cgal::mesh_boolean(VA,FA,VB,FB,type,VC,FC,J);
  422. }
  423. #ifdef IGL_STATIC_LIBRARY
  424. // Explicit template instantiation
  425. // generated by autoexplicit.sh
  426. // generated by autoexplicit.sh
  427. template bool igl::copyleft::cgal::mesh_boolean<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, 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::MatrixBase<Eigen::Matrix<CGAL::Epeck::FT, -1, -1, 1, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<CGAL::Epeck::FT, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, igl::MeshBooleanType const&, Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Epeck::FT, -1, -1, 1, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  428. // generated by autoexplicit.sh
  429. template bool igl::copyleft::cgal::mesh_boolean<Eigen::Matrix<CGAL::Epeck::FT, -1, 3, 0, -1, 3>, 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<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::MatrixBase<Eigen::Matrix<CGAL::Epeck::FT, -1, 3, 0, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<CGAL::Epeck::FT, -1, 3, 0, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, igl::MeshBooleanType const&, Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Epeck::FT, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  430. template bool igl::copyleft::cgal::mesh_boolean<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<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::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, igl::MeshBooleanType const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  431. template bool igl::copyleft::cgal::mesh_boolean<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<double, -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::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, igl::MeshBooleanType const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  432. #ifdef WIN32
  433. #endif
  434. #endif