nrosy.cpp 24 KB

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