nrosy.cpp 22 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2014 Daniele Panozzo <daniele.panozzo@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. #include <igl/comiso/nrosy.h>
  9. #include <igl/triangle_triangle_adjacency.h>
  10. #include <igl/edge_topology.h>
  11. #include <igl/per_face_normals.h>
  12. #include <iostream>
  13. #include <fstream>
  14. #include <Eigen/Geometry>
  15. #include <Eigen/Sparse>
  16. #include <queue>
  17. #include <gmm/gmm.h>
  18. #include <CoMISo/Solver/ConstrainedSolver.hh>
  19. #include <CoMISo/Solver/MISolver.hh>
  20. #include <CoMISo/Solver/GMM_Tools.hh>
  21. namespace igl
  22. {
  23. class NRosyField
  24. {
  25. public:
  26. // Init
  27. IGL_INLINE NRosyField(const Eigen::MatrixXd& _V, const Eigen::MatrixXi& _F);
  28. // Generate the N-rosy field
  29. // N degree of the rosy field
  30. // roundseparately: round the integer variables one at a time, slower but higher quality
  31. IGL_INLINE void solve(const int N = 4);
  32. // Set a hard constraint on fid
  33. // fid: face id
  34. // v: direction to fix (in 3d)
  35. IGL_INLINE void setConstraintHard(const int fid, const Eigen::Vector3d& v);
  36. // Set a soft constraint on fid
  37. // fid: face id
  38. // w: weight of the soft constraint, clipped between 0 and 1
  39. // v: direction to fix (in 3d)
  40. IGL_INLINE void setConstraintSoft(const int fid, const double w, const Eigen::Vector3d& v);
  41. // Set the ratio between smoothness and soft constraints (0 -> smoothness only, 1 -> soft constr only)
  42. IGL_INLINE void setSoftAlpha(double alpha);
  43. // Reset constraints (at least one constraint must be present or solve will fail)
  44. IGL_INLINE void resetConstraints();
  45. // Return the current field
  46. IGL_INLINE Eigen::MatrixXd getFieldPerFace();
  47. // Return the current field (in Ahish's ffield format)
  48. IGL_INLINE Eigen::MatrixXd getFFieldPerFace();
  49. // Compute singularity indexes
  50. IGL_INLINE void findCones(int N);
  51. // Return the singularities
  52. IGL_INLINE Eigen::VectorXd getSingularityIndexPerVertex();
  53. private:
  54. // Compute angle differences between reference frames
  55. IGL_INLINE void computek();
  56. // Remove useless matchings
  57. IGL_INLINE void reduceSpace();
  58. // Prepare the system matrix
  59. IGL_INLINE void prepareSystemMatrix(const int N);
  60. // Solve without roundings
  61. IGL_INLINE void solveNoRoundings();
  62. // Solve with roundings using CoMIso
  63. IGL_INLINE void solveRoundings();
  64. // Round all p to 0 and fix
  65. IGL_INLINE void roundAndFixToZero();
  66. // Round all p and fix
  67. IGL_INLINE void roundAndFix();
  68. // Convert a vector in 3d to an angle wrt the local reference system
  69. IGL_INLINE double convert3DtoLocal(unsigned fid, const Eigen::Vector3d& v);
  70. // Convert an angle wrt the local reference system to a 3d vector
  71. IGL_INLINE Eigen::Vector3d convertLocalto3D(unsigned fid, double a);
  72. // Compute the per vertex angle defect
  73. IGL_INLINE Eigen::VectorXd angleDefect();
  74. // Temporary variable for the field
  75. Eigen::VectorXd angles;
  76. // Hard constraints
  77. Eigen::VectorXd hard;
  78. std::vector<bool> isHard;
  79. // Soft constraints
  80. Eigen::VectorXd soft;
  81. Eigen::VectorXd wSoft;
  82. double softAlpha;
  83. // Face Topology
  84. Eigen::MatrixXi TT, TTi;
  85. // Edge Topology
  86. Eigen::MatrixXi EV, FE, EF;
  87. std::vector<bool> isBorderEdge;
  88. // Per Edge information
  89. // Angle between two reference frames
  90. Eigen::VectorXd k;
  91. // Jumps
  92. Eigen::VectorXi p;
  93. std::vector<bool> pFixed;
  94. // Mesh
  95. Eigen::MatrixXd V;
  96. Eigen::MatrixXi F;
  97. // Normals per face
  98. Eigen::MatrixXd N;
  99. // Singularity index
  100. Eigen::VectorXd singularityIndex;
  101. // Reference frame per triangle
  102. std::vector<Eigen::MatrixXd> TPs;
  103. // System stuff
  104. Eigen::SparseMatrix<double> A;
  105. Eigen::VectorXd b;
  106. Eigen::VectorXi tag_t;
  107. Eigen::VectorXi tag_p;
  108. };
  109. } // NAMESPACE IGL
  110. igl::NRosyField::NRosyField(const Eigen::MatrixXd& _V, const Eigen::MatrixXi& _F)
  111. {
  112. using namespace std;
  113. using namespace Eigen;
  114. V = _V;
  115. F = _F;
  116. assert(V.rows() > 0);
  117. assert(F.rows() > 0);
  118. // Generate topological relations
  119. igl::triangle_triangle_adjacency(V,F,TT,TTi);
  120. igl::edge_topology(V,F, EV, FE, EF);
  121. // Flag border edges
  122. isBorderEdge.resize(EV.rows());
  123. for(unsigned i=0; i<EV.rows(); ++i)
  124. isBorderEdge[i] = (EF(i,0) == -1) || ((EF(i,1) == -1));
  125. // Generate normals per face
  126. igl::per_face_normals(V, F, N);
  127. // Generate reference frames
  128. for(unsigned fid=0; fid<F.rows(); ++fid)
  129. {
  130. // First edge
  131. Vector3d e1 = V.row(F(fid,1)) - V.row(F(fid,0));
  132. e1.normalize();
  133. Vector3d e2 = N.row(fid);
  134. e2 = e2.cross(e1);
  135. e2.normalize();
  136. MatrixXd TP(2,3);
  137. TP << e1.transpose(), e2.transpose();
  138. TPs.push_back(TP);
  139. }
  140. // Alloc internal variables
  141. angles = VectorXd::Zero(F.rows());
  142. p = VectorXi::Zero(EV.rows());
  143. pFixed.resize(EV.rows());
  144. k = VectorXd::Zero(EV.rows());
  145. singularityIndex = VectorXd::Zero(V.rows());
  146. // Reset the constraints
  147. resetConstraints();
  148. // Compute k, differences between reference frames
  149. computek();
  150. softAlpha = 0.5;
  151. }
  152. void igl::NRosyField::setSoftAlpha(double alpha)
  153. {
  154. assert(alpha >= 0 && alpha < 1);
  155. softAlpha = alpha;
  156. }
  157. void igl::NRosyField::prepareSystemMatrix(const int N)
  158. {
  159. using namespace std;
  160. using namespace Eigen;
  161. double Nd = N;
  162. // Minimize the MIQ energy
  163. // Energy on edge ij is
  164. // (t_i - t_j + kij + pij*(2*pi/N))^2
  165. // Partial derivatives:
  166. // t_i: 2 ( t_i - t_j + kij + pij*(2*pi/N)) = 0
  167. // t_j: 2 (-t_i + t_j - kij - pij*(2*pi/N)) = 0
  168. // pij: 4pi/N ( t_i - t_j + kij + pij*(2*pi/N)) = 0
  169. //
  170. // t_i t_j pij kij
  171. // t_i [ 2 -2 4pi/N 2 ]
  172. // t_j [ -2 2 -4pi/N -2 ]
  173. // pij [ 4pi/N -4pi/N 2*(2pi/N)^2 4pi/N ]
  174. // Count and tag the variables
  175. tag_t = VectorXi::Constant(F.rows(),-1);
  176. vector<int> id_t;
  177. int count = 0;
  178. for(unsigned i=0; i<F.rows(); ++i)
  179. if (!isHard[i])
  180. {
  181. tag_t(i) = count++;
  182. id_t.push_back(i);
  183. }
  184. unsigned count_t = id_t.size();
  185. tag_p = VectorXi::Constant(EF.rows(),-1);
  186. vector<int> id_p;
  187. for(unsigned i=0; i<EF.rows(); ++i)
  188. {
  189. if (!pFixed[i])
  190. {
  191. // if it is not fixed then it is a variable
  192. tag_p(i) = count++;
  193. }
  194. // if it is not a border edge,
  195. if (!isBorderEdge[i])
  196. {
  197. // and it is not between two fixed faces
  198. if (!(isHard[EF(i,0)] && isHard[EF(i,1)]))
  199. {
  200. // then it participates in the energy!
  201. id_p.push_back(i);
  202. }
  203. }
  204. }
  205. unsigned count_p = count - count_t;
  206. // System sizes: A (count_t + count_p) x (count_t + count_p)
  207. // b (count_t + count_p)
  208. b = VectorXd::Zero(count_t + count_p);
  209. std::vector<Eigen::Triplet<double> > T;
  210. T.reserve(3 * 4 * count_p);
  211. for(unsigned r=0; r<id_p.size(); ++r)
  212. {
  213. int eid = id_p[r];
  214. int i = EF(eid,0);
  215. int j = EF(eid,1);
  216. bool isFixed_i = isHard[i];
  217. bool isFixed_j = isHard[j];
  218. bool isFixed_p = pFixed[eid];
  219. int row;
  220. // (i)-th row: t_i [ 2 -2 4pi/N 2 ]
  221. if (!isFixed_i)
  222. {
  223. row = tag_t[i];
  224. if (isFixed_i) b(row) += -2 * hard[i]; else T.push_back(Eigen::Triplet<double>(row,tag_t[i] , 2 ));
  225. if (isFixed_j) b(row) += 2 * hard[j]; else T.push_back(Eigen::Triplet<double>(row,tag_t[j] ,-2 ));
  226. if (isFixed_p) b(row) += -((4 * M_PI)/Nd) * p[eid] ; else T.push_back(Eigen::Triplet<double>(row,tag_p[eid],((4 * M_PI)/Nd)));
  227. b(row) += -2 * k[eid];
  228. assert(hard[i] == hard[i]);
  229. assert(hard[j] == hard[j]);
  230. assert(p[eid] == p[eid]);
  231. assert(k[eid] == k[eid]);
  232. assert(b(row) == b(row));
  233. }
  234. // (j)+1 -th row: t_j [ -2 2 -4pi/N -2 ]
  235. if (!isFixed_j)
  236. {
  237. row = tag_t[j];
  238. if (isFixed_i) b(row) += 2 * hard[i]; else T.push_back(Eigen::Triplet<double>(row,tag_t[i] , -2 ));
  239. if (isFixed_j) b(row) += -2 * hard[j]; else T.push_back(Eigen::Triplet<double>(row,tag_t[j] , 2 ));
  240. if (isFixed_p) b(row) += ((4 * M_PI)/Nd) * p[eid] ; else T.push_back(Eigen::Triplet<double>(row,tag_p[eid],-((4 * M_PI)/Nd)));
  241. b(row) += 2 * k[eid];
  242. assert(k[eid] == k[eid]);
  243. assert(b(row) == b(row));
  244. }
  245. // (r*3)+2 -th row: pij [ 4pi/N -4pi/N 2*(2pi/N)^2 4pi/N ]
  246. if (!isFixed_p)
  247. {
  248. row = tag_p[eid];
  249. if (isFixed_i) b(row) += -(4 * M_PI)/Nd * hard[i]; else T.push_back(Eigen::Triplet<double>(row,tag_t[i] , (4 * M_PI)/Nd ));
  250. if (isFixed_j) b(row) += (4 * M_PI)/Nd * hard[j]; else T.push_back(Eigen::Triplet<double>(row,tag_t[j] , -(4 * M_PI)/Nd ));
  251. if (isFixed_p) b(row) += -(2 * pow(((2*M_PI)/Nd),2)) * p[eid] ; else T.push_back(Eigen::Triplet<double>(row,tag_p[eid], (2 * pow(((2*M_PI)/Nd),2))));
  252. b(row) += - (4 * M_PI)/Nd * k[eid];
  253. assert(k[eid] == k[eid]);
  254. assert(b(row) == b(row));
  255. }
  256. }
  257. A = SparseMatrix<double>(count_t + count_p, count_t + count_p);
  258. A.setFromTriplets(T.begin(), T.end());
  259. // Soft constraints
  260. bool addSoft = false;
  261. for(unsigned i=0; i<wSoft.size();++i)
  262. if (wSoft[i] != 0)
  263. addSoft = true;
  264. if (addSoft)
  265. {
  266. cerr << " Adding soft here: " << endl;
  267. cerr << " softAplha: " << softAlpha << endl;
  268. VectorXd bSoft = VectorXd::Zero(count_t + count_p);
  269. std::vector<Eigen::Triplet<double> > TSoft;
  270. TSoft.reserve(2 * count_p);
  271. for(unsigned i=0; i<F.rows(); ++i)
  272. {
  273. int varid = tag_t[i];
  274. if (varid != -1) // if it is a variable in the system
  275. {
  276. TSoft.push_back(Eigen::Triplet<double>(varid,varid,wSoft[i]));
  277. bSoft[varid] += wSoft[i] * soft[i];
  278. }
  279. }
  280. SparseMatrix<double> ASoft(count_t + count_p, count_t + count_p);
  281. ASoft.setFromTriplets(TSoft.begin(), TSoft.end());
  282. // ofstream s("/Users/daniele/As.txt");
  283. // for(unsigned i=0; i<TSoft.size(); ++i)
  284. // s << TSoft[i].row() << " " << TSoft[i].col() << " " << TSoft[i].value() << endl;
  285. // s.close();
  286. // ofstream s2("/Users/daniele/bs.txt");
  287. // for(unsigned i=0; i<bSoft.rows(); ++i)
  288. // s2 << bSoft(i) << endl;
  289. // s2.close();
  290. // Stupid Eigen bug
  291. SparseMatrix<double> Atmp (count_t + count_p, count_t + count_p);
  292. SparseMatrix<double> Atmp2(count_t + count_p, count_t + count_p);
  293. SparseMatrix<double> Atmp3(count_t + count_p, count_t + count_p);
  294. // Merge the two part of the energy
  295. Atmp = (1.0 - softAlpha)*A;
  296. Atmp2 = softAlpha * ASoft;
  297. Atmp3 = Atmp+Atmp2;
  298. A = Atmp3;
  299. b = b*(1.0 - softAlpha) + bSoft * softAlpha;
  300. }
  301. // ofstream s("/Users/daniele/A.txt");
  302. // for (int k=0; k<A.outerSize(); ++k)
  303. // for (SparseMatrix<double>::InnerIterator it(A,k); it; ++it)
  304. // {
  305. // s << it.row() << " " << it.col() << " " << it.value() << endl;
  306. // }
  307. // s.close();
  308. //
  309. // ofstream s2("/Users/daniele/b.txt");
  310. // for(unsigned i=0; i<b.rows(); ++i)
  311. // s2 << b(i) << endl;
  312. // s2.close();
  313. }
  314. void igl::NRosyField::solveNoRoundings()
  315. {
  316. using namespace std;
  317. using namespace Eigen;
  318. // Solve the linear system
  319. SimplicialLDLT<SparseMatrix<double> > solver;
  320. solver.compute(A);
  321. VectorXd x = solver.solve(b);
  322. // Copy the result back
  323. for(unsigned i=0; i<F.rows(); ++i)
  324. if (tag_t[i] != -1)
  325. angles[i] = x(tag_t[i]);
  326. else
  327. angles[i] = hard[i];
  328. for(unsigned i=0; i<EF.rows(); ++i)
  329. if(tag_p[i] != -1)
  330. p[i] = roundl(x[tag_p[i]]);
  331. }
  332. void igl::NRosyField::solveRoundings()
  333. {
  334. using namespace std;
  335. using namespace Eigen;
  336. unsigned n = A.rows();
  337. gmm::col_matrix< gmm::wsvector< double > > gmm_A;
  338. std::vector<double> gmm_b;
  339. std::vector<int> ids_to_round;
  340. std::vector<double> x;
  341. gmm_A.resize(n,n);
  342. gmm_b.resize(n);
  343. x.resize(n);
  344. // Copy A
  345. for (int k=0; k<A.outerSize(); ++k)
  346. for (SparseMatrix<double>::InnerIterator it(A,k); it; ++it)
  347. {
  348. gmm_A(it.row(),it.col()) += it.value();
  349. }
  350. // Copy b
  351. for(unsigned i=0; i<n;++i)
  352. gmm_b[i] = b[i];
  353. // Set variables to round
  354. ids_to_round.clear();
  355. for(unsigned i=0; i<tag_p.size();++i)
  356. if(tag_p[i] != -1)
  357. ids_to_round.push_back(tag_p[i]);
  358. // Empty constraints
  359. gmm::row_matrix< gmm::wsvector< double > > gmm_C(0, n);
  360. COMISO::ConstrainedSolver cs;
  361. //print_miso_settings(cs.misolver());
  362. cs.solve(gmm_C, gmm_A, x, gmm_b, ids_to_round, 0.0, false, true);
  363. // Copy the result back
  364. for(unsigned i=0; i<F.rows(); ++i)
  365. if (tag_t[i] != -1)
  366. angles[i] = x[tag_t[i]];
  367. else
  368. angles[i] = hard[i];
  369. for(unsigned i=0; i<EF.rows(); ++i)
  370. if(tag_p[i] != -1)
  371. p[i] = roundl(x[tag_p[i]]);
  372. }
  373. void igl::NRosyField::roundAndFix()
  374. {
  375. for(unsigned i=0; i<p.rows(); ++i)
  376. pFixed[i] = true;
  377. }
  378. void igl::NRosyField::roundAndFixToZero()
  379. {
  380. for(unsigned i=0; i<p.rows(); ++i)
  381. {
  382. pFixed[i] = true;
  383. p[i] = 0;
  384. }
  385. }
  386. void igl::NRosyField::solve(const int N)
  387. {
  388. // Reduce the search space by fixing matchings
  389. reduceSpace();
  390. // Build the system
  391. prepareSystemMatrix(N);
  392. // Solve with integer roundings
  393. solveRoundings();
  394. // This is a very greedy solving strategy
  395. // // Solve with no roundings
  396. // solveNoRoundings();
  397. //
  398. // // Round all p and fix them
  399. // roundAndFix();
  400. //
  401. // // Build the system
  402. // prepareSystemMatrix(N);
  403. //
  404. // // Solve with no roundings (they are all fixed)
  405. // solveNoRoundings();
  406. // Find the cones
  407. findCones(N);
  408. }
  409. void igl::NRosyField::setConstraintHard(const int fid, const Eigen::Vector3d& v)
  410. {
  411. isHard[fid] = true;
  412. hard(fid) = convert3DtoLocal(fid, v);
  413. }
  414. void igl::NRosyField::setConstraintSoft(const int fid, const double w, const Eigen::Vector3d& v)
  415. {
  416. wSoft(fid) = w;
  417. soft(fid) = convert3DtoLocal(fid, v);
  418. }
  419. void igl::NRosyField::resetConstraints()
  420. {
  421. using namespace std;
  422. using namespace Eigen;
  423. isHard.resize(F.rows());
  424. for(unsigned i=0; i<F.rows(); ++i)
  425. isHard[i] = false;
  426. hard = VectorXd::Zero(F.rows());
  427. wSoft = VectorXd::Zero(F.rows());
  428. soft = VectorXd::Zero(F.rows());
  429. }
  430. Eigen::MatrixXd igl::NRosyField::getFieldPerFace()
  431. {
  432. using namespace std;
  433. using namespace Eigen;
  434. MatrixXd result(F.rows(),3);
  435. for(unsigned i=0; i<F.rows(); ++i)
  436. result.row(i) = convertLocalto3D(i, angles(i));
  437. return result;
  438. }
  439. Eigen::MatrixXd igl::NRosyField::getFFieldPerFace()
  440. {
  441. using namespace std;
  442. using namespace Eigen;
  443. MatrixXd result(F.rows(),6);
  444. for(unsigned i=0; i<F.rows(); ++i)
  445. {
  446. Vector3d v1 = convertLocalto3D(i, angles(i));
  447. Vector3d n = N.row(i);
  448. Vector3d v2 = n.cross(v1);
  449. v1.normalize();
  450. v2.normalize();
  451. result.block(i,0,1,3) = v1.transpose();
  452. result.block(i,3,1,3) = v2.transpose();
  453. }
  454. return result;
  455. }
  456. void igl::NRosyField::computek()
  457. {
  458. using namespace std;
  459. using namespace Eigen;
  460. // For every non-border edge
  461. for (unsigned eid=0; eid<EF.rows(); ++eid)
  462. {
  463. if (!isBorderEdge[eid])
  464. {
  465. int fid0 = EF(eid,0);
  466. int fid1 = EF(eid,1);
  467. Vector3d N0 = N.row(fid0);
  468. Vector3d N1 = N.row(fid1);
  469. // find common edge on triangle 0 and 1
  470. int fid0_vc = -1;
  471. int fid1_vc = -1;
  472. for (unsigned i=0;i<3;++i)
  473. {
  474. if (EV(eid,0) == F(fid0,i))
  475. fid0_vc = i;
  476. if (EV(eid,1) == F(fid1,i))
  477. fid1_vc = i;
  478. }
  479. assert(fid0_vc != -1);
  480. assert(fid1_vc != -1);
  481. Vector3d common_edge = V.row(F(fid0,(fid0_vc+1)%3)) - V.row(F(fid0,fid0_vc));
  482. common_edge.normalize();
  483. // Map the two triangles in a new space where the common edge is the x axis and the N0 the z axis
  484. MatrixXd P(3,3);
  485. VectorXd o = V.row(F(fid0,fid0_vc));
  486. VectorXd tmp = -N0.cross(common_edge);
  487. P << common_edge, tmp, N0;
  488. P.transposeInPlace();
  489. MatrixXd V0(3,3);
  490. V0.row(0) = V.row(F(fid0,0)).transpose() -o;
  491. V0.row(1) = V.row(F(fid0,1)).transpose() -o;
  492. V0.row(2) = V.row(F(fid0,2)).transpose() -o;
  493. V0 = (P*V0.transpose()).transpose();
  494. assert(V0(0,2) < 10e-10);
  495. assert(V0(1,2) < 10e-10);
  496. assert(V0(2,2) < 10e-10);
  497. MatrixXd V1(3,3);
  498. V1.row(0) = V.row(F(fid1,0)).transpose() -o;
  499. V1.row(1) = V.row(F(fid1,1)).transpose() -o;
  500. V1.row(2) = V.row(F(fid1,2)).transpose() -o;
  501. V1 = (P*V1.transpose()).transpose();
  502. assert(V1(fid1_vc,2) < 10e-10);
  503. assert(V1((fid1_vc+1)%3,2) < 10e-10);
  504. // compute rotation R such that R * N1 = N0
  505. // i.e. map both triangles to the same plane
  506. double alpha = -atan2(V1((fid1_vc+2)%3,2),V1((fid1_vc+2)%3,1));
  507. MatrixXd R(3,3);
  508. R << 1, 0, 0,
  509. 0, cos(alpha), -sin(alpha) ,
  510. 0, sin(alpha), cos(alpha);
  511. V1 = (R*V1.transpose()).transpose();
  512. assert(V1(0,2) < 10e-10);
  513. assert(V1(1,2) < 10e-10);
  514. assert(V1(2,2) < 10e-10);
  515. // measure the angle between the reference frames
  516. // k_ij is the angle between the triangle on the left and the one on the right
  517. VectorXd ref0 = V0.row(1) - V0.row(0);
  518. VectorXd ref1 = V1.row(1) - V1.row(0);
  519. ref0.normalize();
  520. ref1.normalize();
  521. double ktemp = atan2(ref1(1),ref1(0)) - atan2(ref0(1),ref0(0));
  522. // just to be sure, rotate ref0 using angle ktemp...
  523. MatrixXd R2(2,2);
  524. R2 << cos(ktemp), -sin(ktemp), sin(ktemp), cos(ktemp);
  525. tmp = R2*ref0.head<2>();
  526. assert(tmp(0) - ref1(0) < 10^10);
  527. assert(tmp(1) - ref1(1) < 10^10);
  528. k[eid] = ktemp;
  529. }
  530. }
  531. }
  532. void igl::NRosyField::reduceSpace()
  533. {
  534. using namespace std;
  535. using namespace Eigen;
  536. // All variables are free in the beginning
  537. for(unsigned i=0; i<EV.rows(); ++i)
  538. pFixed[i] = false;
  539. vector<VectorXd> debug;
  540. // debug
  541. // MatrixXd B(F.rows(),3);
  542. // for(unsigned i=0; i<F.rows(); ++i)
  543. // B.row(i) = 1./3. * (V.row(F(i,0)) + V.row(F(i,1)) + V.row(F(i,2)));
  544. vector<bool> visited(EV.rows());
  545. for(unsigned i=0; i<EV.rows(); ++i)
  546. visited[i] = false;
  547. vector<bool> starting(EV.rows());
  548. for(unsigned i=0; i<EV.rows(); ++i)
  549. starting[i] = false;
  550. queue<int> q;
  551. for(unsigned i=0; i<F.rows(); ++i)
  552. if (isHard[i] || wSoft[i] != 0)
  553. {
  554. q.push(i);
  555. starting[i] = true;
  556. }
  557. // Reduce the search space (see MI paper)
  558. while (!q.empty())
  559. {
  560. int c = q.front();
  561. q.pop();
  562. visited[c] = true;
  563. for(int i=0; i<3; ++i)
  564. {
  565. int eid = FE(c,i);
  566. int fid = TT(c,i);
  567. // skip borders
  568. if (fid != -1)
  569. {
  570. assert((EF(eid,0) == c && EF(eid,1) == fid) || (EF(eid,1) == c && EF(eid,0) == fid));
  571. // for every neighbouring face
  572. if (!visited[fid] && !starting[fid])
  573. {
  574. pFixed[eid] = true;
  575. p[eid] = 0;
  576. visited[fid] = true;
  577. q.push(fid);
  578. }
  579. }
  580. else
  581. {
  582. // fix borders
  583. pFixed[eid] = true;
  584. p[eid] = 0;
  585. }
  586. }
  587. }
  588. // Force matchings between fixed faces
  589. for(unsigned i=0; i<F.rows();++i)
  590. {
  591. if (isHard[i])
  592. {
  593. for(unsigned int j=0; j<3; ++j)
  594. {
  595. int fid = TT(i,j);
  596. if ((fid!=-1) && (isHard[fid]))
  597. {
  598. // i and fid are adjacent and fixed
  599. int eid = FE(i,j);
  600. int fid0 = EF(eid,0);
  601. int fid1 = EF(eid,1);
  602. pFixed[eid] = true;
  603. p[eid] = roundl(2.0/M_PI*(hard(fid1) - hard(fid0) - k(eid)));
  604. }
  605. }
  606. }
  607. }
  608. // std::ofstream s("/Users/daniele/debug.txt");
  609. // for(unsigned i=0; i<debug.size(); i += 2)
  610. // s << debug[i].transpose() << " " << debug[i+1].transpose() << endl;
  611. // s.close();
  612. }
  613. double igl::NRosyField::convert3DtoLocal(unsigned fid, const Eigen::Vector3d& v)
  614. {
  615. using namespace std;
  616. using namespace Eigen;
  617. // Project onto the tangent plane
  618. Vector2d vp = TPs[fid] * v;
  619. // Convert to angle
  620. return atan2(vp(1),vp(0));
  621. }
  622. Eigen::Vector3d igl::NRosyField::convertLocalto3D(unsigned fid, double a)
  623. {
  624. using namespace std;
  625. using namespace Eigen;
  626. Vector2d vp(cos(a),sin(a));
  627. return vp.transpose() * TPs[fid];
  628. }
  629. Eigen::VectorXd igl::NRosyField::angleDefect()
  630. {
  631. Eigen::VectorXd A = Eigen::VectorXd::Constant(V.rows(),-2*M_PI);
  632. for (unsigned i=0; i < F.rows(); ++i)
  633. {
  634. for (int j = 0; j < 3; ++j)
  635. {
  636. Eigen::VectorXd a = V.row(F(i,(j+1)%3)) - V.row(F(i,j));
  637. Eigen::VectorXd b = V.row(F(i,(j+2)%3)) - V.row(F(i,j));
  638. double t = a.transpose()*b;
  639. t /= (a.norm() * b.norm());
  640. A(F(i,j)) += acos(t);
  641. }
  642. }
  643. return A;
  644. }
  645. void igl::NRosyField::findCones(int N)
  646. {
  647. // Compute I0, see http://www.graphics.rwth-aachen.de/media/papers/bommes_zimmer_2009_siggraph_011.pdf for details
  648. Eigen::VectorXd I0 = Eigen::VectorXd::Zero(V.rows());
  649. // first the k
  650. for (unsigned i=0; i < EV.rows(); ++i)
  651. {
  652. if (!isBorderEdge[i])
  653. {
  654. I0(EV(i,0)) -= k(i);
  655. I0(EV(i,1)) += k(i);
  656. }
  657. }
  658. // then the A
  659. Eigen::VectorXd A = angleDefect();
  660. I0 = I0 + A;
  661. // normalize
  662. I0 = I0 / (2*M_PI);
  663. // round to integer (remove numerical noise)
  664. for (unsigned i=0; i < I0.size(); ++i)
  665. I0(i) = round(I0(i));
  666. // compute I
  667. Eigen::VectorXd I = I0;
  668. for (unsigned i=0; i < EV.rows(); ++i)
  669. {
  670. if (!isBorderEdge[i])
  671. {
  672. I(EV(i,0)) -= double(p(i))/double(N);
  673. I(EV(i,1)) += double(p(i))/double(N);
  674. }
  675. }
  676. // Clear the vertices on the edges
  677. for (unsigned i=0; i < EV.rows(); ++i)
  678. {
  679. if (isBorderEdge[i])
  680. {
  681. I0(EV(i,0)) = 0;
  682. I0(EV(i,1)) = 0;
  683. I(EV(i,0)) = 0;
  684. I(EV(i,1)) = 0;
  685. A(EV(i,0)) = 0;
  686. A(EV(i,1)) = 0;
  687. }
  688. }
  689. singularityIndex = I;
  690. }
  691. Eigen::VectorXd igl::NRosyField::getSingularityIndexPerVertex()
  692. {
  693. return singularityIndex;
  694. }
  695. IGL_INLINE void igl::nrosy(
  696. const Eigen::MatrixXd& V,
  697. const Eigen::MatrixXi& F,
  698. const Eigen::VectorXi& b,
  699. const Eigen::MatrixXd& bc,
  700. const Eigen::VectorXi& b_soft,
  701. const Eigen::VectorXd& w_soft,
  702. const Eigen::MatrixXd& bc_soft,
  703. const int N,
  704. const double soft,
  705. Eigen::MatrixXd& R,
  706. Eigen::VectorXd& S
  707. )
  708. {
  709. // Init solver
  710. igl::NRosyField solver(V,F);
  711. // Add hard constraints
  712. for (unsigned i=0; i<b.size();++i)
  713. solver.setConstraintHard(b(i),bc.row(i));
  714. // Add soft constraints
  715. for (unsigned i=0; i<b_soft.size();++i)
  716. solver.setConstraintSoft(b_soft(i),w_soft(i),bc_soft.row(i));
  717. // Set the soft constraints global weight
  718. solver.setSoftAlpha(soft);
  719. // Interpolate
  720. solver.solve(N);
  721. // Copy the result back
  722. R = solver.getFieldPerFace();
  723. // Extract singularity indices
  724. S = solver.getSingularityIndexPerVertex();
  725. }