order_facets_around_edge.cpp 18 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417
  1. #include "order_facets_around_edge.h"
  2. #include <Eigen/Geometry>
  3. #include <CGAL/Exact_predicates_exact_constructions_kernel.h>
  4. #include <stdexcept>
  5. // adj_faces contains signed index starting from +- 1.
  6. template<
  7. typename DerivedV,
  8. typename DerivedF,
  9. typename DerivedI >
  10. void igl::copyleft::cgal::order_facets_around_edge(
  11. const Eigen::PlainObjectBase<DerivedV>& V,
  12. const Eigen::PlainObjectBase<DerivedF>& F,
  13. size_t s,
  14. size_t d,
  15. const std::vector<int>& adj_faces,
  16. Eigen::PlainObjectBase<DerivedI>& order, bool debug)
  17. {
  18. // Although we only need exact predicates in the algorithm,
  19. // exact constructions are needed to avoid degeneracies due to
  20. // casting to double.
  21. typedef CGAL::Exact_predicates_exact_constructions_kernel K;
  22. typedef K::Point_3 Point_3;
  23. typedef K::Plane_3 Plane_3;
  24. auto get_face_index = [&](int adj_f)->size_t
  25. {
  26. return abs(adj_f) - 1;
  27. };
  28. auto get_opposite_vertex = [&](size_t fid)->size_t
  29. {
  30. typedef typename DerivedF::Scalar Index;
  31. if (F(fid, 0) != (Index)s && F(fid, 0) != (Index)d) return F(fid, 0);
  32. if (F(fid, 1) != (Index)s && F(fid, 1) != (Index)d) return F(fid, 1);
  33. if (F(fid, 2) != (Index)s && F(fid, 2) != (Index)d) return F(fid, 2);
  34. assert(false);
  35. return -1;
  36. };
  37. // Handle base cases
  38. if (adj_faces.size() == 0)
  39. {
  40. order.resize(0, 1);
  41. return;
  42. } else if (adj_faces.size() == 1)
  43. {
  44. order.resize(1, 1);
  45. order(0, 0) = 0;
  46. return;
  47. } else if (adj_faces.size() == 2)
  48. {
  49. const size_t o1 = get_opposite_vertex(get_face_index(adj_faces[0]));
  50. const size_t o2 = get_opposite_vertex(get_face_index(adj_faces[1]));
  51. const Point_3 ps(V(s, 0), V(s, 1), V(s, 2));
  52. const Point_3 pd(V(d, 0), V(d, 1), V(d, 2));
  53. const Point_3 p1(V(o1, 0), V(o1, 1), V(o1, 2));
  54. const Point_3 p2(V(o2, 0), V(o2, 1), V(o2, 2));
  55. order.resize(2, 1);
  56. switch (CGAL::orientation(ps, pd, p1, p2))
  57. {
  58. case CGAL::POSITIVE:
  59. order(0, 0) = 1;
  60. order(1, 0) = 0;
  61. break;
  62. case CGAL::NEGATIVE:
  63. order(0, 0) = 0;
  64. order(1, 0) = 1;
  65. break;
  66. case CGAL::COPLANAR:
  67. {
  68. switch (CGAL::coplanar_orientation(ps, pd, p1, p2)) {
  69. case CGAL::POSITIVE:
  70. // Duplicated face, use index to break tie.
  71. order(0, 0) = adj_faces[0] < adj_faces[1] ? 0:1;
  72. order(1, 0) = adj_faces[0] < adj_faces[1] ? 1:0;
  73. break;
  74. case CGAL::NEGATIVE:
  75. // Coplanar faces, one on each side of the edge.
  76. // It is equally valid to order them (0, 1) or (1, 0).
  77. // I cannot think of any reason to prefer one to the
  78. // other. So just use (0, 1) ordering by default.
  79. order(0, 0) = 0;
  80. order(1, 0) = 1;
  81. break;
  82. case CGAL::COLLINEAR:
  83. std::cerr << "Degenerated triangle detected." <<
  84. std::endl;
  85. assert(false);
  86. break;
  87. default:
  88. assert(false);
  89. }
  90. }
  91. break;
  92. default:
  93. assert(false);
  94. }
  95. return;
  96. }
  97. const size_t num_adj_faces = adj_faces.size();
  98. const size_t o = get_opposite_vertex( get_face_index(adj_faces[0]));
  99. const Point_3 p_s(V(s, 0), V(s, 1), V(s, 2));
  100. const Point_3 p_d(V(d, 0), V(d, 1), V(d, 2));
  101. const Point_3 p_o(V(o, 0), V(o, 1), V(o, 2));
  102. const Plane_3 separator(p_s, p_d, p_o);
  103. if (separator.is_degenerate()) {
  104. throw std::runtime_error(
  105. "Cannot order facets around edge due to degenerated facets");
  106. }
  107. std::vector<Point_3> opposite_vertices;
  108. for (size_t i=0; i<num_adj_faces; i++)
  109. {
  110. const size_t o = get_opposite_vertex( get_face_index(adj_faces[i]));
  111. opposite_vertices.emplace_back(
  112. V(o, 0), V(o, 1), V(o, 2));
  113. }
  114. std::vector<int> positive_side;
  115. std::vector<int> negative_side;
  116. std::vector<int> tie_positive_oriented;
  117. std::vector<int> tie_negative_oriented;
  118. std::vector<size_t> positive_side_index;
  119. std::vector<size_t> negative_side_index;
  120. std::vector<size_t> tie_positive_oriented_index;
  121. std::vector<size_t> tie_negative_oriented_index;
  122. for (size_t i=0; i<num_adj_faces; i++)
  123. {
  124. const int f = adj_faces[i];
  125. const Point_3& p_a = opposite_vertices[i];
  126. auto orientation = separator.oriented_side(p_a);
  127. switch (orientation) {
  128. case CGAL::ON_POSITIVE_SIDE:
  129. positive_side.push_back(f);
  130. positive_side_index.push_back(i);
  131. break;
  132. case CGAL::ON_NEGATIVE_SIDE:
  133. negative_side.push_back(f);
  134. negative_side_index.push_back(i);
  135. break;
  136. case CGAL::ON_ORIENTED_BOUNDARY:
  137. {
  138. auto inplane_orientation = CGAL::coplanar_orientation(
  139. p_s, p_d, p_o, p_a);
  140. switch (inplane_orientation) {
  141. case CGAL::POSITIVE:
  142. tie_positive_oriented.push_back(f);
  143. tie_positive_oriented_index.push_back(i);
  144. break;
  145. case CGAL::NEGATIVE:
  146. tie_negative_oriented.push_back(f);
  147. tie_negative_oriented_index.push_back(i);
  148. break;
  149. case CGAL::COLLINEAR:
  150. default:
  151. throw std::runtime_error(
  152. "Degenerated facet detected.");
  153. break;
  154. }
  155. }
  156. break;
  157. default:
  158. // Should not be here.
  159. throw std::runtime_error("Unknown CGAL state detected.");
  160. }
  161. }
  162. if (debug) {
  163. std::cout << "tie positive: " << std::endl;
  164. for (auto& f : tie_positive_oriented) {
  165. std::cout << get_face_index(f) << " ";
  166. }
  167. std::cout << std::endl;
  168. std::cout << "positive side: " << std::endl;
  169. for (auto& f : positive_side) {
  170. std::cout << get_face_index(f) << " ";
  171. }
  172. std::cout << std::endl;
  173. std::cout << "tie negative: " << std::endl;
  174. for (auto& f : tie_negative_oriented) {
  175. std::cout << get_face_index(f) << " ";
  176. }
  177. std::cout << std::endl;
  178. std::cout << "negative side: " << std::endl;
  179. for (auto& f : negative_side) {
  180. std::cout << get_face_index(f) << " ";
  181. }
  182. std::cout << std::endl;
  183. }
  184. auto index_sort = [](std::vector<int>& data) -> std::vector<size_t>{
  185. const size_t len = data.size();
  186. std::vector<size_t> order(len);
  187. for (size_t i=0; i<len; i++) { order[i] = i; }
  188. auto comp = [&](size_t i, size_t j) { return data[i] < data[j]; };
  189. std::sort(order.begin(), order.end(), comp);
  190. return order;
  191. };
  192. DerivedI positive_order, negative_order;
  193. order_facets_around_edge(V, F, s, d, positive_side, positive_order, debug);
  194. order_facets_around_edge(V, F, s, d, negative_side, negative_order, debug);
  195. std::vector<size_t> tie_positive_order = index_sort(tie_positive_oriented);
  196. std::vector<size_t> tie_negative_order = index_sort(tie_negative_oriented);
  197. // Copy results into order vector.
  198. const size_t tie_positive_size = tie_positive_oriented.size();
  199. const size_t tie_negative_size = tie_negative_oriented.size();
  200. const size_t positive_size = positive_order.size();
  201. const size_t negative_size = negative_order.size();
  202. order.resize(
  203. tie_positive_size + positive_size + tie_negative_size + negative_size,1);
  204. size_t count=0;
  205. for (size_t i=0; i<tie_positive_size; i++)
  206. {
  207. order(count+i, 0) = tie_positive_oriented_index[tie_positive_order[i]];
  208. }
  209. count += tie_positive_size;
  210. for (size_t i=0; i<negative_size; i++)
  211. {
  212. order(count+i, 0) = negative_side_index[negative_order(i, 0)];
  213. }
  214. count += negative_size;
  215. for (size_t i=0; i<tie_negative_size; i++)
  216. {
  217. order(count+i, 0) = tie_negative_oriented_index[tie_negative_order[i]];
  218. }
  219. count += tie_negative_size;
  220. for (size_t i=0; i<positive_size; i++)
  221. {
  222. order(count+i, 0) = positive_side_index[positive_order(i, 0)];
  223. }
  224. count += positive_size;
  225. assert(count == num_adj_faces);
  226. // Find the correct start point.
  227. size_t start_idx = 0;
  228. for (size_t i=0; i<num_adj_faces; i++)
  229. {
  230. const Point_3& p_a = opposite_vertices[order(i, 0)];
  231. const Point_3& p_b =
  232. opposite_vertices[order((i+1)%num_adj_faces, 0)];
  233. auto orientation = CGAL::orientation(p_s, p_d, p_a, p_b);
  234. if (orientation == CGAL::POSITIVE)
  235. {
  236. // Angle between triangle (p_s, p_d, p_a) and (p_s, p_d, p_b) is
  237. // more than 180 degrees.
  238. start_idx = (i+1)%num_adj_faces;
  239. break;
  240. } else if (orientation == CGAL::COPLANAR &&
  241. Plane_3(p_s, p_d, p_a).orthogonal_direction() !=
  242. Plane_3(p_s, p_d, p_b).orthogonal_direction())
  243. {
  244. // All 4 points are coplanar, but p_a and p_b are on each side of
  245. // the edge (p_s, p_d). This means the angle between triangle
  246. // (p_s, p_d, p_a) and (p_s, p_d, p_b) is exactly 180 degrees.
  247. start_idx = (i+1)%num_adj_faces;
  248. break;
  249. }
  250. }
  251. DerivedI circular_order = order;
  252. for (size_t i=0; i<num_adj_faces; i++)
  253. {
  254. order(i, 0) = circular_order((start_idx + i)%num_adj_faces, 0);
  255. }
  256. }
  257. template<
  258. typename DerivedV,
  259. typename DerivedF,
  260. typename DerivedI>
  261. IGL_INLINE
  262. void igl::copyleft::cgal::order_facets_around_edge(
  263. const Eigen::PlainObjectBase<DerivedV>& V,
  264. const Eigen::PlainObjectBase<DerivedF>& F,
  265. size_t s,
  266. size_t d,
  267. const std::vector<int>& adj_faces,
  268. const Eigen::PlainObjectBase<DerivedV>& pivot_point,
  269. Eigen::PlainObjectBase<DerivedI>& order)
  270. {
  271. assert(V.cols() == 3);
  272. assert(F.cols() == 3);
  273. assert(pivot_point.cols() == 3);
  274. auto signed_index_to_index = [&](int signed_idx)
  275. {
  276. return abs(signed_idx) -1;
  277. };
  278. auto get_opposite_vertex_index = [&](size_t fid) -> typename DerivedF::Scalar
  279. {
  280. typedef typename DerivedF::Scalar Index;
  281. if (F(fid, 0) != (Index)s && F(fid, 0) != (Index)d) return F(fid, 0);
  282. if (F(fid, 1) != (Index)s && F(fid, 1) != (Index)d) return F(fid, 1);
  283. if (F(fid, 2) != (Index)s && F(fid, 2) != (Index)d) return F(fid, 2);
  284. assert(false);
  285. // avoid warning
  286. return -1;
  287. };
  288. {
  289. // Check if s, d and pivot are collinear.
  290. typedef CGAL::Exact_predicates_exact_constructions_kernel K;
  291. K::Point_3 ps(V(s,0), V(s,1), V(s,2));
  292. K::Point_3 pd(V(d,0), V(d,1), V(d,2));
  293. K::Point_3 pp(pivot_point(0,0), pivot_point(0,1), pivot_point(0,2));
  294. if (CGAL::collinear(ps, pd, pp)) {
  295. throw std::runtime_error(
  296. "Pivot point is collinear with the outer edge!");
  297. }
  298. }
  299. const size_t N = adj_faces.size();
  300. const size_t num_faces = N + 1; // N adj faces + 1 pivot face
  301. // Because face indices are used for tie breaking, the original face indices
  302. // in the new faces array must be ascending.
  303. auto comp = [&](int i, int j)
  304. {
  305. return signed_index_to_index(adj_faces[i]) <
  306. signed_index_to_index(adj_faces[j]);
  307. };
  308. std::vector<size_t> adj_order(N);
  309. for (size_t i=0; i<N; i++) adj_order[i] = i;
  310. std::sort(adj_order.begin(), adj_order.end(), comp);
  311. DerivedV vertices(num_faces + 2, 3);
  312. for (size_t i=0; i<N; i++)
  313. {
  314. const size_t fid = signed_index_to_index(adj_faces[adj_order[i]]);
  315. vertices.row(i) = V.row(get_opposite_vertex_index(fid));
  316. }
  317. vertices.row(N ) = pivot_point;
  318. vertices.row(N+1) = V.row(s);
  319. vertices.row(N+2) = V.row(d);
  320. DerivedF faces(num_faces, 3);
  321. for (size_t i=0; i<N; i++)
  322. {
  323. if (adj_faces[adj_order[i]] < 0)
  324. {
  325. faces(i,0) = N+1; // s
  326. faces(i,1) = N+2; // d
  327. faces(i,2) = i ;
  328. } else
  329. {
  330. faces(i,0) = N+2; // d
  331. faces(i,1) = N+1; // s
  332. faces(i,2) = i ;
  333. }
  334. }
  335. // Last face is the pivot face.
  336. faces(N, 0) = N+1;
  337. faces(N, 1) = N+2;
  338. faces(N, 2) = N;
  339. std::vector<int> adj_faces_with_pivot(num_faces);
  340. for (size_t i=0; i<num_faces; i++)
  341. {
  342. if ((size_t)faces(i,0) == N+1 && (size_t)faces(i,1) == N+2)
  343. {
  344. adj_faces_with_pivot[i] = int(i+1) * -1;
  345. } else
  346. {
  347. adj_faces_with_pivot[i] = int(i+1);
  348. }
  349. }
  350. DerivedI order_with_pivot;
  351. order_facets_around_edge(
  352. vertices, faces, N+1, N+2, adj_faces_with_pivot, order_with_pivot);
  353. assert((size_t)order_with_pivot.size() == num_faces);
  354. order.resize(N);
  355. size_t pivot_index = num_faces + 1;
  356. for (size_t i=0; i<num_faces; i++)
  357. {
  358. if ((size_t)order_with_pivot[i] == N)
  359. {
  360. pivot_index = i;
  361. break;
  362. }
  363. }
  364. assert(pivot_index < num_faces);
  365. for (size_t i=0; i<N; i++)
  366. {
  367. order[i] = adj_order[order_with_pivot[(pivot_index+i+1)%num_faces]];
  368. }
  369. }
  370. #ifdef IGL_STATIC_LIBRARY
  371. // Explicit template instantiation
  372. // generated by autoexplicit.sh
  373. template void igl::copyleft::cgal::order_facets_around_edge<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, -1, 1, -1, -1>, Eigen::Matrix<int, -1, 3, 1, -1, 3>, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, -1, 1, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 1, -1, 3> > const&, unsigned long, unsigned long, std::vector<int, std::allocator<int> > const&, Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, -1, 1, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  374. // generated by autoexplicit.sh
  375. template void igl::copyleft::cgal::order_facets_around_edge<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, -1, 1, -1, -1>, Eigen::Matrix<int, -1, 3, 1, -1, 3>, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, -1, 1, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 1, -1, 3> > const&, unsigned long, unsigned long, std::vector<int, std::allocator<int> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&, bool);
  376. // generated by autoexplicit.sh
  377. template void igl::copyleft::cgal::order_facets_around_edge<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, -1, 1, -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>, -1, -1, 1, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, unsigned long, unsigned long, std::vector<int, std::allocator<int> > const&, Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, -1, 1, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  378. // generated by autoexplicit.sh
  379. template void igl::copyleft::cgal::order_facets_around_edge<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, -1, 1, -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>, -1, -1, 1, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, unsigned long, unsigned long, std::vector<int, std::allocator<int> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&, bool);
  380. template void igl::copyleft::cgal::order_facets_around_edge<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::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&, unsigned long, unsigned long, std::vector<int, std::allocator<int> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&, bool);
  381. template void igl::copyleft::cgal::order_facets_around_edge<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::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, unsigned long, unsigned long, std::vector<int, std::allocator<int> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&, bool);
  382. template void igl::copyleft::cgal::order_facets_around_edge<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::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&, unsigned long, unsigned long, std::vector<int, std::allocator<int> > 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> >&);
  383. template void igl::copyleft::cgal::order_facets_around_edge<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::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, unsigned long, unsigned long, std::vector<int, std::allocator<int> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  384. template void igl::copyleft::cgal::order_facets_around_edge<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&, unsigned long, unsigned long, std::vector<int, std::allocator<int> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  385. template void igl::copyleft::cgal::order_facets_around_edge<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::PlainObjectBase<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, unsigned long, unsigned long, std::vector<int, std::allocator<int> > const&, Eigen::PlainObjectBase<Eigen::Matrix<CGAL::Lazy_exact_nt<CGAL::Gmpq>, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> >&);
  386. #endif