miq.cpp 53 KB

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