mesh_boolean.cpp 24 KB

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