dual_contouring.cpp 29 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575
  1. #include "dual_contouring.h"
  2. #include "quadprog.h"
  3. #include "parallel_for.h"
  4. #include <thread>
  5. #include <mutex>
  6. #include <vector>
  7. #include <unordered_map>
  8. #include <iostream>
  9. #include <cstdint>
  10. namespace igl
  11. {
  12. // These classes not intended to be used directly
  13. class Hash
  14. {
  15. public:
  16. // https://stackoverflow.com/a/26348708/148668
  17. std::uint64_t operator()(const std::tuple<int, int, int> & key) const
  18. {
  19. // Check that conversion is safe. Could use int16_t directly everywhere
  20. // below but it's an uncommon type to expose and grid indices should
  21. // never be more than 2¹⁵-1 in the first place.
  22. assert( std::get<0>(key) == (int)(std::int16_t)std::get<0>(key));
  23. assert( std::get<1>(key) == (int)(std::int16_t)std::get<1>(key));
  24. assert( std::get<2>(key) == (int)(std::int16_t)std::get<2>(key));
  25. std::uint64_t result = std::uint16_t(std::get<0>(key));
  26. result = (result << 16) + std::uint16_t(std::get<1>(key));
  27. result = (result << 16) + std::uint16_t(std::get<2>(key));
  28. return result;
  29. };
  30. };
  31. template <typename Scalar>
  32. class DualContouring
  33. {
  34. // Types
  35. public:
  36. using RowVector3S = Eigen::Matrix<Scalar,1,3>;
  37. using RowVector4S = Eigen::Matrix<Scalar,1,4>;
  38. using Matrix4S = Eigen::Matrix<Scalar,4,4>;
  39. using Matrix3S = Eigen::Matrix<Scalar,3,3>;
  40. using Vector3S = Eigen::Matrix<Scalar,3,1>;
  41. using KeyTriplet = std::tuple<int,int,int>;
  42. public:
  43. // Working variables
  44. // see dual_contouring.h
  45. // f(x) returns >0 outside, <0 inside, and =0 on the surface
  46. std::function<Scalar(const RowVector3S &)> f;
  47. // f_grad(x) returns (df/dx)/‖df/dx‖ (normalization only important when
  48. // f(x) = 0).
  49. std::function<RowVector3S(const RowVector3S &)> f_grad;
  50. bool constrained;
  51. bool triangles;
  52. bool root_finding;
  53. RowVector3S min_corner;
  54. RowVector3S step;
  55. Eigen::Matrix<Scalar,Eigen::Dynamic,3> V;
  56. // Internal variables
  57. // Running number of vertices added during contouring
  58. typename decltype(V)::Index n;
  59. // map from cell subscript to index in V
  60. std::unordered_map< KeyTriplet, typename decltype(V)::Index, Hash > C2V;
  61. // running list of aggregate vertex positions (used for spring
  62. // regularization term)
  63. std::vector<RowVector3S,Eigen::aligned_allocator<RowVector3S>> vV;
  64. // running list of subscripts corresponding to vertices
  65. std::vector<Eigen::RowVector3i,Eigen::aligned_allocator<Eigen::RowVector3i>> vI;
  66. // running list of quadric matrices corresponding to inserted vertices
  67. std::vector<Matrix4S,Eigen::aligned_allocator<Matrix4S>> vH;
  68. // running list of number of faces incident on this vertex (used to
  69. // normalize spring regulatization term)
  70. std::vector<int> vcount;
  71. // running list of output quad faces
  72. Eigen::Matrix<Eigen::Index,Eigen::Dynamic,Eigen::Dynamic> Q;
  73. // running number of real quads in Q (used for dynamic array allocation)
  74. typename decltype(Q)::Index m;
  75. // mutexes used to insert into Q and (vV,vI,vH,vcount)
  76. std::mutex Qmut;
  77. std::mutex Vmut;
  78. public:
  79. DualContouring(
  80. const std::function<Scalar(const RowVector3S &)> & _f,
  81. const std::function<RowVector3S(const RowVector3S &)> & _f_grad,
  82. const bool _constrained = false,
  83. const bool _triangles = false,
  84. const bool _root_finding = true):
  85. f(_f),
  86. f_grad(_f_grad),
  87. constrained(_constrained),
  88. triangles(_triangles),
  89. root_finding(_root_finding),
  90. n(0),
  91. C2V(0),
  92. vV(),vI(),vH(),vcount(),
  93. m(0)
  94. { }
  95. // Side effects: new entry in vV,vI,vH,vcount, increment n
  96. // Returns index of new vertex
  97. typename decltype(V)::Index new_vertex()
  98. {
  99. const auto v = n;
  100. n++;
  101. vcount.resize(n);
  102. vV.resize(n);
  103. vI.resize(n);
  104. vH.resize(n);
  105. vcount[v] = 0;
  106. vV[v].setZero();
  107. vH[v].setZero();
  108. return v;
  109. };
  110. // Inputs:
  111. // kc 3-long vector of {x,y,z} index of primal grid **cell**
  112. // Returns index to corresponding dual vertex
  113. // Side effects: if vertex for this cell does not yet exist, creates it
  114. typename decltype(V)::Index sub2dual(const Eigen::RowVector3i & kc)
  115. {
  116. const KeyTriplet key = {kc(0),kc(1),kc(2)};
  117. const auto it = C2V.find(key);
  118. typename decltype(V)::Index v = -1;
  119. if(it == C2V.end())
  120. {
  121. v = new_vertex();
  122. C2V[key] = v;
  123. vI[v] = kc;
  124. }else
  125. {
  126. v = it->second;
  127. }
  128. return v;
  129. };
  130. RowVector3S primal(const Eigen::RowVector3i & ic) const
  131. {
  132. return min_corner + (ic.cast<Scalar>().array() * step.array()).matrix();
  133. }
  134. Eigen::RowVector3i inverse_primal(const RowVector3S & x) const
  135. {
  136. // x = min_corner + (ic.cast<Scalar>().array() * step.array()).matrix();
  137. // x-min_corner = (ic.cast<Scalar>().array() * step.array())
  138. // (x-min_corner).array() / step.array() = ic.cast<Scalar>().array()
  139. // ((x-min_corner).array() / step.array()).round() = ic
  140. return
  141. ((x-min_corner).array()/step.array()).round().template cast<int>();
  142. }
  143. // Inputs:
  144. // x x-index of vertex on primal grid
  145. // y y-index of vertex on primal grid
  146. // z z-index of vertex on primal grid
  147. // o which edge are we looking back on? o=0->x,o=1->y,o=2->z
  148. // Side effects: may insert new vertices into vV,vI,vH,vcount, new faces
  149. // into Q
  150. bool single_edge(const int & x, const int & y, const int & z, const int & o)
  151. {
  152. const RowVector3S e0 = primal(Eigen::RowVector3i(x,y,z));
  153. const Scalar f0 = f(e0);
  154. return single_edge(x,y,z,o,e0,f0);
  155. }
  156. bool single_edge(
  157. const int & x,
  158. const int & y,
  159. const int & z,
  160. const int & o,
  161. const RowVector3S & e0,
  162. const Scalar & f0)
  163. {
  164. //e1 computed here needs to precisely agree with e0 when called with
  165. //correspond x,y,z. So, don't do this:
  166. //Eigen::RowVector3d e1 = e0;
  167. //e1(o) -= step(o);
  168. Eigen::RowVector3i jc(x,y,z);
  169. jc(o) -= 1;
  170. const RowVector3S e1 = primal(jc);
  171. const Scalar f1 = f(e1);
  172. return single_edge(x,y,z,o,e0,f0,e1,f1);
  173. }
  174. bool single_edge(
  175. const int & x,
  176. const int & y,
  177. const int & z,
  178. const int & o,
  179. const RowVector3S & e0,
  180. const Scalar & f0,
  181. const RowVector3S & e1,
  182. const Scalar & f1)
  183. {
  184. const Scalar isovalue = 0;
  185. if((f0>isovalue) == (f1>isovalue)) { return false; }
  186. // Position of crossing point along edge
  187. RowVector3S p;
  188. Scalar t = -1;
  189. if(root_finding)
  190. {
  191. Scalar tl = 0;
  192. bool gl = f0>0;
  193. Scalar tu = 1;
  194. bool gu = f1>0;
  195. assert(gu ^ gl);
  196. int riter = 0;
  197. const int max_riter = 7;
  198. while(true)
  199. {
  200. t = 0.5*(tu + tl);
  201. p = e0+t*(e1-e0);
  202. riter++;
  203. if(riter > max_riter) { break;}
  204. const Scalar ft = f(p);
  205. if( (ft>0) == gu) { tu = t; }
  206. else if( (ft>0) == gl){ tl = t; }
  207. else { break; }
  208. }
  209. }else
  210. {
  211. // inverse lerp
  212. const Scalar delta = f1-f0;
  213. if(delta == 0) { t = 0.5; }
  214. t = (isovalue - f0)/delta;
  215. p = e0+t*(e1-e0);
  216. }
  217. typename decltype(V)::Index ev;
  218. {
  219. const std::lock_guard<std::mutex> lock(Vmut);
  220. // insert vertex at this point to triangulate quad face
  221. ev = triangles ? new_vertex() : -1;
  222. if (triangles)
  223. {
  224. vV[ev] = p;
  225. vcount[ev] = 1;
  226. vI[ev] = Eigen::RowVector3i(-1, -1, -1);
  227. }
  228. }
  229. // edge normal from function handle (could use grid finite
  230. // differences/interpolation gradients)
  231. const RowVector3S dfdx = f_grad(p);
  232. // homogenous plane equation
  233. const RowVector4S P = (RowVector4S()<<dfdx,-dfdx.dot(p)).finished();
  234. // quadric contribution
  235. const Matrix4S H = P.transpose() * P;
  236. // Build quad face from dual vertices of 4 cells around this edge
  237. Eigen::RowVector4i face;
  238. int k = 0;
  239. for(int i = -1;i<=0;i++)
  240. {
  241. for(int j = -1;j<=0;j++)
  242. {
  243. Eigen::RowVector3i kc(x,y,z);
  244. kc(o)--;
  245. kc((o+1)%3)+=i;
  246. kc((o+2)%3)+=j;
  247. const std::lock_guard<std::mutex> lock(Vmut);
  248. const typename decltype(V)::Index v = sub2dual(kc);
  249. vV[v] += p;
  250. vcount[v]++;
  251. vH[v] += H;
  252. face(k++) = v;
  253. }
  254. }
  255. {
  256. const std::lock_guard<std::mutex> lock(Qmut);
  257. if(triangles)
  258. {
  259. if(m+4 >= Q.rows()){ Q.conservativeResize(2*m+4,Q.cols()); }
  260. if(f0>f1)
  261. {
  262. Q.row(m+0)<< ev,face(3),face(1) ;
  263. Q.row(m+1)<< ev,face(1),face(0);
  264. Q.row(m+2)<< face(2), ev,face(0);
  265. Q.row(m+3)<< face(2),face(3), ev;
  266. }else
  267. {
  268. Q.row(m+0)<< ev,face(1),face(3) ;
  269. Q.row(m+1)<< ev,face(3),face(2);
  270. Q.row(m+2)<< face(0), ev,face(2);
  271. Q.row(m+3)<< face(0),face(1), ev;
  272. }
  273. m+=4;
  274. }else
  275. {
  276. if(m+1 >= Q.rows()){ Q.conservativeResize(2*m+1,Q.cols()); }
  277. if(f0>f1)
  278. {
  279. Q.row(m)<< face(2),face(3),face(1),face(0);
  280. }else
  281. {
  282. Q.row(m)<< face(0),face(1),face(3),face(2);
  283. }
  284. m++;
  285. }
  286. }
  287. return true;
  288. }
  289. // Side effects: Q resized to fit m, V constructed to fit n and
  290. // reconstruct data in vH,vI,vV,vcount
  291. void dual_vertex_positions()
  292. {
  293. Q.conservativeResize(m,Q.cols());
  294. V.resize(n,3);
  295. igl::parallel_for(n,[&](const Eigen::Index v)
  296. {
  297. RowVector3S mid = vV[v] / Scalar(vcount[v]);
  298. if(triangles && vI[v](0)<0 ){ V.row(v) = mid; return; }
  299. const Scalar w = 1e-2*(0.01+vcount[v]);
  300. Matrix3S A = vH[v].block(0,0,3,3) + w*Matrix3S::Identity();
  301. RowVector3S b = -vH[v].block(3,0,1,3) + w*mid;
  302. // Replace with solver
  303. //RowVector3S p = b * A.inverse();
  304. //
  305. // min_p ½ pᵀ A p - pᵀb
  306. //
  307. // let p = p₀ + x
  308. //
  309. // min ½ (p₀ + x )ᵀ A (p₀ + x ) - (p₀ + x )ᵀb
  310. // step≥x≥0
  311. const RowVector3S p0 =
  312. min_corner + ((vI[v].template cast<Scalar>().array()) * step.array()).matrix();
  313. const RowVector3S x =
  314. constrained ?
  315. igl::quadprog<Scalar,3>(A,(p0*A-b).transpose(),Vector3S(0,0,0),step.transpose()) :
  316. Eigen::LLT<Matrix3S>(A).solve(-(p0*A-b).transpose());
  317. V.row(v) = p0+x;
  318. },1000ul);
  319. }
  320. // Inputs:
  321. // _min_corner minimum (bottomLeftBack) corner of primal grid
  322. // max_corner maximum (topRightFront) corner of primal grid
  323. // nx number of primal grid vertices along x-axis
  324. // ny number of primal grid vertices along y-ayis
  325. // nz number of primal grid vertices along z-azis
  326. // Side effects: prepares vV,vI,vH,vcount, Q for vertex_positions()
  327. void dense(
  328. const RowVector3S & _min_corner,
  329. const RowVector3S & max_corner,
  330. const int nx,
  331. const int ny,
  332. const int nz)
  333. {
  334. min_corner = _min_corner;
  335. step =
  336. (max_corner-min_corner).array()/(RowVector3S(nx,ny,nz).array()-1);
  337. // Should do some reasonable reserves for C2V,vV,vI,vH,vcount
  338. Q.resize(std::pow(nx*ny*nz,2./3.),triangles?3:4);
  339. // loop over grid
  340. igl::parallel_for(nx,[&](const int x)
  341. {
  342. for(int y = 0;y<ny;y++)
  343. {
  344. for(int z = 0;z<nz;z++)
  345. {
  346. const RowVector3S e0 = primal(Eigen::RowVector3i(x,y,z));
  347. const Scalar f0 = f(e0);
  348. // we'll consider the edges going "back" from this vertex
  349. for(int o = 0;o<3;o++)
  350. {
  351. // off-by-one boundary cases
  352. if(((o==0)&&x==0)||((o==1)&&y==0)||((o==2)&&z==0)){ continue;}
  353. single_edge(x,y,z,o,e0,f0);
  354. }
  355. }
  356. }
  357. },10ul);
  358. dual_vertex_positions();
  359. }
  360. template <typename DerivedGf, typename DerivedGV>
  361. void dense(
  362. const Eigen::MatrixBase<DerivedGf> & Gf,
  363. const Eigen::MatrixBase<DerivedGV> & GV,
  364. const int nx,
  365. const int ny,
  366. const int nz)
  367. {
  368. min_corner = GV.colwise().minCoeff();
  369. const RowVector3S max_corner = GV.colwise().maxCoeff();
  370. step =
  371. (max_corner-min_corner).array()/(RowVector3S(nx,ny,nz).array()-1);
  372. // Should do some reasonable reserves for C2V,vV,vI,vH,vcount
  373. Q.resize(std::pow(nx*ny*nz,2./3.),triangles?3:4);
  374. const auto xyz2i = [&nx,&ny]
  375. (const int & x, const int & y, const int & z)->Eigen::Index
  376. {
  377. return x+nx*(y+ny*(z));
  378. };
  379. // loop over grid
  380. igl::parallel_for(nz,[&](const int z)
  381. {
  382. for(int y = 0;y<ny;y++)
  383. {
  384. for(int x = 0;x<nx;x++)
  385. {
  386. //const Scalar f0 = f(e0);
  387. const Eigen::Index k0 = xyz2i(x,y,z);
  388. const RowVector3S e0 = GV.row(k0);
  389. const Scalar f0 = Gf(k0);
  390. // we'll consider the edges going "back" from this vertex
  391. for(int o = 0;o<3;o++)
  392. {
  393. Eigen::RowVector3i jc(x,y,z);
  394. jc(o) -= 1;
  395. if(jc(o)<0) { continue;} // off-by-one boundary cases
  396. const int k1 = xyz2i(jc(0),jc(1),jc(2));
  397. const RowVector3S e1 = GV.row(k1);
  398. const Scalar f1 = Gf(k1);
  399. single_edge(x,y,z,o,e0,f0,e1,f1);
  400. }
  401. }
  402. }
  403. },10ul);
  404. dual_vertex_positions();
  405. }
  406. void sparse(
  407. const RowVector3S & _step,
  408. const Eigen::Matrix<Scalar,Eigen::Dynamic,1> & Gf,
  409. const Eigen::Matrix<Scalar,Eigen::Dynamic,3> & GV,
  410. const Eigen::Matrix<int,Eigen::Dynamic,2> & GI)
  411. {
  412. step = _step;
  413. Q.resize((triangles?4:1)*GI.rows(),triangles?3:4);
  414. // in perfect world doesn't matter where min_corner is so long as it is
  415. // _on_ the grid. For models very far from origin, centering grid near
  416. // model avoids possible rounding error in hash()/inverse_primal()
  417. // [still very unlikely, but let's be safe]
  418. min_corner = GV.colwise().minCoeff();
  419. // igl::parallel_for here made things worse. Probably need to do proper
  420. // map-reduce rather than locks on mutexes.
  421. for(Eigen::Index i = 0;i<GI.rows();i++)
  422. {
  423. RowVector3S e0 = GV.row(GI(i,0));
  424. RowVector3S e1 = GV.row(GI(i,1));
  425. Eigen::RowVector3i ic0 = inverse_primal(e0);
  426. Eigen::RowVector3i ic1 = inverse_primal(e1);
  427. #ifndef NDEBUG
  428. RowVector3S p0 = primal(ic0);
  429. RowVector3S p1 = primal(ic1);
  430. assert( (p0-e0).norm() < 1e-10);
  431. assert( (p1-e1).norm() < 1e-10);
  432. #endif
  433. Scalar f0 = Gf(GI(i,0)); //f(e0);
  434. Scalar f1 = Gf(GI(i,1)); //f(e1);
  435. // should differ in just one coordinate. Find that coordinate.
  436. int o = -1;
  437. for(int j = 0;j<3;j++)
  438. {
  439. if(ic0(j) == ic1(j)){ continue;}
  440. if(ic0(j) - ic1(j) == 1)
  441. {
  442. assert(o == -1 && "Edges should differ in just one coordinate");
  443. o = j;
  444. continue; // rather than break so assertions fire
  445. }
  446. if(ic1(j) - ic0(j) == 1)
  447. {
  448. assert(o == -1 && "Edges should differ in just one coordinate");
  449. std::swap(e0,e1);
  450. std::swap(f0,f1);
  451. std::swap(ic0,ic1);
  452. o = j;
  453. continue; // rather than break so assertions fire
  454. } else
  455. {
  456. assert(false && "Edges should differ in just one coordinate");
  457. }
  458. }
  459. assert(o>=0 && "Edges should differ in just one coordinate");
  460. // i0 is the larger subscript location and ic1 is backward in the o
  461. // direction.
  462. for(int j = 0;j<3;j++){ assert(ic0(j) == ic1(j)+(o==j)); }
  463. const int x = ic0(0);
  464. const int y = ic0(1);
  465. const int z = ic0(2);
  466. single_edge(x,y,z,o,e0,f0,e1,f1);
  467. }
  468. dual_vertex_positions();
  469. }
  470. };
  471. }
  472. template <
  473. typename DerivedV,
  474. typename DerivedQ>
  475. IGL_INLINE void igl::dual_contouring(
  476. const std::function<
  477. typename DerivedV::Scalar(const Eigen::Matrix<typename DerivedV::Scalar,1,3> &)> & f,
  478. const std::function<
  479. Eigen::Matrix<typename DerivedV::Scalar,1,3>(
  480. const Eigen::Matrix<typename DerivedV::Scalar,1,3> &)> & f_grad,
  481. const Eigen::Matrix<typename DerivedV::Scalar,1,3> & min_corner,
  482. const Eigen::Matrix<typename DerivedV::Scalar,1,3> & max_corner,
  483. const int nx,
  484. const int ny,
  485. const int nz,
  486. const bool constrained,
  487. const bool triangles,
  488. const bool root_finding,
  489. Eigen::PlainObjectBase<DerivedV> & V,
  490. Eigen::PlainObjectBase<DerivedQ> & Q)
  491. {
  492. typedef typename DerivedV::Scalar Scalar;
  493. DualContouring<Scalar> DC(f,f_grad,constrained,triangles,root_finding);
  494. DC.dense(min_corner,max_corner,nx,ny,nz);
  495. V = DC.V;
  496. Q = DC.Q.template cast<typename DerivedQ::Scalar>();
  497. }
  498. template <
  499. typename DerivedGf,
  500. typename DerivedGV,
  501. typename DerivedV,
  502. typename DerivedQ>
  503. IGL_INLINE void igl::dual_contouring(
  504. const std::function<
  505. typename DerivedV::Scalar(const Eigen::Matrix<typename DerivedV::Scalar,1,3> &)> & f,
  506. const std::function<
  507. Eigen::Matrix<typename DerivedV::Scalar,1,3>(
  508. const Eigen::Matrix<typename DerivedV::Scalar,1,3> &)> & f_grad,
  509. const Eigen::MatrixBase<DerivedGf> & Gf,
  510. const Eigen::MatrixBase<DerivedGV> & GV,
  511. const int nx,
  512. const int ny,
  513. const int nz,
  514. const bool constrained,
  515. const bool triangles,
  516. const bool root_finding,
  517. Eigen::PlainObjectBase<DerivedV> & V,
  518. Eigen::PlainObjectBase<DerivedQ> & Q)
  519. {
  520. typedef typename DerivedV::Scalar Scalar;
  521. DualContouring<Scalar> DC(f,f_grad,constrained,triangles,root_finding);
  522. DC.dense(Gf,GV,nx,ny,nz);
  523. V = DC.V;
  524. Q = DC.Q.template cast<typename DerivedQ::Scalar>();
  525. }
  526. template <
  527. typename DerivedGf,
  528. typename DerivedGV,
  529. typename DerivedGI,
  530. typename DerivedV,
  531. typename DerivedQ>
  532. IGL_INLINE void igl::dual_contouring(
  533. const std::function<typename DerivedV::Scalar(const Eigen::Matrix<typename DerivedV::Scalar,1,3> &)> & f,
  534. const std::function<Eigen::Matrix<typename DerivedV::Scalar,1,3>(const Eigen::Matrix<typename DerivedV::Scalar,1,3> &)> & f_grad,
  535. const Eigen::Matrix<typename DerivedV::Scalar,1,3> & step,
  536. const Eigen::MatrixBase<DerivedGf> & Gf,
  537. const Eigen::MatrixBase<DerivedGV> & GV,
  538. const Eigen::MatrixBase<DerivedGI> & GI,
  539. const bool constrained,
  540. const bool triangles,
  541. const bool root_finding,
  542. Eigen::PlainObjectBase<DerivedV> & V,
  543. Eigen::PlainObjectBase<DerivedQ> & Q)
  544. {
  545. if(GI.rows() == 0){ return;}
  546. assert(GI.cols() == 2);
  547. typedef typename DerivedV::Scalar Scalar;
  548. DualContouring<Scalar> DC(f,f_grad,constrained,triangles,root_finding);
  549. DC.sparse(step,Gf,GV,GI);
  550. V = DC.V;
  551. Q = DC.Q.template cast<typename DerivedQ::Scalar>();
  552. }
  553. #ifdef IGL_STATIC_LIBRARY
  554. // Explicit template instantiation
  555. // generated by autoexplicit.sh
  556. template void igl::dual_contouring<Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 3, 1, -1, 3>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(std::function<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar (Eigen::Matrix<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&)> const&, std::function<Eigen::Matrix<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> (Eigen::Matrix<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&)> const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 3, 1, -1, 3> > const&, int, int, int, bool, bool, bool, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  557. template void igl::dual_contouring<Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(std::function<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar (Eigen::Matrix<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&)> const&, std::function<Eigen::Matrix<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> (Eigen::Matrix<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&)> const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, int, int, int, bool, bool, bool, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  558. template void igl::dual_contouring<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(std::function<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar (Eigen::Matrix<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&)> const&, std::function<Eigen::Matrix<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> (Eigen::Matrix<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&)> const&, Eigen::Matrix<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&, Eigen::Matrix<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&, int, int, int, bool, bool, bool, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  559. template void igl::dual_contouring<Eigen::Matrix<double, -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> >(std::function<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar (Eigen::Matrix<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&)> const&, std::function<Eigen::Matrix<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> (Eigen::Matrix<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&)> const&, Eigen::Matrix<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, bool, bool, bool, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  560. template void igl::dual_contouring<Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 2, 0, -1, 2>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(std::function<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar (Eigen::Matrix<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&)> const&, std::function<Eigen::Matrix<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> (Eigen::Matrix<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&)> const&, Eigen::Matrix<Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 2, 0, -1, 2> > const&, bool, bool, bool, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  561. template void igl::dual_contouring<Eigen::Matrix<float, -1, 1, 0, -1, 1>, Eigen::Matrix<float, -1, 3, 1, -1, 3>, Eigen::Matrix<float, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(std::function<Eigen::Matrix<float, -1, -1, 0, -1, -1>::Scalar (Eigen::Matrix<Eigen::Matrix<float, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&)> const&, std::function<Eigen::Matrix<Eigen::Matrix<float, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> (Eigen::Matrix<Eigen::Matrix<float, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&)> const&, Eigen::MatrixBase<Eigen::Matrix<float, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<float, -1, 3, 1, -1, 3> > const&, int, int, int, bool, bool, bool, Eigen::PlainObjectBase<Eigen::Matrix<float, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  562. template void igl::dual_contouring<Eigen::Matrix<float, -1, 1, 0, -1, 1>, Eigen::Matrix<float, -1, -1, 0, -1, -1>, Eigen::Matrix<float, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(std::function<Eigen::Matrix<float, -1, -1, 0, -1, -1>::Scalar (Eigen::Matrix<Eigen::Matrix<float, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&)> const&, std::function<Eigen::Matrix<Eigen::Matrix<float, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> (Eigen::Matrix<Eigen::Matrix<float, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&)> const&, Eigen::MatrixBase<Eigen::Matrix<float, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<float, -1, -1, 0, -1, -1> > const&, int, int, int, bool, bool, bool, Eigen::PlainObjectBase<Eigen::Matrix<float, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  563. template void igl::dual_contouring<Eigen::Matrix<float, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(std::function<Eigen::Matrix<float, -1, -1, 0, -1, -1>::Scalar (Eigen::Matrix<Eigen::Matrix<float, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&)> const&, std::function<Eigen::Matrix<Eigen::Matrix<float, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> (Eigen::Matrix<Eigen::Matrix<float, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&)> const&, Eigen::Matrix<Eigen::Matrix<float, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&, Eigen::Matrix<Eigen::Matrix<float, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&, int, int, int, bool, bool, bool, Eigen::PlainObjectBase<Eigen::Matrix<float, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  564. template void igl::dual_contouring<Eigen::Matrix<float, -1, 1, 0, -1, 1>, Eigen::Matrix<float, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<float, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(std::function<Eigen::Matrix<float, -1, -1, 0, -1, -1>::Scalar (Eigen::Matrix<Eigen::Matrix<float, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&)> const&, std::function<Eigen::Matrix<Eigen::Matrix<float, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> (Eigen::Matrix<Eigen::Matrix<float, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&)> const&, Eigen::Matrix<Eigen::Matrix<float, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&, Eigen::MatrixBase<Eigen::Matrix<float, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<float, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, bool, bool, bool, Eigen::PlainObjectBase<Eigen::Matrix<float, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  565. template void igl::dual_contouring<Eigen::Matrix<float, -1, 1, 0, -1, 1>, Eigen::Matrix<float, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 2, 0, -1, 2>, Eigen::Matrix<float, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(std::function<Eigen::Matrix<float, -1, -1, 0, -1, -1>::Scalar (Eigen::Matrix<Eigen::Matrix<float, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&)> const&, std::function<Eigen::Matrix<Eigen::Matrix<float, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> (Eigen::Matrix<Eigen::Matrix<float, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&)> const&, Eigen::Matrix<Eigen::Matrix<float, -1, -1, 0, -1, -1>::Scalar, 1, 3, 1, 1, 3> const&, Eigen::MatrixBase<Eigen::Matrix<float, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<float, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 2, 0, -1, 2> > const&, bool, bool, bool, Eigen::PlainObjectBase<Eigen::Matrix<float, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> >&);
  566. #endif