miq.cpp 52 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210121112121213121412151216121712181219122012211222122312241225122612271228122912301231123212331234123512361237123812391240124112421243124412451246124712481249125012511252125312541255125612571258125912601261126212631264126512661267126812691270127112721273127412751276127712781279128012811282128312841285128612871288128912901291129212931294129512961297129812991300130113021303130413051306130713081309131013111312131313141315131613171318131913201321132213231324132513261327132813291330133113321333133413351336133713381339134013411342134313441345134613471348134913501351135213531354135513561357135813591360136113621363136413651366136713681369137013711372137313741375137613771378137913801381138213831384138513861387138813891390139113921393139413951396139713981399140014011402140314041405140614071408140914101411141214131414141514161417141814191420142114221423142414251426142714281429143014311432143314341435143614371438143914401441144214431444144514461447144814491450145114521453145414551456145714581459146014611462146314641465146614671468146914701471147214731474147514761477147814791480148114821483148414851486148714881489149014911492149314941495149614971498149915001501150215031504150515061507150815091510151115121513151415151516151715181519152015211522152315241525152615271528152915301531153215331534153515361537153815391540154115421543
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2014 Daniele Panozzo <daniele.panozzo@gmail.com>, Olga Diamanti <olga.diam@gmail.com>, Kevin Walliman <wkevin@student.ethz.ch>
  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 "miq.h"
  9. #include "../../local_basis.h"
  10. #include "../../triangle_triangle_adjacency.h"
  11. #include "../../cut_mesh.h"
  12. #include "../../LinSpaced.h"
  13. // includes for VertexIndexing
  14. #include "../../HalfEdgeIterator.h"
  15. #include "../../is_border_vertex.h"
  16. #include "../../vertex_triangle_adjacency.h"
  17. // includes for PoissonSolver
  18. #include "../../slice_into.h"
  19. #include "../../grad.h"
  20. #include "../../cotmatrix.h"
  21. #include "../../doublearea.h"
  22. #include <gmm/gmm.h>
  23. #include <CoMISo/Solver/ConstrainedSolver.hh>
  24. #include <CoMISo/Solver/MISolver.hh>
  25. #include <CoMISo/Solver/GMM_Tools.hh>
  26. //
  27. #include "igl/cross_field_mismatch.h"
  28. #include "../../comb_frame_field.h"
  29. #include "../../comb_cross_field.h"
  30. #include "../../cut_mesh_from_singularities.h"
  31. #include "../../find_cross_field_singularities.h"
  32. #include "../../compute_frame_field_bisectors.h"
  33. #include "../../rotate_vectors.h"
  34. #ifndef NDEBUG
  35. #include <fstream>
  36. #endif
  37. #include <iostream>
  38. #include "../../matlab_format.h"
  39. #define DEBUGPRINT 0
  40. namespace igl {
  41. namespace copyleft {
  42. namespace comiso {
  43. struct SeamInfo
  44. {
  45. int v0,v0p;
  46. unsigned int integerVar;
  47. int mismatch;
  48. IGL_INLINE SeamInfo(int _v0,
  49. int _v0p,
  50. int _mismatch,
  51. unsigned int _integerVar);
  52. IGL_INLINE SeamInfo(const SeamInfo &S1);
  53. };
  54. struct MeshSystemInfo
  55. {
  56. MeshSystemInfo()
  57. {
  58. num_vert_variables = 0;
  59. num_integer_cuts = 0;
  60. }
  61. ////number of vertices variables
  62. unsigned int num_vert_variables;
  63. ///num of integer for cuts
  64. unsigned int num_integer_cuts;
  65. ///this are used for drawing purposes
  66. std::vector<SeamInfo> edgeSeamInfo;
  67. };
  68. template <typename DerivedV, typename DerivedF>
  69. class VertexIndexing
  70. {
  71. public:
  72. // Input:
  73. const Eigen::PlainObjectBase<DerivedV> &V;
  74. const Eigen::PlainObjectBase<DerivedF> &F;
  75. const Eigen::PlainObjectBase<DerivedV> &Vcut;
  76. const Eigen::PlainObjectBase<DerivedF> &Fcut;
  77. const Eigen::PlainObjectBase<DerivedF> &TT;
  78. const Eigen::PlainObjectBase<DerivedF> &TTi;
  79. const Eigen::Matrix<int, Eigen::Dynamic, 3> &mismatch;
  80. const Eigen::Matrix<int, Eigen::Dynamic, 1> &singular; // bool
  81. const Eigen::Matrix<int, Eigen::Dynamic, 3> &seams; // 3 bool
  82. ///this handle for mesh TODO: move with the other global variables
  83. MeshSystemInfo systemInfo;
  84. IGL_INLINE VertexIndexing(const Eigen::PlainObjectBase<DerivedV> &_V,
  85. const Eigen::PlainObjectBase<DerivedF> &_F,
  86. const Eigen::PlainObjectBase<DerivedV> &_Vcut,
  87. const Eigen::PlainObjectBase<DerivedF> &_Fcut,
  88. const Eigen::PlainObjectBase<DerivedF> &_TT,
  89. const Eigen::PlainObjectBase<DerivedF> &_TTi,
  90. const Eigen::Matrix<int, Eigen::Dynamic, 3> &_mismatch,
  91. const Eigen::Matrix<int, Eigen::Dynamic, 1> &_singular,
  92. const Eigen::Matrix<int, Eigen::Dynamic, 3> &_seams
  93. );
  94. // provide information about every vertex per seam
  95. IGL_INLINE void initSeamInfo();
  96. private:
  97. struct VertexInfo{
  98. int v; // vertex index (according to V)
  99. int f0, k0; // face and local edge information of the edge that connects this vertex to the previous vertex (previous in the vector)
  100. int f1, k1; // face and local edge information of the other face corresponding to the same edge
  101. VertexInfo(int _v, int _f0, int _k0, int _f1, int _k1) :
  102. v(_v), f0(_f0), k0(_k0), f1(_f1), k1(_k1){}
  103. bool operator==(VertexInfo const& other){
  104. return other.v == v;
  105. }
  106. };
  107. IGL_INLINE void getSeamInfo(int f0,
  108. int f1,
  109. int indexE,
  110. int &v0, int &v1,
  111. int &v0p, int &v1p,
  112. int &_mismatch);
  113. IGL_INLINE std::vector<std::vector<VertexInfo> > getVerticesPerSeam();
  114. };
  115. template <typename DerivedV, typename DerivedF>
  116. class PoissonSolver
  117. {
  118. private:
  119. // Penalization term for integer variables used in mixedIntegerSolve
  120. const double PENALIZATION = 0.000001;
  121. public:
  122. IGL_INLINE void solvePoisson(Eigen::VectorXd stiffness,
  123. double gradientSize = 0.1,
  124. double gridResolution = 1.,
  125. bool directRound = true,
  126. unsigned int localIter = 0,
  127. bool doRound = true,
  128. bool singularityRound = true,
  129. const std::vector<int> &roundVertices = std::vector<int>(),
  130. const std::vector<std::vector<int>> &hardFeatures = std::vector<std::vector<int> >());
  131. IGL_INLINE PoissonSolver(const Eigen::PlainObjectBase<DerivedV> &_V,
  132. const Eigen::PlainObjectBase<DerivedF> &_F,
  133. const Eigen::PlainObjectBase<DerivedV> &_Vcut,
  134. const Eigen::PlainObjectBase<DerivedF> &_Fcut,
  135. const Eigen::PlainObjectBase<DerivedF> &_TT,
  136. const Eigen::PlainObjectBase<DerivedF> &_TTi,
  137. const Eigen::PlainObjectBase<DerivedV> &_PD1,
  138. const Eigen::PlainObjectBase<DerivedV> &_PD2,
  139. const Eigen::Matrix<int, Eigen::Dynamic, 1>&_singular,
  140. const MeshSystemInfo &_systemInfo
  141. );
  142. const Eigen::PlainObjectBase<DerivedV> &V;
  143. const Eigen::PlainObjectBase<DerivedF> &F;
  144. const Eigen::PlainObjectBase<DerivedV> &Vcut;
  145. const Eigen::PlainObjectBase<DerivedF> &Fcut;
  146. const Eigen::PlainObjectBase<DerivedF> &TT;
  147. const Eigen::PlainObjectBase<DerivedF> &TTi;
  148. const Eigen::PlainObjectBase<DerivedV> &PD1;
  149. const Eigen::PlainObjectBase<DerivedV> &PD2;
  150. const Eigen::Matrix<int, Eigen::Dynamic, 1> &singular; // bool
  151. const MeshSystemInfo &systemInfo;
  152. // Internal:
  153. Eigen::VectorXd Handle_Stiffness;
  154. std::vector<std::vector<int> > VF;
  155. std::vector<std::vector<int> > VFi;
  156. Eigen::MatrixXd UV; // this is probably useless
  157. // Output:
  158. // per wedge UV coordinates, 6 coordinates (1 face) per row
  159. Eigen::MatrixXd WUV;
  160. // per vertex UV coordinates, Vcut.rows() x 2
  161. Eigen::MatrixXd UV_out;
  162. // Matrices
  163. Eigen::SparseMatrix<double> Lhs;
  164. Eigen::SparseMatrix<double> Constraints;
  165. Eigen::VectorXd rhs;
  166. Eigen::VectorXd constraints_rhs;
  167. ///vector of unknowns
  168. std::vector< double > X;
  169. ////REAL PART
  170. ///number of fixed vertex
  171. unsigned int n_fixed_vars;
  172. ///the number of REAL variables for vertices
  173. unsigned int n_vert_vars;
  174. ///total number of variables of the system,
  175. ///do not consider constraints, but consider integer vars
  176. unsigned int num_total_vars;
  177. //////INTEGER PART
  178. ///the total number of integer variables
  179. unsigned int n_integer_vars;
  180. ///CONSTRAINT PART
  181. ///number of cuts constraints
  182. unsigned int num_cut_constraint;
  183. // number of user-defined constraints
  184. unsigned int num_userdefined_constraint;
  185. ///total number of constraints equations
  186. unsigned int num_constraint_equations;
  187. ///vector of blocked vertices
  188. std::vector<int> Hard_constraints;
  189. ///vector of indexes to round
  190. std::vector<int> ids_to_round;
  191. ///vector of indexes to round
  192. std::vector<std::vector<int > > userdefined_constraints;
  193. ///boolean that is true if rounding to integer is needed
  194. bool integer_rounding;
  195. ///START COMMON MATH FUNCTIONS
  196. ///return the complex encoding the rotation
  197. ///for a given mismatch interval
  198. IGL_INLINE std::complex<double> getRotationComplex(int interval);
  199. ///END COMMON MATH FUNCTIONS
  200. ///START FIXING VERTICES
  201. ///set a given vertex as fixed
  202. IGL_INLINE void addFixedVertex(int v);
  203. ///find vertex to fix in case we're using
  204. ///a vector field NB: multiple components not handled
  205. IGL_INLINE void findFixedVertField();
  206. ///find hard constraint depending if using or not
  207. ///a vector field
  208. IGL_INLINE void findFixedVert();
  209. IGL_INLINE int getFirstVertexIndex(int v);
  210. ///fix the vertices which are flagged as fixed
  211. IGL_INLINE void fixBlockedVertex();
  212. ///END FIXING VERTICES
  213. ///HANDLING SINGULARITY
  214. //set the singularity round to integer location
  215. IGL_INLINE void addSingularityRound();
  216. IGL_INLINE void addToRoundVertices(std::vector<int> ids);
  217. ///START GENERIC SYSTEM FUNCTIONS
  218. //build the Laplacian matrix cycling over all range maps
  219. //and over all faces
  220. IGL_INLINE void buildLaplacianMatrix(double vfscale = 1);
  221. ///find different sized of the system
  222. IGL_INLINE void findSizes();
  223. IGL_INLINE void allocateSystem();
  224. ///intitialize the whole matrix
  225. IGL_INLINE void initMatrix();
  226. ///map back coordinates after that
  227. ///the system has been solved
  228. IGL_INLINE void mapCoords();
  229. ///END GENERIC SYSTEM FUNCTIONS
  230. ///set the constraints for the inter-range cuts
  231. IGL_INLINE void buildSeamConstraintsExplicitTranslation();
  232. ///set the constraints for the inter-range cuts
  233. IGL_INLINE void buildUserDefinedConstraints();
  234. ///call of the mixed integer solver
  235. IGL_INLINE void mixedIntegerSolve(double coneGridRes = 1,
  236. bool directRound = true,
  237. unsigned int localIter = 0);
  238. IGL_INLINE void clearUserConstraint();
  239. IGL_INLINE void addSharpEdgeConstraint(int fid, int vid);
  240. };
  241. template <typename DerivedV, typename DerivedF, typename DerivedU>
  242. class MIQ_class
  243. {
  244. private:
  245. const Eigen::PlainObjectBase<DerivedV> &V;
  246. const Eigen::PlainObjectBase<DerivedF> &F;
  247. DerivedV Vcut;
  248. DerivedF Fcut;
  249. Eigen::MatrixXd UV_out;
  250. DerivedF FUV_out;
  251. // internal
  252. DerivedF TT;
  253. DerivedF TTi;
  254. // Stiffness per face
  255. Eigen::VectorXd stiffnessVector;
  256. DerivedV B1, B2, B3;
  257. public:
  258. IGL_INLINE MIQ_class(const Eigen::PlainObjectBase<DerivedV> &V_,
  259. const Eigen::PlainObjectBase<DerivedF> &F_,
  260. const Eigen::PlainObjectBase<DerivedV> &PD1_combed,
  261. const Eigen::PlainObjectBase<DerivedV> &PD2_combed,
  262. const Eigen::Matrix<int, Eigen::Dynamic, 3> &mismatch,
  263. const Eigen::Matrix<int, Eigen::Dynamic, 1> &singular,
  264. const Eigen::Matrix<int, Eigen::Dynamic, 3> &seams,
  265. Eigen::PlainObjectBase<DerivedU> &UV,
  266. Eigen::PlainObjectBase<DerivedF> &FUV,
  267. double gradientSize = 30.0,
  268. double stiffness = 5.0,
  269. bool directRound = false,
  270. unsigned int iter = 5,
  271. unsigned int localIter = 5,
  272. bool doRound = true,
  273. bool singularityRound = true,
  274. std::vector<int> roundVertices = std::vector<int>(),
  275. std::vector<std::vector<int> > hardFeatures = std::vector<std::vector<int> >());
  276. IGL_INLINE void extractUV(Eigen::PlainObjectBase<DerivedU> &UV_out,
  277. Eigen::PlainObjectBase<DerivedF> &FUV_out);
  278. private:
  279. IGL_INLINE int NumFlips(const Eigen::MatrixXd& WUV);
  280. IGL_INLINE double Distortion(int f, double h, const Eigen::MatrixXd& WUV);
  281. IGL_INLINE double LaplaceDistortion(int f, double h, const Eigen::MatrixXd& WUV);
  282. IGL_INLINE bool updateStiffeningJacobianDistorsion(double grad_size, const Eigen::MatrixXd& WUV);
  283. IGL_INLINE bool IsFlipped(const Eigen::Vector2d &uv0,
  284. const Eigen::Vector2d &uv1,
  285. const Eigen::Vector2d &uv2);
  286. IGL_INLINE bool IsFlipped(int i, const Eigen::MatrixXd& WUV);
  287. };
  288. };
  289. };
  290. }
  291. IGL_INLINE igl::copyleft::comiso::SeamInfo::SeamInfo(int _v0,
  292. int _v0p,
  293. int _mismatch,
  294. unsigned int _integerVar)
  295. {
  296. v0=_v0;
  297. v0p=_v0p;
  298. integerVar=_integerVar;
  299. mismatch=_mismatch;
  300. }
  301. IGL_INLINE igl::copyleft::comiso::SeamInfo::SeamInfo(const SeamInfo &S1)
  302. {
  303. v0=S1.v0;
  304. v0p=S1.v0p;
  305. integerVar=S1.integerVar;
  306. mismatch=S1.mismatch;
  307. }
  308. template <typename DerivedV, typename DerivedF>
  309. IGL_INLINE igl::copyleft::comiso::VertexIndexing<DerivedV, DerivedF>::VertexIndexing(const Eigen::PlainObjectBase<DerivedV> &_V,
  310. const Eigen::PlainObjectBase<DerivedF> &_F,
  311. const Eigen::PlainObjectBase<DerivedV> &_Vcut,
  312. const Eigen::PlainObjectBase<DerivedF> &_Fcut,
  313. const Eigen::PlainObjectBase<DerivedF> &_TT,
  314. const Eigen::PlainObjectBase<DerivedF> &_TTi,
  315. const Eigen::Matrix<int, Eigen::Dynamic, 3> &_mismatch,
  316. const Eigen::Matrix<int, Eigen::Dynamic, 1> &_singular,
  317. const Eigen::Matrix<int, Eigen::Dynamic, 3> &_seams
  318. ):
  319. V(_V),
  320. F(_F),
  321. Vcut(_Vcut),
  322. Fcut(_Fcut),
  323. TT(_TT),
  324. TTi(_TTi),
  325. mismatch(_mismatch),
  326. singular(_singular),
  327. seams(_seams)
  328. {
  329. #ifdef DEBUG_PRINT
  330. cerr<<igl::matlab_format(Handle_Seams,"Handle_Seams");
  331. #endif
  332. systemInfo.num_vert_variables=Vcut.rows();
  333. systemInfo.num_integer_cuts=0;
  334. }
  335. template <typename DerivedV, typename DerivedF>
  336. IGL_INLINE void igl::copyleft::comiso::VertexIndexing<DerivedV, DerivedF>::getSeamInfo(const int f0,
  337. const int f1,
  338. const int indexE,
  339. int &v0, int &v1,
  340. int &v0p, int &v1p,
  341. int &_mismatch)
  342. {
  343. int edgef0 = indexE;
  344. v0 = Fcut(f0,edgef0);
  345. v1 = Fcut(f0,(edgef0+1)%3);
  346. ////get the index on opposite side
  347. assert(TT(f0,edgef0) == f1);
  348. int edgef1 = TTi(f0,edgef0);
  349. v1p = Fcut(f1,edgef1);
  350. v0p = Fcut(f1,(edgef1+1)%3);
  351. _mismatch = mismatch(f0,edgef0);
  352. assert(F(f0,edgef0) == F(f1,((edgef1+1)%3)));
  353. assert(F(f0,((edgef0+1)%3)) == F(f1,edgef1));
  354. }
  355. template <typename DerivedV, typename DerivedF>
  356. IGL_INLINE std::vector<std::vector<typename igl::copyleft::comiso::VertexIndexing<DerivedV, DerivedF>::VertexInfo> > igl::copyleft::comiso::VertexIndexing<DerivedV, DerivedF>::getVerticesPerSeam()
  357. {
  358. // Return value
  359. std::vector<std::vector<VertexInfo> >verticesPerSeam;
  360. // for every vertex, keep track of their adjacent vertices on seams.
  361. // regular vertices have two neighbors on a seam, start- and endvertices may have any other numbers of neighbors (e.g. 1 or 3)
  362. std::vector<std::list<VertexInfo> > VVSeam(V.rows());
  363. Eigen::MatrixXi F_hit = Eigen::MatrixXi::Zero(F.rows(), 3);
  364. for (unsigned int f=0; f<F.rows();f++)
  365. {
  366. int f0 = f;
  367. for(int k0=0; k0<3; k0++){
  368. int f1 = TT(f0,k0);
  369. if(f1 == -1)
  370. continue;
  371. if (seams(f0, k0) != 0 && F_hit(f0,k0) == 0)
  372. {
  373. int v0 = F(f0, k0);
  374. int v1 = F(f0, (k0+1)%3);
  375. int k1 = TTi(f0,k0);
  376. VVSeam[v0].push_back(VertexInfo(v1, f0, k0, f1, k1));
  377. VVSeam[v1].push_back(VertexInfo(v0, f0, k0, f1, k1));
  378. F_hit(f0, k0) = 1;
  379. F_hit(f1, k1) = 1;
  380. }
  381. }
  382. }
  383. // Find start vertices, i.e. vertices that start or end a seam branch
  384. std::vector<int> startVertexIndices;
  385. std::vector<bool> isStartVertex(V.rows());
  386. for (unsigned int i=0;i<V.rows();i++)
  387. {
  388. isStartVertex[i] = false;
  389. // vertices with two neighbors are regular vertices, unless the vertex is a singularity, in which case it qualifies as a start vertex
  390. if ((!VVSeam[i].empty() && VVSeam[i].size() != 2) || singular(i) != 0)
  391. {
  392. startVertexIndices.push_back(i);
  393. isStartVertex[i] = true;
  394. }
  395. }
  396. // For each startVertex, walk along its seam
  397. for (auto element : startVertexIndices)
  398. {
  399. auto startVertexNeighbors = &VVSeam[element];
  400. size_t neighborSize = startVertexNeighbors->size();
  401. // explore every seam to which this vertex is a start vertex
  402. // note: a vertex can never be a start vertex and a regular vertex simultaneously
  403. for (size_t j = 0; j < neighborSize; j++)
  404. {
  405. std::vector<VertexInfo> thisSeam; // temporary container
  406. // Create vertexInfo struct for start vertex
  407. VertexInfo startVertex = VertexInfo(element, -1, -1, -1, -1);// -1 values are arbitrary (will never be used)
  408. VertexInfo currentVertex = startVertex;
  409. // Add start vertex to the seam
  410. thisSeam.push_back(currentVertex);
  411. // advance on the seam
  412. auto currentVertexNeighbors = startVertexNeighbors;
  413. auto nextVertex = currentVertexNeighbors->front();
  414. currentVertexNeighbors->pop_front();
  415. // bogus initialization due to lack of def. constructor
  416. VertexInfo prevVertex = startVertex;
  417. while (true)
  418. {
  419. // move to the next vertex
  420. prevVertex = currentVertex;
  421. currentVertex = nextVertex;
  422. currentVertexNeighbors = &VVSeam[nextVertex.v];
  423. // add current vertex to this seam
  424. thisSeam.push_back(currentVertex);
  425. // remove the previous vertex
  426. auto it = std::find(currentVertexNeighbors->begin(), currentVertexNeighbors->end(), prevVertex);
  427. assert(it != currentVertexNeighbors->end());
  428. currentVertexNeighbors->erase(it);
  429. if (currentVertexNeighbors->size() == 1 && !isStartVertex[currentVertex.v])
  430. {
  431. nextVertex = currentVertexNeighbors->front();
  432. currentVertexNeighbors->pop_front();
  433. }
  434. else
  435. break;
  436. }
  437. verticesPerSeam.push_back(thisSeam);
  438. }
  439. }
  440. return verticesPerSeam;
  441. }
  442. template <typename DerivedV, typename DerivedF>
  443. IGL_INLINE void igl::copyleft::comiso::VertexIndexing<DerivedV, DerivedF>::initSeamInfo()
  444. {
  445. auto verticesPerSeam = getVerticesPerSeam();
  446. systemInfo.edgeSeamInfo.clear();
  447. unsigned int integerVar = 0;
  448. // Loop over each seam
  449. for(auto seam : verticesPerSeam){
  450. //choose initial side of the seam such that the start vertex corresponds to Fcut(f, k) and the end vertex corresponds to Fcut(f, (k+1)%3) and not vice versa.
  451. int priorVertexIdx;
  452. if(seam.size() > 2){
  453. auto v1 = seam[1];
  454. auto v2 = seam[2];
  455. if(Fcut(v1.f0, (v1.k0+1) % 3) == Fcut(v2.f0, v2.k0) || Fcut(v1.f0, (v1.k0+1) % 3) == Fcut(v2.f1, v2.k1)){
  456. priorVertexIdx = Fcut(v1.f0, v1.k0);
  457. }
  458. else{
  459. priorVertexIdx = Fcut(v1.f1, v1.k1);
  460. assert(Fcut(v1.f1, (v1.k1+1) % 3) == Fcut(v2.f0, v2.k0) || Fcut(v1.f1, (v1.k1+1) % 3) == Fcut(v2.f1, v2.k1));
  461. }
  462. }
  463. else{
  464. auto v1 = seam[1];
  465. priorVertexIdx = Fcut(v1.f0, v1.k0);
  466. }
  467. // Loop over each vertex of the seam
  468. for(auto it=seam.begin()+1; it != seam.end(); ++it){
  469. auto vertex = *it;
  470. // choose the correct side of the seam
  471. int f,k,ff;
  472. if(priorVertexIdx == Fcut(vertex.f0, vertex.k0)){
  473. f = vertex.f0; ff = vertex.f1;
  474. k = vertex.k0;
  475. }
  476. else{
  477. f = vertex.f1; ff = vertex.f0;
  478. k = vertex.k1;
  479. assert(priorVertexIdx == Fcut(vertex.f1, vertex.k1));
  480. }
  481. int vtx0,vtx0p,vtx1,vtx1p;
  482. int MM;
  483. getSeamInfo(f, ff, k, vtx0, vtx1, vtx0p, vtx1p, MM);
  484. systemInfo.edgeSeamInfo.push_back(SeamInfo(vtx0,vtx0p,MM,integerVar));
  485. if(it == seam.end() -1){
  486. systemInfo.edgeSeamInfo.push_back(SeamInfo(vtx1,vtx1p,MM,integerVar));
  487. }
  488. priorVertexIdx = vtx1;
  489. }
  490. // use the same integer for each seam
  491. integerVar++;
  492. }
  493. systemInfo.num_integer_cuts = integerVar;
  494. #ifndef NDEBUG
  495. int totalNVerticesOnSeams = 0;
  496. for(auto const & seam : verticesPerSeam){
  497. totalNVerticesOnSeams += seam.size();
  498. }
  499. assert(systemInfo.edgeSeamInfo.size() == totalNVerticesOnSeams);
  500. #endif
  501. }
  502. template <typename DerivedV, typename DerivedF>
  503. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::solvePoisson(Eigen::VectorXd stiffness,
  504. double gradientSize,
  505. double gridResolution,
  506. bool directRound,
  507. unsigned int localIter,
  508. bool doRound,
  509. bool singularityRound,
  510. const std::vector<int> &roundVertices,
  511. const std::vector<std::vector<int>> &hardFeatures)
  512. {
  513. Handle_Stiffness = stiffness;
  514. //initialization of flags and data structures
  515. integer_rounding=doRound;
  516. ids_to_round.clear();
  517. clearUserConstraint();
  518. // copy the user constraints number
  519. for (const auto & element : hardFeatures)
  520. {
  521. addSharpEdgeConstraint(element[0], element[1]);
  522. }
  523. ///Initializing Matrix
  524. clock_t t0 = clock();
  525. ///initialize the matrix ALLOCATING SPACE
  526. initMatrix();
  527. if (DEBUGPRINT)
  528. printf("\n ALLOCATED THE MATRIX \n");
  529. ///build the Laplacian system
  530. buildLaplacianMatrix(gradientSize);
  531. // add seam constraints
  532. buildSeamConstraintsExplicitTranslation();
  533. // add user defined constraints
  534. buildUserDefinedConstraints();
  535. ////add the Lagrange multiplier
  536. fixBlockedVertex();
  537. if (DEBUGPRINT)
  538. printf("\n BUILT THE MATRIX \n");
  539. if (integer_rounding)
  540. addToRoundVertices(roundVertices);
  541. if (singularityRound)
  542. addSingularityRound();
  543. clock_t t1 = clock();
  544. if (DEBUGPRINT) printf("\n time:%ld \n",t1-t0);
  545. if (DEBUGPRINT) printf("\n SOLVING \n");
  546. mixedIntegerSolve(gridResolution, directRound, localIter);
  547. clock_t t2 = clock();
  548. if (DEBUGPRINT) printf("\n time:%ld \n",t2-t1);
  549. if (DEBUGPRINT) printf("\n ASSIGNING COORDS \n");
  550. mapCoords();
  551. clock_t t3 = clock();
  552. if (DEBUGPRINT) printf("\n time:%ld \n",t3-t2);
  553. if (DEBUGPRINT) printf("\n FINISHED \n");
  554. }
  555. template <typename DerivedV, typename DerivedF>
  556. IGL_INLINE igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>
  557. ::PoissonSolver(const Eigen::PlainObjectBase<DerivedV> &_V,
  558. const Eigen::PlainObjectBase<DerivedF> &_F,
  559. const Eigen::PlainObjectBase<DerivedV> &_Vcut,
  560. const Eigen::PlainObjectBase<DerivedF> &_Fcut,
  561. const Eigen::PlainObjectBase<DerivedF> &_TT,
  562. const Eigen::PlainObjectBase<DerivedF> &_TTi,
  563. const Eigen::PlainObjectBase<DerivedV> &_PD1,
  564. const Eigen::PlainObjectBase<DerivedV> &_PD2,
  565. const Eigen::Matrix<int, Eigen::Dynamic, 1>&_singular,
  566. const MeshSystemInfo &_systemInfo
  567. ):
  568. V(_V),
  569. F(_F),
  570. Vcut(_Vcut),
  571. Fcut(_Fcut),
  572. TT(_TT),
  573. TTi(_TTi),
  574. PD1(_PD1),
  575. PD2(_PD2),
  576. singular(_singular),
  577. systemInfo(_systemInfo)
  578. {
  579. n_fixed_vars = 0;
  580. n_vert_vars = 0;
  581. num_total_vars = 0;
  582. n_integer_vars = 0;
  583. num_cut_constraint = 0;
  584. num_userdefined_constraint = 0;
  585. num_constraint_equations = 0;
  586. integer_rounding = false;
  587. UV = Eigen::MatrixXd(V.rows(),2);
  588. WUV = Eigen::MatrixXd(F.rows(),6);
  589. UV_out = Eigen::MatrixXd(Vcut.rows(),2);
  590. igl::vertex_triangle_adjacency(V,F,VF,VFi);
  591. }
  592. ///START COMMON MATH FUNCTIONS
  593. ///return the complex encoding the rotation
  594. ///for a given mismatch interval
  595. template <typename DerivedV, typename DerivedF>
  596. IGL_INLINE std::complex<double> igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::getRotationComplex(
  597. int interval)
  598. {
  599. assert((interval>=0)&&(interval<4));
  600. switch(interval)
  601. {
  602. case 0:return {1,0};
  603. case 1:return {0,1};
  604. case 2: return {-1,0};
  605. default:return {0,-1};
  606. }
  607. }
  608. ///END COMMON MATH FUNCTIONS
  609. ///START FIXING VERTICES
  610. ///set a given vertex as fixed
  611. template <typename DerivedV, typename DerivedF>
  612. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::addFixedVertex(int v)
  613. {
  614. n_fixed_vars++;
  615. Hard_constraints.push_back(v);
  616. }
  617. ///find vertex to fix in case we're using
  618. ///a vector field NB: multiple components not handled
  619. template <typename DerivedV, typename DerivedF>
  620. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::findFixedVertField()
  621. {
  622. Hard_constraints.clear();
  623. n_fixed_vars=0;
  624. //fix the first singularity
  625. for (unsigned int v=0;v<V.rows();v++)
  626. {
  627. if (singular(v))
  628. {
  629. addFixedVertex(v);
  630. UV.row(v) << 0,0;
  631. return;
  632. }
  633. }
  634. ///if anything fixed fix the first
  635. addFixedVertex(0);
  636. UV.row(0) << 0,0;
  637. std::cerr << "No vertices to fix, I am fixing the first vertex to the origin!" << std::endl;
  638. }
  639. ///find hard constraint depending if using or not
  640. ///a vector field
  641. template <typename DerivedV, typename DerivedF>
  642. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::findFixedVert()
  643. {
  644. Hard_constraints.clear();
  645. findFixedVertField();
  646. }
  647. template <typename DerivedV, typename DerivedF>
  648. IGL_INLINE int igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::getFirstVertexIndex(int v)
  649. {
  650. return Fcut(VF[v][0],VFi[v][0]);
  651. }
  652. ///fix the vertices which are flagged as fixed
  653. template <typename DerivedV, typename DerivedF>
  654. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::fixBlockedVertex()
  655. {
  656. int offset_row = num_cut_constraint*2;
  657. unsigned int constr_num = 0;
  658. for (unsigned int i=0;i<Hard_constraints.size();i++)
  659. {
  660. int v = Hard_constraints[i];
  661. ///get first index of the vertex that must blocked
  662. //int index=v->vertex_index[0];
  663. int index = getFirstVertexIndex(v);
  664. ///multiply times 2 because of uv
  665. int indexvert = index*2;
  666. ///find the first free row to add the constraint
  667. int indexRow = offset_row + constr_num * 2;
  668. int indexCol = indexRow;
  669. ///add fixing constraint LHS
  670. Constraints.coeffRef(indexRow, indexvert) += 1;
  671. Constraints.coeffRef(indexRow+1,indexvert+1) += 1;
  672. ///add fixing constraint RHS
  673. constraints_rhs[indexCol] = UV(v,0);
  674. constraints_rhs[indexCol+1] = UV(v,1);
  675. constr_num++;
  676. }
  677. assert(constr_num==n_fixed_vars);
  678. }
  679. ///END FIXING VERTICES
  680. ///HANDLING SINGULARITY
  681. //set the singularity round to integer location
  682. template <typename DerivedV, typename DerivedF>
  683. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::addSingularityRound()
  684. {
  685. for (unsigned int v=0;v<V.rows();v++)
  686. {
  687. if (singular(v))
  688. {
  689. int index0= getFirstVertexIndex(v);
  690. ids_to_round.push_back( index0*2 );
  691. ids_to_round.push_back((index0*2)+1);
  692. }
  693. }
  694. }
  695. template <typename DerivedV, typename DerivedF>
  696. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::addToRoundVertices(std::vector<int> ids)
  697. {
  698. for(auto index : ids)
  699. {
  700. if (index < 0 || index >= V.rows())
  701. std::cerr << "WARNING: Ignored round vertex constraint, vertex " << index << " does not exist in the mesh." << std::endl;
  702. int index0 = getFirstVertexIndex(index);
  703. ids_to_round.push_back( index0 * 2 );
  704. ids_to_round.push_back((index0 * 2)+1);
  705. }
  706. }
  707. ///START GENERIC SYSTEM FUNCTIONS
  708. template <typename DerivedV, typename DerivedF>
  709. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::buildLaplacianMatrix(double vfscale)
  710. {
  711. Eigen::VectorXi idx = igl::LinSpaced<Eigen::VectorXi >(Vcut.rows(), 0, 2*Vcut.rows()-2);
  712. Eigen::VectorXi idx2 = igl::LinSpaced<Eigen::VectorXi >(Vcut.rows(), 1, 2*Vcut.rows()-1);
  713. // get gradient matrix
  714. Eigen::SparseMatrix<double> G(Fcut.rows() * 3, Vcut.rows());
  715. igl::grad(Vcut, Fcut, G);
  716. // get triangle weights
  717. Eigen::VectorXd dblA(Fcut.rows());
  718. igl::doublearea(Vcut, Fcut, dblA);
  719. // compute intermediate result
  720. Eigen::SparseMatrix<double> G2;
  721. G2 = G.transpose() * dblA.replicate<3,1>().asDiagonal() * Handle_Stiffness.replicate<3,1>().asDiagonal();
  722. /// Compute LHS
  723. Eigen::SparseMatrix<double> Cotmatrix;
  724. Cotmatrix = 0.5 * G2 * G;
  725. igl::slice_into(Cotmatrix, idx, idx, Lhs);
  726. igl::slice_into(Cotmatrix, idx2, idx2, Lhs);
  727. /// Compute RHS
  728. // reshape nrosy vectors
  729. const Eigen::MatrixXd u = Eigen::Map<const Eigen::MatrixXd>(PD1.data(),Fcut.rows()*3,1); // this mimics a reshape at the cost of a copy.
  730. const Eigen::MatrixXd v = Eigen::Map<const Eigen::MatrixXd>(PD2.data(),Fcut.rows()*3,1); // this mimics a reshape at the cost of a copy.
  731. // multiply with weights
  732. Eigen::VectorXd rhs1 = G2 * u * 0.5 * vfscale;
  733. Eigen::VectorXd rhs2 = -G2 * v * 0.5 * vfscale;
  734. igl::slice_into(rhs1, idx, 1, rhs);
  735. igl::slice_into(rhs2, idx2, 1, rhs);
  736. }
  737. ///find different sized of the system
  738. template <typename DerivedV, typename DerivedF>
  739. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::findSizes()
  740. {
  741. ///find the vertex that need to be fixed
  742. findFixedVert();
  743. ///REAL PART
  744. n_vert_vars = systemInfo.num_vert_variables;
  745. ///INTEGER PART
  746. ///the total number of integer variables
  747. n_integer_vars = systemInfo.num_integer_cuts;
  748. ///CONSTRAINT PART
  749. num_cut_constraint = systemInfo.edgeSeamInfo.size();
  750. num_constraint_equations = num_cut_constraint * 2 + n_fixed_vars * 2 + num_userdefined_constraint;
  751. ///total variable of the system
  752. num_total_vars = (n_vert_vars+n_integer_vars) * 2;
  753. ///initialize matrix size
  754. if (DEBUGPRINT) printf("\n*** SYSTEM VARIABLES *** \n");
  755. if (DEBUGPRINT) printf("* NUM REAL VERTEX VARIABLES %ud \n",n_vert_vars);
  756. if (DEBUGPRINT) printf("\n*** INTEGER VARIABLES *** \n");
  757. if (DEBUGPRINT) printf("* NUM INTEGER VARIABLES %ud \n",n_integer_vars);
  758. if (DEBUGPRINT) printf("\n*** CONSTRAINTS *** \n ");
  759. if (DEBUGPRINT) printf("* NUM FIXED CONSTRAINTS %ud\n",n_fixed_vars);
  760. if (DEBUGPRINT) printf("* NUM CUTS CONSTRAINTS %ud\n",num_cut_constraint);
  761. if (DEBUGPRINT) printf("* NUM USER DEFINED CONSTRAINTS %ud\n",num_userdefined_constraint);
  762. if (DEBUGPRINT) printf("\n*** TOTAL SIZE *** \n");
  763. if (DEBUGPRINT) printf("* TOTAL VARIABLE SIZE (WITH INTEGER TRASL) %ud \n",num_total_vars);
  764. if (DEBUGPRINT) printf("* TOTAL CONSTRAINTS %ud \n",num_constraint_equations);
  765. }
  766. template <typename DerivedV, typename DerivedF>
  767. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::allocateSystem()
  768. {
  769. Lhs.resize(n_vert_vars * 2, n_vert_vars * 2);
  770. Constraints.resize(num_constraint_equations, num_total_vars);
  771. rhs.resize(n_vert_vars * 2);
  772. constraints_rhs.resize(num_constraint_equations);
  773. printf("\n INITIALIZED SPARSE MATRIX OF %ud x %ud \n",n_vert_vars*2, n_vert_vars*2);
  774. printf("\n INITIALIZED SPARSE MATRIX OF %ud x %ud \n",num_constraint_equations, num_total_vars);
  775. printf("\n INITIALIZED VECTOR OF %ud x 1 \n",n_vert_vars*2);
  776. printf("\n INITIALIZED VECTOR OF %ud x 1 \n",num_constraint_equations);
  777. }
  778. ///intitialize the whole matrix
  779. template <typename DerivedV, typename DerivedF>
  780. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::initMatrix()
  781. {
  782. findSizes();
  783. allocateSystem();
  784. }
  785. ///map back coordinates after that
  786. ///the system has been solved
  787. template <typename DerivedV, typename DerivedF>
  788. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::mapCoords()
  789. {
  790. ///map coords to faces
  791. for (unsigned int f=0;f<Fcut.rows();f++)
  792. {
  793. for (int k=0;k<3;k++)
  794. {
  795. //get the index of the variable in the system
  796. int indexUV = Fcut(f,k);
  797. ///then get U and V coords
  798. double U=X[indexUV*2];
  799. double V=X[indexUV*2+1];
  800. WUV(f,k*2 + 0) = U;
  801. WUV(f,k*2 + 1) = V;
  802. }
  803. }
  804. for(unsigned int i = 0; i < Vcut.rows(); i++){
  805. UV_out(i,0) = X[i*2];
  806. UV_out(i,1) = X[i*2+1];
  807. }
  808. }
  809. ///END GENERIC SYSTEM FUNCTIONS
  810. ///set the constraints for the inter-range cuts
  811. template <typename DerivedV, typename DerivedF>
  812. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::buildSeamConstraintsExplicitTranslation()
  813. {
  814. ///current constraint row
  815. int constr_row = 0;
  816. for (unsigned int i=0; i<num_cut_constraint; i++)
  817. {
  818. int interval = systemInfo.edgeSeamInfo[i].mismatch;
  819. if (interval==1)
  820. interval=3;
  821. else
  822. if(interval==3)
  823. interval=1;
  824. int p0 = systemInfo.edgeSeamInfo[i].v0;
  825. int p0p = systemInfo.edgeSeamInfo[i].v0p;
  826. std::complex<double> rot = getRotationComplex(interval);
  827. ///get the integer variable
  828. unsigned int integerVar = n_vert_vars + systemInfo.edgeSeamInfo[i].integerVar;
  829. if (integer_rounding)
  830. {
  831. ids_to_round.push_back(integerVar*2);
  832. ids_to_round.push_back(integerVar*2+1);
  833. }
  834. // cross boundary compatibility conditions
  835. Constraints.coeffRef(constr_row, 2*p0) += rot.real();
  836. Constraints.coeffRef(constr_row, 2*p0+1) += -rot.imag();
  837. Constraints.coeffRef(constr_row+1, 2*p0) += rot.imag();
  838. Constraints.coeffRef(constr_row+1, 2*p0+1) += rot.real();
  839. Constraints.coeffRef(constr_row, 2*p0p) += -1;
  840. Constraints.coeffRef(constr_row+1, 2*p0p+1) += -1;
  841. Constraints.coeffRef(constr_row, 2*integerVar) += 1;
  842. Constraints.coeffRef(constr_row+1, 2*integerVar+1) += 1;
  843. constraints_rhs[constr_row] = 0;
  844. constraints_rhs[constr_row+1] = 0;
  845. constr_row += 2;
  846. }
  847. }
  848. ///set the constraints for the inter-range cuts
  849. template <typename DerivedV, typename DerivedF>
  850. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::buildUserDefinedConstraints()
  851. {
  852. /// the user defined constraints are at the end
  853. unsigned int constr_row = num_cut_constraint*2 + n_fixed_vars*2;
  854. assert(num_userdefined_constraint == userdefined_constraints.size());
  855. for (unsigned int i = 0; i < num_userdefined_constraint; i++)
  856. {
  857. for (unsigned int j = 0; j < userdefined_constraints[i].size()-1; ++j)
  858. {
  859. Constraints.coeffRef(constr_row, j) = userdefined_constraints[i][j];
  860. }
  861. constraints_rhs[constr_row] = userdefined_constraints[i][userdefined_constraints[i].size()-1];
  862. constr_row +=1;
  863. }
  864. }
  865. ///call of the mixed integer solver
  866. template <typename DerivedV, typename DerivedF>
  867. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::mixedIntegerSolve(double coneGridRes,
  868. bool directRound,
  869. unsigned int localIter)
  870. {
  871. X = std::vector<double>((n_vert_vars+n_integer_vars)*2);
  872. if (DEBUGPRINT)
  873. printf("\n ALLOCATED X \n");
  874. ///variables part
  875. const int sizeMatrix = (n_vert_vars + n_integer_vars) * 2;
  876. const int scalarSize = n_vert_vars * 2;
  877. ///matrix A
  878. gmm::col_matrix< gmm::wsvector< double > > A(sizeMatrix,sizeMatrix); // lhs matrix variables
  879. ///constraints part
  880. int CsizeX = num_constraint_equations;
  881. int CsizeY = sizeMatrix+1;
  882. gmm::row_matrix< gmm::wsvector< double > > C(CsizeX,CsizeY); // constraints
  883. if (DEBUGPRINT)
  884. printf("\n ALLOCATED QMM STRUCTURES \n");
  885. std::vector<double> B(sizeMatrix,0); // rhs
  886. if (DEBUGPRINT)
  887. printf("\n ALLOCATED RHS STRUCTURES \n");
  888. //// copy LHS
  889. for (int k=0; k < Lhs.outerSize(); ++k){
  890. for (Eigen::SparseMatrix<double>::InnerIterator it(Lhs,k); it; ++it){
  891. int row = it.row();
  892. int col = it.col();
  893. A(row, col) += it.value();
  894. }
  895. }
  896. //// copy Constraints
  897. for (int k=0; k < Constraints.outerSize(); ++k){
  898. for (Eigen::SparseMatrix<double>::InnerIterator it(Constraints,k); it; ++it){
  899. int row = it.row();
  900. int col = it.col();
  901. C(row, col) += it.value();
  902. }
  903. }
  904. if (DEBUGPRINT)
  905. printf("\n SET %d INTEGER VALUES \n",n_integer_vars);
  906. int offline_index = scalarSize;
  907. for(unsigned int i = 0; i < n_integer_vars*2; ++i)
  908. {
  909. int index=offline_index+i;
  910. A(index, index) = PENALIZATION;
  911. }
  912. if (DEBUGPRINT)
  913. printf("\n SET RHS \n");
  914. // copy RHS
  915. for(unsigned int i = 0; i < scalarSize; ++i)
  916. {
  917. B[i] = rhs[i] * coneGridRes;
  918. }
  919. // copy constraint RHS
  920. if (DEBUGPRINT)
  921. printf("\n SET %d CONSTRAINTS \n",num_constraint_equations);
  922. for(unsigned int i = 0; i < num_constraint_equations; ++i)
  923. {
  924. C(i, sizeMatrix) = -constraints_rhs[i] * coneGridRes;
  925. }
  926. COMISO::ConstrainedSolver solver;
  927. solver.misolver().set_local_iters(localIter);
  928. solver.misolver().set_direct_rounding(directRound);
  929. std::sort(ids_to_round.begin(),ids_to_round.end());
  930. auto new_end=std::unique(ids_to_round.begin(),ids_to_round.end());
  931. long int dist = distance(ids_to_round.begin(),new_end);
  932. ids_to_round.resize(dist);
  933. solver.solve( C, A, X, B, ids_to_round, 0.0, false, false);
  934. }
  935. template <typename DerivedV, typename DerivedF>
  936. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::clearUserConstraint()
  937. {
  938. num_userdefined_constraint = 0;
  939. userdefined_constraints.clear();
  940. }
  941. template <typename DerivedV, typename DerivedF>
  942. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::addSharpEdgeConstraint(int fid, int vid)
  943. {
  944. // prepare constraint
  945. std::vector<int> c(systemInfo.num_vert_variables*2 + 1, 0);
  946. int v1 = Fcut(fid,vid);
  947. int v2 = Fcut(fid,(vid+1)%3);
  948. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> e = Vcut.row(v2) - Vcut.row(v1);
  949. e = e.normalized();
  950. double d1 = fabs(e.dot(PD1.row(fid).normalized()));
  951. double d2 = fabs(e.dot(PD2.row(fid).normalized()));
  952. int offset = 0;
  953. if (d1>d2)
  954. offset = 1;
  955. ids_to_round.push_back((v1 * 2) + offset);
  956. ids_to_round.push_back((v2 * 2) + offset);
  957. // add constraint
  958. c[(v1 * 2) + offset] = 1;
  959. c[(v2 * 2) + offset] = -1;
  960. // add to the user-defined constraints
  961. num_userdefined_constraint++;
  962. userdefined_constraints.push_back(c);
  963. }
  964. template <typename DerivedV, typename DerivedF, typename DerivedU>
  965. IGL_INLINE igl::copyleft::comiso::MIQ_class<DerivedV, DerivedF, DerivedU>::MIQ_class(
  966. const Eigen::PlainObjectBase<DerivedV> &V_,
  967. const Eigen::PlainObjectBase<DerivedF> &F_,
  968. const Eigen::PlainObjectBase<DerivedV> &PD1_combed,
  969. const Eigen::PlainObjectBase<DerivedV> &PD2_combed,
  970. const Eigen::Matrix<int, Eigen::Dynamic, 3> &mismatch,
  971. const Eigen::Matrix<int, Eigen::Dynamic, 1> &singular,
  972. const Eigen::Matrix<int, Eigen::Dynamic, 3> &seams,
  973. Eigen::PlainObjectBase<DerivedU> &UV,
  974. Eigen::PlainObjectBase<DerivedF> &FUV,
  975. double gradientSize,
  976. double stiffness,
  977. bool directRound,
  978. unsigned int iter,
  979. unsigned int localIter,
  980. bool doRound,
  981. bool singularityRound,
  982. std::vector<int> roundVertices,
  983. std::vector<std::vector<int> > hardFeatures):
  984. V(V_),
  985. F(F_)
  986. {
  987. igl::cut_mesh(V, F, seams, Vcut, Fcut);
  988. igl::local_basis(V,F,B1,B2,B3);
  989. igl::triangle_triangle_adjacency(F,TT,TTi);
  990. // Prepare indexing for the linear system
  991. VertexIndexing<DerivedV, DerivedF> VInd(V, F, Vcut, Fcut, TT, TTi, mismatch, singular, seams);
  992. VInd.initSeamInfo();
  993. // Assemble the system and solve
  994. PoissonSolver<DerivedV, DerivedF> PSolver(V,
  995. F,
  996. Vcut,
  997. Fcut,
  998. TT,
  999. TTi,
  1000. PD1_combed,
  1001. PD2_combed,
  1002. singular,
  1003. VInd.systemInfo);
  1004. stiffnessVector = Eigen::VectorXd::Constant(F.rows(),1);
  1005. if (iter > 0) // do stiffening
  1006. {
  1007. for (unsigned int i=0;i<iter;i++)
  1008. {
  1009. PSolver.solvePoisson(stiffnessVector, gradientSize, 1.f, directRound, localIter, doRound, singularityRound,
  1010. roundVertices, hardFeatures);
  1011. int nflips=NumFlips(PSolver.WUV);
  1012. bool folded = updateStiffeningJacobianDistorsion(gradientSize,PSolver.WUV);
  1013. printf("ITERATION %d FLIPS %d \n",i,nflips);
  1014. if (!folded)break;
  1015. }
  1016. }
  1017. else
  1018. {
  1019. PSolver.solvePoisson(stiffnessVector, gradientSize, 1.f, directRound, localIter, doRound, singularityRound,
  1020. roundVertices, hardFeatures);
  1021. }
  1022. int nflips=NumFlips(PSolver.WUV);
  1023. printf("**** END OPTIMIZING #FLIPS %d ****\n",nflips);
  1024. UV_out = PSolver.UV_out;
  1025. FUV_out = PSolver.Fcut;
  1026. fflush(stdout);
  1027. }
  1028. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1029. IGL_INLINE void igl::copyleft::comiso::MIQ_class<DerivedV, DerivedF, DerivedU>::extractUV(Eigen::PlainObjectBase<DerivedU> &UV_out,
  1030. Eigen::PlainObjectBase<DerivedF> &FUV_out)
  1031. {
  1032. UV_out = this->UV_out;
  1033. FUV_out = this->FUV_out;
  1034. }
  1035. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1036. IGL_INLINE int igl::copyleft::comiso::MIQ_class<DerivedV, DerivedF, DerivedU>::NumFlips(const Eigen::MatrixXd& WUV)
  1037. {
  1038. int numFl=0;
  1039. for (unsigned int i=0;i<F.rows();i++)
  1040. {
  1041. if (IsFlipped(i, WUV))
  1042. numFl++;
  1043. }
  1044. return numFl;
  1045. }
  1046. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1047. IGL_INLINE double igl::copyleft::comiso::MIQ_class<DerivedV, DerivedF, DerivedU>::Distortion(int f, double h, const Eigen::MatrixXd& WUV)
  1048. {
  1049. assert(h > 0);
  1050. Eigen::Vector2d uv0,uv1,uv2;
  1051. uv0 << WUV(f,0), WUV(f,1);
  1052. uv1 << WUV(f,2), WUV(f,3);
  1053. uv2 << WUV(f,4), WUV(f,5);
  1054. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> p0 = Vcut.row(Fcut(f,0));
  1055. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> p1 = Vcut.row(Fcut(f,1));
  1056. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> p2 = Vcut.row(Fcut(f,2));
  1057. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> norm = (p1 - p0).cross(p2 - p0);
  1058. double area2 = norm.norm();
  1059. double area2_inv = 1.0 / area2;
  1060. norm *= area2_inv;
  1061. if (area2 > 0)
  1062. {
  1063. // Singular values of the Jacobian
  1064. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> neg_t0 = norm.cross(p2 - p1);
  1065. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> neg_t1 = norm.cross(p0 - p2);
  1066. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> neg_t2 = norm.cross(p1 - p0);
  1067. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> diffu = (neg_t0 * uv0(0) +neg_t1 *uv1(0) + neg_t2 * uv2(0) )*area2_inv;
  1068. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> diffv = (neg_t0 * uv0(1) + neg_t1*uv1(1) + neg_t2*uv2(1) )*area2_inv;
  1069. // first fundamental form
  1070. double I00 = diffu.dot(diffu); // guaranteed non-neg
  1071. double I01 = diffu.dot(diffv); // I01 = I10
  1072. double I11 = diffv.dot(diffv); // guaranteed non-neg
  1073. // eigenvalues of a 2x2 matrix
  1074. // [a00 a01]
  1075. // [a10 a11]
  1076. // 1/2 * [ (a00 + a11) +/- sqrt((a00 - a11)^2 + 4 a01 a10) ]
  1077. double trI = I00 + I11; // guaranteed non-neg
  1078. double diffDiag = I00 - I11; // guaranteed non-neg
  1079. double sqrtDet = sqrt(std::max(0.0, diffDiag*diffDiag +
  1080. 4 * I01 * I01)); // guaranteed non-neg
  1081. double sig1 = 0.5 * (trI + sqrtDet); // higher singular value
  1082. double sig2 = 0.5 * (trI - sqrtDet); // lower singular value
  1083. // Avoid sig2 < 0 due to numerical error
  1084. if (fabs(sig2) < 1.0e-8)
  1085. sig2 = 0;
  1086. assert(sig1 >= 0);
  1087. assert(sig2 >= 0);
  1088. if (sig2 < 0) {
  1089. printf("Distortion will be NaN! sig1^2 is negative (%lg)\n",
  1090. sig2);
  1091. }
  1092. // The singular values of the Jacobian are the sqrts of the
  1093. // eigenvalues of the first fundamental form.
  1094. sig1 = sqrt(sig1);
  1095. sig2 = sqrt(sig2);
  1096. // distortion
  1097. double tao = IsFlipped(f,WUV) ? -1 : 1;
  1098. double factor = tao / h;
  1099. double lam = fabs(factor * sig1 - 1) + fabs(factor * sig2 - 1);
  1100. return lam;
  1101. }
  1102. else {
  1103. return 10; // something "large"
  1104. }
  1105. }
  1106. ////////////////////////////////////////////////////////////////////////////
  1107. // Approximate the distortion Laplacian using a uniform Laplacian on
  1108. // the dual mesh:
  1109. // ___________
  1110. // \-1 / \-1 /
  1111. // \ / 3 \ /
  1112. // \-----/
  1113. // \-1 /
  1114. // \ /
  1115. //
  1116. // @param[in] f facet on which to compute distortion Laplacian
  1117. // @param[in] h scaling factor applied to cross field
  1118. // @return distortion Laplacian for f
  1119. ///////////////////////////////////////////////////////////////////////////
  1120. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1121. IGL_INLINE double igl::copyleft::comiso::MIQ_class<DerivedV, DerivedF, DerivedU>::LaplaceDistortion(const int f, double h, const Eigen::MatrixXd& WUV)
  1122. {
  1123. double mydist = Distortion(f, h, WUV);
  1124. double lapl=0;
  1125. for (int i=0;i<3;i++)
  1126. {
  1127. if (TT(f,i) != -1)
  1128. lapl += (mydist - Distortion(TT(f,i), h, WUV));
  1129. }
  1130. return lapl;
  1131. }
  1132. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1133. IGL_INLINE bool igl::copyleft::comiso::MIQ_class<DerivedV, DerivedF, DerivedU>::updateStiffeningJacobianDistorsion(double grad_size, const Eigen::MatrixXd& WUV)
  1134. {
  1135. bool flipped = NumFlips(WUV)>0;
  1136. if (!flipped)
  1137. return false;
  1138. double maxL=0;
  1139. double maxD=0;
  1140. if (flipped)
  1141. {
  1142. const double c = 1.0;
  1143. const double d = 5.0;
  1144. for (unsigned int i = 0; i < Fcut.rows(); ++i)
  1145. {
  1146. double dist=Distortion(i,grad_size,WUV);
  1147. if (dist > maxD)
  1148. maxD=dist;
  1149. double absLap=fabs(LaplaceDistortion(i, grad_size,WUV));
  1150. if (absLap > maxL)
  1151. maxL = absLap;
  1152. double stiffDelta = std::min(c * absLap, d);
  1153. stiffnessVector[i]+=stiffDelta;
  1154. }
  1155. }
  1156. printf("Maximum Distorsion %4.4f \n",maxD);
  1157. printf("Maximum Laplacian %4.4f \n",maxL);
  1158. return flipped;
  1159. }
  1160. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1161. IGL_INLINE bool igl::copyleft::comiso::MIQ_class<DerivedV, DerivedF, DerivedU>::IsFlipped(const Eigen::Vector2d &uv0,
  1162. const Eigen::Vector2d &uv1,
  1163. const Eigen::Vector2d &uv2)
  1164. {
  1165. Eigen::Vector2d e0 = (uv1-uv0);
  1166. Eigen::Vector2d e1 = (uv2-uv0);
  1167. double Area = e0(0)*e1(1) - e0(1)*e1(0);
  1168. return (Area<=0);
  1169. }
  1170. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1171. IGL_INLINE bool igl::copyleft::comiso::MIQ_class<DerivedV, DerivedF, DerivedU>::IsFlipped(
  1172. const int i, const Eigen::MatrixXd& WUV)
  1173. {
  1174. Eigen::Vector2d uv0,uv1,uv2;
  1175. uv0 << WUV(i,0), WUV(i,1);
  1176. uv1 << WUV(i,2), WUV(i,3);
  1177. uv2 << WUV(i,4), WUV(i,5);
  1178. return (IsFlipped(uv0,uv1,uv2));
  1179. }
  1180. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1181. IGL_INLINE void igl::copyleft::comiso::miq(
  1182. const Eigen::PlainObjectBase<DerivedV> &V,
  1183. const Eigen::PlainObjectBase<DerivedF> &F,
  1184. const Eigen::PlainObjectBase<DerivedV> &PD1_combed,
  1185. const Eigen::PlainObjectBase<DerivedV> &PD2_combed,
  1186. const Eigen::Matrix<int, Eigen::Dynamic, 3> &mismatch,
  1187. const Eigen::Matrix<int, Eigen::Dynamic, 1> &singular,
  1188. const Eigen::Matrix<int, Eigen::Dynamic, 3> &seams,
  1189. Eigen::PlainObjectBase<DerivedU> &UV,
  1190. Eigen::PlainObjectBase<DerivedF> &FUV,
  1191. double gradientSize,
  1192. double stiffness,
  1193. bool directRound,
  1194. unsigned int iter,
  1195. unsigned int localIter,
  1196. bool doRound,
  1197. bool singularityRound,
  1198. const std::vector<int> &roundVertices,
  1199. const std::vector<std::vector<int>> &hardFeatures)
  1200. {
  1201. gradientSize = gradientSize/(V.colwise().maxCoeff()-V.colwise().minCoeff()).norm();
  1202. igl::copyleft::comiso::MIQ_class<DerivedV, DerivedF, DerivedU> miq(V,
  1203. F,
  1204. PD1_combed,
  1205. PD2_combed,
  1206. mismatch,
  1207. singular,
  1208. seams,
  1209. UV,
  1210. FUV,
  1211. gradientSize,
  1212. stiffness,
  1213. directRound,
  1214. iter,
  1215. localIter,
  1216. doRound,
  1217. singularityRound,
  1218. roundVertices,
  1219. hardFeatures);
  1220. miq.extractUV(UV,FUV);
  1221. }
  1222. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1223. IGL_INLINE void igl::copyleft::comiso::miq(
  1224. const Eigen::PlainObjectBase<DerivedV> &V,
  1225. const Eigen::PlainObjectBase<DerivedF> &F,
  1226. const Eigen::PlainObjectBase<DerivedV> &PD1,
  1227. const Eigen::PlainObjectBase<DerivedV> &PD2,
  1228. Eigen::PlainObjectBase<DerivedU> &UV,
  1229. Eigen::PlainObjectBase<DerivedF> &FUV,
  1230. double gradientSize,
  1231. double stiffness,
  1232. bool directRound,
  1233. unsigned int iter,
  1234. unsigned int localIter,
  1235. bool doRound,
  1236. bool singularityRound,
  1237. const std::vector<int> &roundVertices,
  1238. const std::vector<std::vector<int>> &hardFeatures)
  1239. {
  1240. DerivedV BIS1, BIS2;
  1241. igl::compute_frame_field_bisectors(V, F, PD1, PD2, BIS1, BIS2);
  1242. DerivedV BIS1_combed, BIS2_combed;
  1243. igl::comb_cross_field(V, F, BIS1, BIS2, BIS1_combed, BIS2_combed);
  1244. DerivedF Handle_MMatch;
  1245. igl::cross_field_mismatch(V, F, BIS1_combed, BIS2_combed, true, Handle_MMatch);
  1246. Eigen::Matrix<int, Eigen::Dynamic, 1> isSingularity, singularityIndex;
  1247. igl::find_cross_field_singularities(V, F, Handle_MMatch, isSingularity, singularityIndex);
  1248. Eigen::Matrix<int, Eigen::Dynamic, 3> Handle_Seams;
  1249. igl::cut_mesh_from_singularities(V, F, Handle_MMatch, Handle_Seams);
  1250. DerivedV PD1_combed, PD2_combed;
  1251. igl::comb_frame_field(V, F, PD1, PD2, BIS1_combed, BIS2_combed, PD1_combed, PD2_combed);
  1252. igl::copyleft::comiso::miq(V,
  1253. F,
  1254. PD1_combed,
  1255. PD2_combed,
  1256. Handle_MMatch,
  1257. isSingularity,
  1258. Handle_Seams,
  1259. UV,
  1260. FUV,
  1261. gradientSize,
  1262. stiffness,
  1263. directRound,
  1264. iter,
  1265. localIter,
  1266. doRound,
  1267. singularityRound,
  1268. roundVertices,
  1269. hardFeatures);
  1270. }
  1271. #ifdef IGL_STATIC_LIBRARY
  1272. // Explicit template instantiation
  1273. template void igl::copyleft::comiso::miq<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, Eigen::Matrix<double, -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 &, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const &, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const &, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > &, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > &, double, double, bool, unsigned int, unsigned int, bool, bool, const std::vector<int> &, const std::vector<std::vector<int>> &);
  1274. template void igl::copyleft::comiso::miq<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::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const &, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const &, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const &, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const &, Eigen::Matrix<int, -1, 3, 0, -1, 3> const &, Eigen::Matrix<int, -1, 1, 0, -1, 1> const &, Eigen::Matrix<int, -1, 3, 0, -1, 3> const &, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > &, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > &, double, double, bool, unsigned int, unsigned int, bool, bool, const std::vector<int> &, const std::vector<std::vector<int>> &);
  1275. template void igl::copyleft::comiso::miq<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::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const &, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const &, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const &, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const &, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > &, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > &, double, double, bool, unsigned int, unsigned int, bool, bool, const std::vector<int> &, const std::vector<std::vector<int>> &);
  1276. #endif