remesh_intersections.cpp 56 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2015 Qingnan Zhou <qnzhou@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. //
  9. #include "remesh_intersections.h"
  10. #include "assign_scalar.h"
  11. #include "../../get_seconds.h"
  12. #include "../../unique.h"
  13. #include <vector>
  14. #include <map>
  15. #include <queue>
  16. #include <unordered_map>
  17. #include <iostream>
  18. //#define REMESH_INTERSECTIONS_TIMING
  19. template <
  20. typename DerivedV,
  21. typename DerivedF,
  22. typename Kernel,
  23. typename DerivedVV,
  24. typename DerivedFF,
  25. typename DerivedJ,
  26. typename DerivedIM>
  27. IGL_INLINE void igl::copyleft::cgal::remesh_intersections(
  28. const Eigen::PlainObjectBase<DerivedV> & V,
  29. const Eigen::PlainObjectBase<DerivedF> & F,
  30. const std::vector<CGAL::Triangle_3<Kernel> > & T,
  31. const std::map<
  32. typename DerivedF::Index,
  33. std::vector<
  34. std::pair<typename DerivedF::Index, CGAL::Object> > > & offending,
  35. const std::map<
  36. std::pair<typename DerivedF::Index,typename DerivedF::Index>,
  37. std::vector<typename DerivedF::Index> > & /*edge2faces*/,
  38. Eigen::PlainObjectBase<DerivedVV> & VV,
  39. Eigen::PlainObjectBase<DerivedFF> & FF,
  40. Eigen::PlainObjectBase<DerivedJ> & J,
  41. Eigen::PlainObjectBase<DerivedIM> & IM) {
  42. #ifdef REMESH_INTERSECTIONS_TIMING
  43. const auto & tictoc = []() -> double
  44. {
  45. static double t_start = igl::get_seconds();
  46. double diff = igl::get_seconds()-t_start;
  47. t_start += diff;
  48. return diff;
  49. };
  50. const auto log_time = [&](const std::string& label) -> void {
  51. std::cout << "remesh_intersections." << label << ": "
  52. << tictoc() << std::endl;
  53. };
  54. tictoc();
  55. #endif
  56. typedef CGAL::Point_3<Kernel> Point_3;
  57. typedef CGAL::Segment_3<Kernel> Segment_3;
  58. typedef CGAL::Triangle_3<Kernel> Triangle_3;
  59. typedef CGAL::Plane_3<Kernel> Plane_3;
  60. typedef CGAL::Triangulation_vertex_base_2<Kernel> TVB_2;
  61. typedef CGAL::Constrained_triangulation_face_base_2<Kernel> CTFB_2;
  62. typedef CGAL::Triangulation_data_structure_2<TVB_2,CTFB_2> TDS_2;
  63. typedef CGAL::Exact_intersections_tag Itag;
  64. typedef CGAL::Constrained_Delaunay_triangulation_2<Kernel,TDS_2,Itag>
  65. CDT_2;
  66. typedef CGAL::Constrained_triangulation_plus_2<CDT_2> CDT_plus_2;
  67. typedef typename DerivedF::Index Index;
  68. typedef std::pair<Index, Index> Edge;
  69. struct EdgeHash {
  70. size_t operator()(const Edge& e) const {
  71. return (e.first * 805306457) ^ (e.second * 201326611);
  72. }
  73. };
  74. typedef std::unordered_map<Edge, std::vector<Index>, EdgeHash > EdgeMap;
  75. const size_t num_faces = F.rows();
  76. const size_t num_base_vertices = V.rows();
  77. assert(num_faces == T.size());
  78. std::vector<bool> is_offending(num_faces, false);
  79. for (const auto itr : offending) {
  80. const auto& fid = itr.first;
  81. is_offending[fid] = true;
  82. }
  83. std::unordered_map<Index, std::vector<Index> > intersecting_and_coplanar;
  84. for (const auto itr : offending) {
  85. const auto& fi = itr.first;
  86. const auto P = T[fi].supporting_plane();
  87. assert(!P.is_degenerate());
  88. for (const auto jtr : itr.second) {
  89. const auto& fj = jtr.first;
  90. const auto& tj = T[fj];
  91. if (P.has_on(tj[0]) && P.has_on(tj[1]) && P.has_on(tj[2])) {
  92. auto loc = intersecting_and_coplanar.find(fi);
  93. if (loc == intersecting_and_coplanar.end()) {
  94. intersecting_and_coplanar[fi] = {fj};
  95. } else {
  96. loc->second.push_back(fj);
  97. }
  98. }
  99. }
  100. }
  101. #ifdef REMESH_INTERSECTIONS_TIMING
  102. log_time("overlap_analysis");
  103. #endif
  104. std::vector<std::vector<Index> > resolved_faces;
  105. std::vector<Index> source_faces;
  106. std::vector<Point_3> new_vertices;
  107. EdgeMap edge_vertices;
  108. /**
  109. * Run constraint Delaunay triangulation on the plane.
  110. */
  111. auto run_delaunay_triangulation = [&offending, &T](const Plane_3& P,
  112. const std::vector<Index>& involved_faces,
  113. std::vector<Point_3>& vertices,
  114. std::vector<std::vector<Index> >& faces) -> void {
  115. CDT_plus_2 cdt;
  116. for (const auto& fid : involved_faces) {
  117. const auto itr = offending.find(fid);
  118. const auto& triangle = T[fid];
  119. cdt.insert_constraint(P.to_2d(triangle[0]), P.to_2d(triangle[1]));
  120. cdt.insert_constraint(P.to_2d(triangle[1]), P.to_2d(triangle[2]));
  121. cdt.insert_constraint(P.to_2d(triangle[2]), P.to_2d(triangle[0]));
  122. if (itr == offending.end()) continue;
  123. for (const auto& index_obj : itr->second) {
  124. const auto& ofid = index_obj.first;
  125. const auto& obj = index_obj.second;
  126. if(const Segment_3 *iseg = CGAL::object_cast<Segment_3 >(&obj)) {
  127. // Add segment constraint
  128. cdt.insert_constraint(
  129. P.to_2d(iseg->vertex(0)),P.to_2d(iseg->vertex(1)));
  130. }else if(const Point_3 *ipoint = CGAL::object_cast<Point_3 >(&obj)) {
  131. // Add point
  132. cdt.insert(P.to_2d(*ipoint));
  133. } else if(const Triangle_3 *itri = CGAL::object_cast<Triangle_3 >(&obj)) {
  134. // Add 3 segment constraints
  135. cdt.insert_constraint(
  136. P.to_2d(itri->vertex(0)),P.to_2d(itri->vertex(1)));
  137. cdt.insert_constraint(
  138. P.to_2d(itri->vertex(1)),P.to_2d(itri->vertex(2)));
  139. cdt.insert_constraint(
  140. P.to_2d(itri->vertex(2)),P.to_2d(itri->vertex(0)));
  141. } else if(const std::vector<Point_3 > *polyp =
  142. CGAL::object_cast< std::vector<Point_3 > >(&obj)) {
  143. //cerr<<REDRUM("Poly...")<<endl;
  144. const std::vector<Point_3 > & poly = *polyp;
  145. const Index m = poly.size();
  146. assert(m>=2);
  147. for(Index p = 0;p<m;p++)
  148. {
  149. const Index np = (p+1)%m;
  150. cdt.insert_constraint(P.to_2d(poly[p]),P.to_2d(poly[np]));
  151. }
  152. }else {
  153. throw std::runtime_error("Unknown intersection object!");
  154. }
  155. }
  156. }
  157. std::map<typename CDT_plus_2::Vertex_handle,Index> v2i;
  158. size_t count=0;
  159. for (auto itr = cdt.finite_vertices_begin();
  160. itr != cdt.finite_vertices_end(); itr++) {
  161. vertices.push_back(P.to_3d(itr->point()));
  162. v2i[itr] = count;
  163. count++;
  164. }
  165. for (auto itr = cdt.finite_faces_begin();
  166. itr != cdt.finite_faces_end(); itr++) {
  167. faces.push_back( {
  168. v2i[itr->vertex(0)],
  169. v2i[itr->vertex(1)],
  170. v2i[itr->vertex(2)] });
  171. }
  172. };
  173. /**
  174. * Given p on triangle indexed by ori_f, determine the index of p.
  175. */
  176. auto determine_point_index = [&](
  177. const Point_3& p, const size_t ori_f) -> Index {
  178. const auto& triangle = T[ori_f];
  179. const auto& f = F.row(ori_f).eval();
  180. // Check if p is one of the triangle corners.
  181. for (size_t i=0; i<3; i++) {
  182. if (p == triangle[i]) return f[i];
  183. }
  184. // Check if p is on one of the edges.
  185. for (size_t i=0; i<3; i++) {
  186. const Point_3 curr_corner = triangle[i];
  187. const Point_3 next_corner = triangle[(i+1)%3];
  188. const Segment_3 edge(curr_corner, next_corner);
  189. if (edge.has_on(p)) {
  190. const Index curr = f[i];
  191. const Index next = f[(i+1)%3];
  192. Edge key;
  193. key.first = curr<next?curr:next;
  194. key.second = curr<next?next:curr;
  195. auto itr = edge_vertices.find(key);
  196. if (itr == edge_vertices.end()) {
  197. const Index index =
  198. num_base_vertices + new_vertices.size();
  199. edge_vertices.insert({key, {index}});
  200. new_vertices.push_back(p);
  201. return index;
  202. } else {
  203. for (const auto vid : itr->second) {
  204. if (p == new_vertices[vid - num_base_vertices]) {
  205. return vid;
  206. }
  207. }
  208. const size_t index = num_base_vertices + new_vertices.size();
  209. new_vertices.push_back(p);
  210. itr->second.push_back(index);
  211. return index;
  212. }
  213. }
  214. }
  215. // p must be in the middle of the triangle.
  216. const size_t index = num_base_vertices + new_vertices.size();
  217. new_vertices.push_back(p);
  218. return index;
  219. };
  220. /**
  221. * Determine the vertex indices for each corner of each output triangle.
  222. */
  223. auto post_triangulation_process = [&](
  224. const std::vector<Point_3>& vertices,
  225. const std::vector<std::vector<Index> >& faces,
  226. const std::vector<Index>& involved_faces) -> void {
  227. for (const auto& f : faces) {
  228. const Point_3& v0 = vertices[f[0]];
  229. const Point_3& v1 = vertices[f[1]];
  230. const Point_3& v2 = vertices[f[2]];
  231. Point_3 center(
  232. (v0[0] + v1[0] + v2[0]) / 3.0,
  233. (v0[1] + v1[1] + v2[1]) / 3.0,
  234. (v0[2] + v1[2] + v2[2]) / 3.0);
  235. for (const auto& ori_f : involved_faces) {
  236. const auto& triangle = T[ori_f];
  237. const Plane_3 P = triangle.supporting_plane();
  238. if (triangle.has_on(center)) {
  239. std::vector<Index> corners(3);
  240. corners[0] = determine_point_index(v0, ori_f);
  241. corners[1] = determine_point_index(v1, ori_f);
  242. corners[2] = determine_point_index(v2, ori_f);
  243. if (CGAL::orientation(
  244. P.to_2d(v0), P.to_2d(v1), P.to_2d(v2))
  245. == CGAL::RIGHT_TURN) {
  246. std::swap(corners[0], corners[1]);
  247. }
  248. resolved_faces.emplace_back(corners);
  249. source_faces.push_back(ori_f);
  250. }
  251. }
  252. }
  253. };
  254. // Process un-touched faces.
  255. for (size_t i=0; i<num_faces; i++) {
  256. if (!is_offending[i] && !T[i].is_degenerate()) {
  257. resolved_faces.push_back(
  258. { F(i,0), F(i,1), F(i,2) } );
  259. source_faces.push_back(i);
  260. }
  261. }
  262. // Process self-intersecting faces.
  263. std::vector<bool> processed(num_faces, false);
  264. std::vector<std::pair<Plane_3, std::vector<Index> > > cdt_inputs;
  265. for (const auto itr : offending) {
  266. const auto fid = itr.first;
  267. if (processed[fid]) continue;
  268. processed[fid] = true;
  269. const auto loc = intersecting_and_coplanar.find(fid);
  270. std::vector<Index> involved_faces;
  271. if (loc == intersecting_and_coplanar.end()) {
  272. involved_faces.push_back(fid);
  273. } else {
  274. std::queue<Index> Q;
  275. Q.push(fid);
  276. while (!Q.empty()) {
  277. const auto index = Q.front();
  278. involved_faces.push_back(index);
  279. Q.pop();
  280. const auto overlapping_faces =
  281. intersecting_and_coplanar.find(index);
  282. assert(overlapping_faces != intersecting_and_coplanar.end());
  283. for (const auto other_index : overlapping_faces->second) {
  284. if (processed[other_index]) continue;
  285. processed[other_index] = true;
  286. Q.push(other_index);
  287. }
  288. }
  289. }
  290. Plane_3 P = T[fid].supporting_plane();
  291. cdt_inputs.emplace_back(P, involved_faces);
  292. }
  293. #ifdef REMESH_INTERSECTIONS_TIMING
  294. log_time("preprocess");
  295. #endif
  296. const size_t num_cdts = cdt_inputs.size();
  297. std::vector<std::vector<Point_3> > cdt_vertices(num_cdts);
  298. std::vector<std::vector<std::vector<Index> > > cdt_faces(num_cdts);
  299. const auto run_cdt = [&](const size_t first, const size_t last) {
  300. for (size_t i=first; i<last; i++) {
  301. auto& vertices = cdt_vertices[i];
  302. auto& faces = cdt_faces[i];
  303. const auto& P = cdt_inputs[i].first;
  304. const auto& involved_faces = cdt_inputs[i].second;
  305. run_delaunay_triangulation(P, involved_faces, vertices, faces);
  306. }
  307. };
  308. run_cdt(0, num_cdts);
  309. #ifdef REMESH_INTERSECTIONS_TIMING
  310. log_time("cdt");
  311. #endif
  312. // Post process
  313. for (size_t i=0; i<num_cdts; i++) {
  314. const auto& vertices = cdt_vertices[i];
  315. const auto& faces = cdt_faces[i];
  316. const auto& involved_faces = cdt_inputs[i].second;
  317. post_triangulation_process(vertices, faces, involved_faces);
  318. }
  319. #ifdef REMESH_INTERSECTIONS_TIMING
  320. log_time("stitching");
  321. #endif
  322. // Output resolved mesh.
  323. const size_t num_out_vertices = new_vertices.size() + num_base_vertices;
  324. VV.resize(num_out_vertices, 3);
  325. for (size_t i=0; i<num_base_vertices; i++) {
  326. assign_scalar(V(i,0), VV(i,0));
  327. assign_scalar(V(i,1), VV(i,1));
  328. assign_scalar(V(i,2), VV(i,2));
  329. }
  330. for (size_t i=num_base_vertices; i<num_out_vertices; i++) {
  331. assign_scalar(new_vertices[i-num_base_vertices][0], VV(i,0));
  332. assign_scalar(new_vertices[i-num_base_vertices][1], VV(i,1));
  333. assign_scalar(new_vertices[i-num_base_vertices][2], VV(i,2));
  334. }
  335. const size_t num_out_faces = resolved_faces.size();
  336. FF.resize(num_out_faces, 3);
  337. for (size_t i=0; i<num_out_faces; i++) {
  338. FF(i,0) = resolved_faces[i][0];
  339. FF(i,1) = resolved_faces[i][1];
  340. FF(i,2) = resolved_faces[i][2];
  341. }
  342. J.resize(num_out_faces);
  343. std::copy(source_faces.begin(), source_faces.end(), J.data());
  344. struct PointHash {
  345. size_t operator()(const Point_3& p) const {
  346. static const std::hash<double> hasher{};
  347. double x,y,z;
  348. assign_scalar(p.x(), x);
  349. assign_scalar(p.y(), y);
  350. assign_scalar(p.z(), z);
  351. return hasher(x) ^ hasher(y) ^ hasher(z);
  352. }
  353. };
  354. // TODO: use igl::unique()
  355. // Extract unique vertex indices.
  356. const size_t VV_size = VV.rows();
  357. IM.resize(VV_size,1);
  358. DerivedVV unique_vv;
  359. Eigen::VectorXi unique_to_vv, vv_to_unique;
  360. igl::unique_rows(VV, unique_vv, unique_to_vv, vv_to_unique);
  361. for (Index v=0; v<VV_size; v++) {
  362. IM(v) = unique_to_vv[vv_to_unique[v]];
  363. }
  364. #ifdef REMESH_INTERSECTIONS_TIMING
  365. log_time("store_results");
  366. #endif
  367. }
  368. #ifdef IGL_STATIC_LIBRARY
  369. template void igl::copyleft::cgal::remesh_intersections<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, CGAL::Epeck, Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<long, -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, 3, 0, -1, 3> > const&, std::__1::vector<CGAL::Triangle_3<CGAL::Epeck>, std::__1::allocator<CGAL::Triangle_3<CGAL::Epeck> > > const&, std::__1::map<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object> > >, std::__1::less<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index const, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index,CGAL::Object> > > > > > const&, std::__1::map<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index>, std::__1::vector<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::allocator<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> >, std::__1::less<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> >, std::__1::allocator<std::__1::pair<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> const, std::__1::vector<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::allocator<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> > > > > const&, Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<long, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  370. template void igl::copyleft::cgal::remesh_intersections<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, CGAL::Epick, Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<long, -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, 3, 0, -1, 3> > const&, std::__1::vector<CGAL::Triangle_3<CGAL::Epick>, std::__1::allocator<CGAL::Triangle_3<CGAL::Epick> > > const&, std::__1::map<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object> > >, std::__1::less<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index const, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object> > > > > > const&, std::__1::map<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index>, std::__1::vector<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::allocator<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> >, std::__1::less<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> >, std::__1::allocator<std::__1::pair<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> const, std::__1::vector<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::allocator<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> > > > > const&, Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<long, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int , -1, 1, 0, -1, 1> >&);
  371. template void igl::copyleft::cgal::remesh_intersections<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, CGAL::Epeck, Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<long, -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&, std::__1::vector<CGAL::Triangle_3<CGAL::Epeck>, std::__1::allocator<CGAL::Triangle_3<CGAL::Epeck> > > const&, std::__1::map<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object> > >, std::__1::less<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index const, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object> > > > > > const&, std::__1::map<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index>, std::__1::vector<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, std::__1::allocator<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> >, std::__1::less<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> >, std::__1::allocator<std::__1::pair<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> const, std::__1::vector<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, std::__1::allocator<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> > > > > const&, Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<long, -1, 1, 0, -1, 1 > >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  372. template void igl::copyleft::cgal::remesh_intersections<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, CGAL::Epick, Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<long, -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&, std::__1::vector<CGAL::Triangle_3<CGAL::Epick>, std::__1::allocator<CGAL::Triangle_3<CGAL::Epick> > > const&, std::__1::map<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object> > >, std::__1::less<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index const, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object> > > > > > const&, std::__1::map<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index>, std::__1::vector<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, std::__1::allocator<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> >, std::__1::less<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> >, std::__1::allocator<std::__1::pair<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> const, std::__1::vector<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, std::__1::allocator<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> > > > > const&, Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<long, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  373. template void igl::copyleft::cgal::remesh_intersections<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, CGAL::Epeck, Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -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<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, std::__1::vector<CGAL::Triangle_3<CGAL::Epeck>, std::__1::allocator<CGAL::Triangle_3<CGAL::Epeck> > > const&, std::__1::map<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object> > >, std::__1::less<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index const, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object> > > > > > const&, std::__1::map<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index>, std::__1::vector<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::allocator<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> >, std::__1::less<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> >, std::__1::allocator<std::__1::pair<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0,-1, 3>::Index> const, std::__1::vector<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::allocator<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> > > > > const&, Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1,3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -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> >&);
  374. template void igl::copyleft::cgal::remesh_intersections<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, CGAL::Epeck, Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<long, -1, 1, 0, -1, 1>, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, std::__1::vector<CGAL::Triangle_3<CGAL::Epeck>, std::__1::allocator<CGAL::Triangle_3<CGAL::Epeck> > > const&, std::__1::map<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object> > >, std::__1::less<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index const, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object> > > > > > const&, std::__1::map<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index>, std::__1::vector<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::allocator<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> >, std::__1::less<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> >, std::__1::allocator<std::__1::pair<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> const, std::__1::vector<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::allocator<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> > > > > const&, Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<long, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  375. template void igl::copyleft::cgal::remesh_intersections<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, CGAL::Epeck, Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<long, -1, 1, 0, -1, 1>, Eigen::Matrix<long, -1, 1, 0, -1, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, std::__1::vector<CGAL::Triangle_3<CGAL::Epeck>, std::__1::allocator<CGAL::Triangle_3<CGAL::Epeck> > > const&, std::__1::map<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object> > >, std::__1::less<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1,3>::Index const, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object> > > > > > const&, std::__1::map<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index>, std::__1::vector<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::allocator<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> >, std::__1::less<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> >, std::__1::allocator<std::__1::pair<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> const, std::__1::vector<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::allocator<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> > > > > const&, Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<long, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<long, -1, 1, 0, -1, 1> >&);
  376. template void igl::copyleft::cgal::remesh_intersections<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, CGAL::Epeck, Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -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<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, std::__1::vector<CGAL::Triangle_3<CGAL::Epeck>, std::__1::allocator<CGAL::Triangle_3<CGAL::Epeck> > > const&, std::__1::map<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object> > >, std::__1::less<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index const, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object> > > > > > const&, std::__1::map<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index>, std::__1::vector<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::allocator<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> >, std::__1::less<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> >, std::__1::allocator<std::__1::pair<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> const, std::__1::vector<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::allocator<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> > > > > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -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> >&);
  377. template void igl::copyleft::cgal::remesh_intersections<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, CGAL::Epeck, Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<long, -1, 1, 0, -1, 1>, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, std::__1::vector<CGAL::Triangle_3<CGAL::Epeck>, std::__1::allocator<CGAL::Triangle_3<CGAL::Epeck> > > const&, std::__1::map<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object> > >, std::__1::less<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index const, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object> > > > > > const&, std::__1::map<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index>, std::__1::vector<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::allocator<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> >, std::__1::less<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> >, std::__1::allocator<std::__1::pair<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> const, std::__1::vector<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::allocator<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> > > > > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<long, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  378. template void igl::copyleft::cgal::remesh_intersections<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, CGAL::Epick, Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -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<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, std::__1::vector<CGAL::Triangle_3<CGAL::Epick>, std::__1::allocator<CGAL::Triangle_3<CGAL::Epick> > > const&, std::__1::map<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object> > >, std::__1::less<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index const, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object> > > > > > const&, std::__1::map<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index>, std::__1::vector<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::allocator<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> >, std::__1::less<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> >, std::__1::allocator<std::__1::pair<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> const, std::__1::vector<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::allocator<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> > > > > const&, Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -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> >&);
  379. template void igl::copyleft::cgal::remesh_intersections<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, CGAL::Epick, Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<long, -1, 1, 0, -1, 1>, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, std::__1::vector<CGAL::Triangle_3<CGAL::Epick>, std::__1::allocator<CGAL::Triangle_3<CGAL::Epick> > > const&, std::__1::map<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object> > >, std::__1::less<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1,3>::Index const, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object> > > > > > const&, std::__1::map<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index>, std::__1::vector<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::allocator<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> >, std::__1::less<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> >, std::__1::allocator<std::__1::pair<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> const, std::__1::vector<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::allocator<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> > > > > const&, Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<long, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  380. template void igl::copyleft::cgal::remesh_intersections<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, CGAL::Epick, Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<long, -1, 1, 0, -1, 1>, Eigen::Matrix<long, -1, 1, 0, -1, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, std::__1::vector<CGAL::Triangle_3<CGAL::Epick>, std::__1::allocator<CGAL::Triangle_3<CGAL::Epick> > > const&, std::__1::map<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object> > >, std::__1::less<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index const, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object> > > > > > const&, std::__1::map<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index>, std::__1::vector<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::allocator<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> >, std::__1::less<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> >, std::__1::allocator<std::__1::pair<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> const, std::__1::vector<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::allocator<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> > > > > const&, Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<long, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<long, -1, 1, 0, -1, 1> >&);
  381. template void igl::copyleft::cgal::remesh_intersections<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, CGAL::Epick, Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -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<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, std::__1::vector<CGAL::Triangle_3<CGAL::Epick>, std::__1::allocator<CGAL::Triangle_3<CGAL::Epick> > > const&, std::__1::map<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object> > >, std::__1::less<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index const, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object> > > > > > const&, std::__1::map<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index>, std::__1::vector<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::allocator<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> >, std::__1::less<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> >, std::__1::allocator<std::__1::pair<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> const, std::__1::vector<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::allocator<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> > > > > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -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> >&);
  382. template void igl::copyleft::cgal::remesh_intersections<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, CGAL::Epick, Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<long, -1, 1, 0, -1, 1>, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, std::__1::vector<CGAL::Triangle_3<CGAL::Epick>, std::__1::allocator<CGAL::Triangle_3<CGAL::Epick> > > const&, std::__1::map<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object> > >, std::__1::less<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index const, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, CGAL::Object> > > > > > const&, std::__1::map<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index>, std::__1::vector<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::allocator<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> >, std::__1::less<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> >, std::__1::allocator<std::__1::pair<std::__1::pair<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> const, std::__1::vector<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index, std::__1::allocator<Eigen::Matrix<int, -1, 3, 0, -1, 3>::Index> > > > > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> >&, Eigen::PlainObjectBase<Eigen::Matrix<long, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  383. template void igl::copyleft::cgal::remesh_intersections<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, CGAL::Epeck, Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, -1, 0, -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> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, std::__1::vector<CGAL::Triangle_3<CGAL::Epeck>, std::__1::allocator<CGAL::Triangle_3<CGAL::Epeck> > > const&, std::__1::map<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object> > >, std::__1::less<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index const, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object> > > > > > const&, std::__1::map<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index>, std::__1::vector<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, std::__1::allocator<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> >, std::__1::less<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> >, std::__1::allocator<std::__1::pair<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> const, std::__1::vector<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, std::__1::allocator<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> > > > > const&, Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -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> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  384. template void igl::copyleft::cgal::remesh_intersections<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, CGAL::Epeck, Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<long, -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&, std::__1::vector<CGAL::Triangle_3<CGAL::Epeck>, std::__1::allocator<CGAL::Triangle_3<CGAL::Epeck> > > const&, std::__1::map<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object> > >, std::__1::less<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index const, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object> > > > > > const&, std::__1::map<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index>, std::__1::vector<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, std::__1::allocator<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> >, std::__1::less<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> >, std::__1::allocator<std::__1::pair<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> const, std::__1::vector<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, std::__1::allocator<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> > > > > const&, Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<long, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  385. template void igl::copyleft::cgal::remesh_intersections<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, CGAL::Epeck, 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<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&, std::__1::vector<CGAL::Triangle_3<CGAL::Epeck>, std::__1::allocator<CGAL::Triangle_3<CGAL::Epeck> > > const&, std::__1::map<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object> > >, std::__1::less<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index const, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object> > > > > > const&, std::__1::map<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index>, std::__1::vector<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, std::__1::allocator<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> >, std::__1::less<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> >, std::__1::allocator<std::__1::pair<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> const, std::__1::vector<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, std::__1::allocator<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> > > > > 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> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  386. template void igl::copyleft::cgal::remesh_intersections<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, CGAL::Epick, Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, -1, 0, -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> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, std::__1::vector<CGAL::Triangle_3<CGAL::Epick>, std::__1::allocator<CGAL::Triangle_3<CGAL::Epick> > > const&, std::__1::map<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object> > >, std::__1::less<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index const, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object> > > > > > const&, std::__1::map<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index>, std::__1::vector<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, std::__1::allocator<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> >, std::__1::less<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> >, std::__1::allocator<std::__1::pair<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> const, std::__1::vector<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, std::__1::allocator<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> > > > > const&, Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -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> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  387. template void igl::copyleft::cgal::remesh_intersections<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, CGAL::Epick, Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<long, -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&, std::__1::vector<CGAL::Triangle_3<CGAL::Epick>, std::__1::allocator<CGAL::Triangle_3<CGAL::Epick> > > const&, std::__1::map<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object> > >, std::__1::less<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index const, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object> > > > > > const&, std::__1::map<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index>, std::__1::vector<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, std::__1::allocator<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> >, std::__1::less<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> >, std::__1::allocator<std::__1::pair<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> const, std::__1::vector<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, std::__1::allocator<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> > > > > const&, Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<long, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  388. template void igl::copyleft::cgal::remesh_intersections<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, CGAL::Epick, 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<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&, std::__1::vector<CGAL::Triangle_3<CGAL::Epick>, std::__1::allocator<CGAL::Triangle_3<CGAL::Epick> > > const&, std::__1::map<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object> > >, std::__1::less<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index const, std::__1::vector<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object>, std::__1::allocator<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, CGAL::Object> > > > > > const&, std::__1::map<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index>, std::__1::vector<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, std::__1::allocator<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> >, std::__1::less<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> >, std::__1::allocator<std::__1::pair<std::__1::pair<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> const, std::__1::vector<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index, std::__1::allocator<Eigen::Matrix<int, -1, -1, 0, -1, -1>::Index> > > > > 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> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  389. #endif