mesh_boolean.cpp 17 KB

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