miq.cpp 53 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365136613671368136913701371137213731374137513761377137813791380138113821383138413851386138713881389139013911392139313941395139613971398139914001401140214031404140514061407140814091410141114121413141414151416141714181419142014211422142314241425142614271428142914301431143214331434143514361437143814391440144114421443144414451446144714481449145014511452145314541455145614571458145914601461146214631464146514661467146814691470147114721473147414751476147714781479148014811482148314841485148614871488148914901491149214931494149514961497149814991500150115021503150415051506150715081509151015111512151315141515151615171518151915201521152215231524152515261527152815291530153115321533
  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 <igl/copyleft/comiso/miq.h>
  9. #include <igl/local_basis.h>
  10. #include <igl/triangle_triangle_adjacency.h>
  11. #include <igl/cut_mesh.h>
  12. // includes for VertexIndexing
  13. #include <igl/HalfEdgeIterator.h>
  14. #include <igl/is_border_vertex.h>
  15. #include <igl/vertex_triangle_adjacency.h>
  16. // includes for PoissonSolver
  17. #include <igl/slice_into.h>
  18. #include <igl/grad.h>
  19. #include <igl/cotmatrix.h>
  20. #include <igl/doublearea.h>
  21. #include <gmm/gmm.h>
  22. #include <CoMISo/Solver/ConstrainedSolver.hh>
  23. #include <CoMISo/Solver/MISolver.hh>
  24. #include <CoMISo/Solver/GMM_Tools.hh>
  25. //
  26. #include <igl/cross_field_missmatch.h>
  27. #include <igl/comb_frame_field.h>
  28. #include <igl/comb_cross_field.h>
  29. #include <igl/cut_mesh_from_singularities.h>
  30. #include <igl/find_cross_field_singularities.h>
  31. #include <igl/compute_frame_field_bisectors.h>
  32. #include <igl/rotate_vectors.h>
  33. #ifndef NDEBUG
  34. #include <fstream>
  35. #endif
  36. #include <iostream>
  37. #include <igl/matlab_format.h>
  38. #define DEBUGPRINT 0
  39. namespace igl {
  40. namespace copyleft {
  41. namespace comiso {
  42. struct SeamInfo
  43. {
  44. int v0,v0p;
  45. int integerVar;
  46. unsigned char MMatch;
  47. IGL_INLINE SeamInfo(int _v0,
  48. int _v0p,
  49. int _MMatch,
  50. int _integerVar);
  51. IGL_INLINE SeamInfo(const SeamInfo &S1);
  52. };
  53. struct MeshSystemInfo
  54. {
  55. ////number of vertices variables
  56. int num_vert_variables;
  57. ///num of integer for cuts
  58. int num_integer_cuts;
  59. ///this are used for drawing purposes
  60. std::vector<SeamInfo> EdgeSeamInfo;
  61. };
  62. template <typename DerivedV, typename DerivedF>
  63. class VertexIndexing
  64. {
  65. public:
  66. // Input:
  67. const Eigen::PlainObjectBase<DerivedV> &V;
  68. const Eigen::PlainObjectBase<DerivedF> &F;
  69. const Eigen::PlainObjectBase<DerivedV> &Vcut;
  70. const Eigen::PlainObjectBase<DerivedF> &Fcut;
  71. const Eigen::PlainObjectBase<DerivedF> &TT;
  72. const Eigen::PlainObjectBase<DerivedF> &TTi;
  73. const Eigen::Matrix<int, Eigen::Dynamic, 3> &Handle_MMatch;
  74. const Eigen::Matrix<int, Eigen::Dynamic, 1> &Handle_Singular; // bool
  75. const Eigen::Matrix<int, Eigen::Dynamic, 3> &Handle_Seams; // 3 bool
  76. ///this handle for mesh TODO: move with the other global variables
  77. MeshSystemInfo Handle_SystemInfo;
  78. IGL_INLINE VertexIndexing(const Eigen::PlainObjectBase<DerivedV> &_V,
  79. const Eigen::PlainObjectBase<DerivedF> &_F,
  80. const Eigen::PlainObjectBase<DerivedV> &_Vcut,
  81. const Eigen::PlainObjectBase<DerivedF> &_Fcut,
  82. const Eigen::PlainObjectBase<DerivedF> &_TT,
  83. const Eigen::PlainObjectBase<DerivedF> &_TTi,
  84. const Eigen::Matrix<int, Eigen::Dynamic, 3> &_Handle_MMatch,
  85. const Eigen::Matrix<int, Eigen::Dynamic, 1> &_Handle_Singular,
  86. const Eigen::Matrix<int, Eigen::Dynamic, 3> &_Handle_Seams
  87. );
  88. // provide information about every vertex per seam
  89. IGL_INLINE void InitSeamInfo();
  90. private:
  91. struct VertexInfo{
  92. int v; // vertex index (according to V)
  93. int f0, k0; // face and local edge information of the edge that connects this vertex to the previous vertex (previous in the vector)
  94. int f1, k1; // face and local edge information of the other face corresponding to the same edge
  95. VertexInfo(int _v, int _f0, int _k0, int _f1, int _k1) :
  96. v(_v), f0(_f0), k0(_k0), f1(_f1), k1(_k1){}
  97. bool operator==(VertexInfo const& other){
  98. return other.v == v;
  99. }
  100. };
  101. IGL_INLINE void GetSeamInfo(const int f0,
  102. const int f1,
  103. const int indexE,
  104. int &v0,int &v1,
  105. int &v0p,int &v1p,
  106. unsigned char &_MMatch);
  107. IGL_INLINE std::vector<std::vector<VertexInfo> > GetVerticesPerSeam();
  108. };
  109. template <typename DerivedV, typename DerivedF>
  110. class PoissonSolver
  111. {
  112. public:
  113. IGL_INLINE void SolvePoisson(Eigen::VectorXd Stiffness,
  114. double vector_field_scale=0.1f,
  115. double grid_res=1.f,
  116. bool direct_round=true,
  117. int localIter=0,
  118. bool _integer_rounding=true,
  119. bool _singularity_rounding=true,
  120. std::vector<int> roundVertices = std::vector<int>(),
  121. std::vector<std::vector<int> > hardFeatures = std::vector<std::vector<int> >());
  122. IGL_INLINE PoissonSolver(const Eigen::PlainObjectBase<DerivedV> &_V,
  123. const Eigen::PlainObjectBase<DerivedF> &_F,
  124. const Eigen::PlainObjectBase<DerivedV> &_Vcut,
  125. const Eigen::PlainObjectBase<DerivedF> &_Fcut,
  126. const Eigen::PlainObjectBase<DerivedF> &_TT,
  127. const Eigen::PlainObjectBase<DerivedF> &_TTi,
  128. const Eigen::PlainObjectBase<DerivedV> &_PD1,
  129. const Eigen::PlainObjectBase<DerivedV> &_PD2,
  130. const Eigen::Matrix<int, Eigen::Dynamic, 1>&_Handle_Singular,
  131. const MeshSystemInfo &_Handle_SystemInfo
  132. );
  133. const Eigen::PlainObjectBase<DerivedV> &V;
  134. const Eigen::PlainObjectBase<DerivedF> &F;
  135. const Eigen::PlainObjectBase<DerivedV> &Vcut;
  136. const Eigen::PlainObjectBase<DerivedF> &Fcut;
  137. const Eigen::PlainObjectBase<DerivedF> &TT;
  138. const Eigen::PlainObjectBase<DerivedF> &TTi;
  139. const Eigen::PlainObjectBase<DerivedV> &PD1;
  140. const Eigen::PlainObjectBase<DerivedV> &PD2;
  141. const Eigen::Matrix<int, Eigen::Dynamic, 1> &Handle_Singular; // bool
  142. const MeshSystemInfo &Handle_SystemInfo;
  143. // Internal:
  144. Eigen::VectorXd Handle_Stiffness;
  145. std::vector<std::vector<int> > VF;
  146. std::vector<std::vector<int> > VFi;
  147. Eigen::MatrixXd UV; // this is probably useless
  148. // Output:
  149. // per wedge UV coordinates, 6 coordinates (1 face) per row
  150. Eigen::MatrixXd WUV;
  151. // per vertex UV coordinates, Vcut.rows() x 2
  152. Eigen::MatrixXd UV_out;
  153. // Matrices
  154. Eigen::SparseMatrix<double> Lhs;
  155. Eigen::SparseMatrix<double> Constraints;
  156. Eigen::VectorXd rhs;
  157. Eigen::VectorXd constraints_rhs;
  158. ///vector of unknowns
  159. std::vector< double > X;
  160. ////REAL PART
  161. ///number of fixed vertex
  162. unsigned int n_fixed_vars;
  163. ///the number of REAL variables for vertices
  164. unsigned int n_vert_vars;
  165. ///total number of variables of the system,
  166. ///do not consider constraints, but consider integer vars
  167. unsigned int num_total_vars;
  168. //////INTEGER PART
  169. ///the total number of integer variables
  170. unsigned int n_integer_vars;
  171. ///CONSTRAINT PART
  172. ///number of cuts constraints
  173. unsigned int num_cut_constraint;
  174. // number of user-defined constraints
  175. unsigned int num_userdefined_constraint;
  176. ///total number of constraints equations
  177. unsigned int num_constraint_equations;
  178. ///vector of blocked vertices
  179. std::vector<int> Hard_constraints;
  180. ///vector of indexes to round
  181. std::vector<int> ids_to_round;
  182. ///vector of indexes to round
  183. std::vector<std::vector<int > > userdefined_constraints;
  184. ///boolean that is true if rounding to integer is needed
  185. bool integer_rounding;
  186. ///START COMMON MATH FUNCTIONS
  187. ///return the complex encoding the rotation
  188. ///for a given missmatch interval
  189. IGL_INLINE std::complex<double> GetRotationComplex(int interval);
  190. ///END COMMON MATH FUNCTIONS
  191. ///START FIXING VERTICES
  192. ///set a given vertex as fixed
  193. IGL_INLINE void AddFixedVertex(int v);
  194. ///find vertex to fix in case we're using
  195. ///a vector field NB: multiple components not handled
  196. IGL_INLINE void FindFixedVertField();
  197. ///find hard constraint depending if using or not
  198. ///a vector field
  199. IGL_INLINE void FindFixedVert();
  200. IGL_INLINE int GetFirstVertexIndex(int v);
  201. ///fix the vertices which are flagged as fixed
  202. IGL_INLINE void FixBlockedVertex();
  203. ///END FIXING VERTICES
  204. ///HANDLING SINGULARITY
  205. //set the singularity round to integer location
  206. IGL_INLINE void AddSingularityRound();
  207. IGL_INLINE void AddToRoundVertices(std::vector<int> ids);
  208. ///START GENERIC SYSTEM FUNCTIONS
  209. //build the laplacian matrix cyclyng over all rangemaps
  210. //and over all faces
  211. IGL_INLINE void BuildLaplacianMatrix(double vfscale=1);
  212. ///find different sized of the system
  213. IGL_INLINE void FindSizes();
  214. IGL_INLINE void AllocateSystem();
  215. ///intitialize the whole matrix
  216. IGL_INLINE void InitMatrix();
  217. ///map back coordinates after that
  218. ///the system has been solved
  219. IGL_INLINE void MapCoords();
  220. ///END GENERIC SYSTEM FUNCTIONS
  221. ///set the constraints for the inter-range cuts
  222. IGL_INLINE void BuildSeamConstraintsExplicitTranslation();
  223. ///set the constraints for the inter-range cuts
  224. IGL_INLINE void BuildUserDefinedConstraints();
  225. ///call of the mixed integer solver
  226. IGL_INLINE void MixedIntegerSolve(double cone_grid_res=1,
  227. bool direct_round=true,
  228. int localIter=0);
  229. IGL_INLINE void clearUserConstraint();
  230. IGL_INLINE void addSharpEdgeConstraint(int fid, int vid);
  231. };
  232. template <typename DerivedV, typename DerivedF, typename DerivedU>
  233. class MIQ_class
  234. {
  235. private:
  236. const Eigen::PlainObjectBase<DerivedV> &V;
  237. const Eigen::PlainObjectBase<DerivedF> &F;
  238. Eigen::PlainObjectBase<DerivedV> Vcut;
  239. Eigen::PlainObjectBase<DerivedF> Fcut;
  240. Eigen::MatrixXd UV_out;
  241. Eigen::PlainObjectBase<DerivedF> FUV_out;
  242. // internal
  243. Eigen::PlainObjectBase<DerivedF> TT;
  244. Eigen::PlainObjectBase<DerivedF> TTi;
  245. // Stiffness per face
  246. Eigen::VectorXd Handle_Stiffness;
  247. Eigen::PlainObjectBase<DerivedV> B1, B2, B3;
  248. public:
  249. IGL_INLINE MIQ_class(const Eigen::PlainObjectBase<DerivedV> &V_,
  250. const Eigen::PlainObjectBase<DerivedF> &F_,
  251. const Eigen::PlainObjectBase<DerivedV> &PD1_combed,
  252. const Eigen::PlainObjectBase<DerivedV> &PD2_combed,
  253. const Eigen::Matrix<int, Eigen::Dynamic, 3> &Handle_MMatch,
  254. const Eigen::Matrix<int, Eigen::Dynamic, 1> &Handle_Singular,
  255. const Eigen::Matrix<int, Eigen::Dynamic, 3> &Handle_Seams,
  256. Eigen::PlainObjectBase<DerivedU> &UV,
  257. Eigen::PlainObjectBase<DerivedF> &FUV,
  258. double GradientSize = 30.0,
  259. double Stiffness = 5.0,
  260. bool DirectRound = false,
  261. int iter = 5,
  262. int localIter = 5,
  263. bool DoRound = true,
  264. bool SingularityRound=true,
  265. std::vector<int> roundVertices = std::vector<int>(),
  266. std::vector<std::vector<int> > hardFeatures = std::vector<std::vector<int> >());
  267. IGL_INLINE void extractUV(Eigen::PlainObjectBase<DerivedU> &UV_out,
  268. Eigen::PlainObjectBase<DerivedF> &FUV_out);
  269. private:
  270. IGL_INLINE int NumFlips(const Eigen::MatrixXd& WUV);
  271. IGL_INLINE double Distortion(int f, double h, const Eigen::MatrixXd& WUV);
  272. IGL_INLINE double LaplaceDistortion(const int f, double h, const Eigen::MatrixXd& WUV);
  273. IGL_INLINE bool updateStiffeningJacobianDistorsion(double grad_size, const Eigen::MatrixXd& WUV);
  274. IGL_INLINE bool IsFlipped(const Eigen::Vector2d &uv0,
  275. const Eigen::Vector2d &uv1,
  276. const Eigen::Vector2d &uv2);
  277. IGL_INLINE bool IsFlipped(const int i, const Eigen::MatrixXd& WUV);
  278. };
  279. };
  280. };
  281. }
  282. IGL_INLINE igl::copyleft::comiso::SeamInfo::SeamInfo(int _v0,
  283. int _v0p,
  284. int _MMatch,
  285. int _integerVar)
  286. {
  287. v0=_v0;
  288. v0p=_v0p;
  289. integerVar=_integerVar;
  290. MMatch=_MMatch;
  291. }
  292. IGL_INLINE igl::copyleft::comiso::SeamInfo::SeamInfo(const SeamInfo &S1)
  293. {
  294. v0=S1.v0;
  295. v0p=S1.v0p;
  296. integerVar=S1.integerVar;
  297. MMatch=S1.MMatch;
  298. }
  299. template <typename DerivedV, typename DerivedF>
  300. IGL_INLINE igl::copyleft::comiso::VertexIndexing<DerivedV, DerivedF>::VertexIndexing(const Eigen::PlainObjectBase<DerivedV> &_V,
  301. const Eigen::PlainObjectBase<DerivedF> &_F,
  302. const Eigen::PlainObjectBase<DerivedV> &_Vcut,
  303. const Eigen::PlainObjectBase<DerivedF> &_Fcut,
  304. const Eigen::PlainObjectBase<DerivedF> &_TT,
  305. const Eigen::PlainObjectBase<DerivedF> &_TTi,
  306. const Eigen::Matrix<int, Eigen::Dynamic, 3> &_Handle_MMatch,
  307. const Eigen::Matrix<int, Eigen::Dynamic, 1> &_Handle_Singular,
  308. const Eigen::Matrix<int, Eigen::Dynamic, 3> &_Handle_Seams
  309. ):
  310. V(_V),
  311. F(_F),
  312. Vcut(_Vcut),
  313. Fcut(_Fcut),
  314. TT(_TT),
  315. TTi(_TTi),
  316. Handle_MMatch(_Handle_MMatch),
  317. Handle_Singular(_Handle_Singular),
  318. Handle_Seams(_Handle_Seams)
  319. {
  320. #ifdef DEBUG_PRINT
  321. cerr<<igl::matlab_format(Handle_Seams,"Handle_Seams");
  322. #endif
  323. Handle_SystemInfo.num_vert_variables=Vcut.rows();
  324. Handle_SystemInfo.num_integer_cuts=0;
  325. }
  326. template <typename DerivedV, typename DerivedF>
  327. IGL_INLINE void igl::copyleft::comiso::VertexIndexing<DerivedV, DerivedF>::GetSeamInfo(const int f0,
  328. const int f1,
  329. const int indexE,
  330. int &v0,int &v1,
  331. int &v0p,int &v1p,
  332. unsigned char &_MMatch)
  333. {
  334. int edgef0 = indexE;
  335. v0 = Fcut(f0,edgef0);
  336. v1 = Fcut(f0,(edgef0+1)%3);
  337. ////get the index on opposite side
  338. assert(TT(f0,edgef0) == f1);
  339. int edgef1 = TTi(f0,edgef0);
  340. v1p = Fcut(f1,edgef1);
  341. v0p = Fcut(f1,(edgef1+1)%3);
  342. _MMatch = Handle_MMatch(f0,edgef0);
  343. assert(F(f0,edgef0) == F(f1,((edgef1+1)%3)));
  344. assert(F(f0,((edgef0+1)%3)) == F(f1,edgef1));
  345. }
  346. template <typename DerivedV, typename DerivedF>
  347. IGL_INLINE std::vector<std::vector<typename igl::copyleft::comiso::VertexIndexing<DerivedV, DerivedF>::VertexInfo> > igl::copyleft::comiso::VertexIndexing<DerivedV, DerivedF>::GetVerticesPerSeam()
  348. {
  349. // Return value
  350. std::vector<std::vector<VertexInfo> >verticesPerSeam;
  351. // for every vertex, keep track of their adjacent vertices on seams.
  352. // regular vertices have two neighbors on a seam, start- and endvertices may have any other numbers of neighbors (e.g. 1 or 3)
  353. std::vector<std::list<VertexInfo> > VVSeam(V.rows());
  354. Eigen::MatrixXi F_hit = Eigen::MatrixXi::Zero(F.rows(), 3);
  355. for (unsigned int f=0; f<F.rows();f++)
  356. {
  357. int f0 = f;
  358. for(int k0=0; k0<3; k0++){
  359. int f1 = TT(f0,k0);
  360. if(f1 == -1)
  361. continue;
  362. bool seam = Handle_Seams(f0,k0);
  363. if (seam && F_hit(f0,k0) == 0)
  364. {
  365. int v0 = F(f0, k0);
  366. int v1 = F(f0, (k0+1)%3);
  367. int k1 = TTi(f0,k0);
  368. VVSeam[v0].push_back(VertexInfo(v1, f0, k0, f1, k1));
  369. VVSeam[v1].push_back(VertexInfo(v0, f0, k0, f1, k1));
  370. F_hit(f0, k0) = 1;
  371. F_hit(f1, k1) = 1;
  372. }
  373. }
  374. }
  375. // Find start vertices, i.e. vertices that start or end a seam branch
  376. std::vector<int> startVertexIndices;
  377. std::vector<bool> isStartVertex(V.rows());
  378. for (unsigned int i=0;i<V.rows();i++)
  379. {
  380. isStartVertex[i] = false;
  381. // vertices with two neighbors are regular vertices, unless the vertex is a singularity, in which case it qualifies as a start vertex
  382. if (VVSeam[i].size() > 0 && VVSeam[i].size() != 2 || Handle_Singular(i) == true)
  383. {
  384. startVertexIndices.push_back(i);
  385. isStartVertex[i] = true;
  386. }
  387. }
  388. // For each startVertex, walk along its seam
  389. for (unsigned int i=0;i<startVertexIndices.size();i++)
  390. {
  391. auto startVertexNeighbors = &VVSeam[startVertexIndices[i]];
  392. const int neighborSize = startVertexNeighbors->size();
  393. // explore every seam to which this vertex is a start vertex
  394. // note: a vertex can never be a start vertex and a regular vertex simultaneously
  395. for (unsigned int j=0;j<neighborSize;j++)
  396. {
  397. std::vector<VertexInfo> thisSeam; // temporary container
  398. // Create vertexInfo struct for start vertex
  399. auto startVertex = VertexInfo(startVertexIndices[i], -1, -1, -1, -1);// -1 values are arbitrary (will never be used)
  400. auto currentVertex = startVertex;
  401. // Add start vertex to the seam
  402. thisSeam.push_back(currentVertex);
  403. // advance on the seam
  404. auto currentVertexNeighbors = startVertexNeighbors;
  405. auto nextVertex = currentVertexNeighbors->front();
  406. currentVertexNeighbors->pop_front();
  407. auto prevVertex = startVertex; // bogus initialization to get the type
  408. while (true)
  409. {
  410. // move to the next vertex
  411. prevVertex = currentVertex;
  412. currentVertex = nextVertex;
  413. currentVertexNeighbors = &VVSeam[nextVertex.v];
  414. // add current vertex to this seam
  415. thisSeam.push_back(currentVertex);
  416. // remove the previous vertex
  417. auto it = std::find(currentVertexNeighbors->begin(), currentVertexNeighbors->end(), prevVertex);
  418. assert(it != currentVertexNeighbors->end());
  419. currentVertexNeighbors->erase(it);
  420. if (currentVertexNeighbors->size() == 1 && !isStartVertex[currentVertex.v])
  421. {
  422. nextVertex = currentVertexNeighbors->front();
  423. currentVertexNeighbors->pop_front();
  424. }
  425. else
  426. break;
  427. }
  428. verticesPerSeam.push_back(thisSeam);
  429. }
  430. }
  431. return verticesPerSeam;
  432. }
  433. template <typename DerivedV, typename DerivedF>
  434. IGL_INLINE void igl::copyleft::comiso::VertexIndexing<DerivedV, DerivedF>::InitSeamInfo()
  435. {
  436. auto verticesPerSeam = GetVerticesPerSeam();
  437. Handle_SystemInfo.EdgeSeamInfo.clear();
  438. int integerVar = 0;
  439. // Loop over each seam
  440. for(auto seam : verticesPerSeam){
  441. //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.
  442. int priorVertexIdx;
  443. if(seam.size() > 2){
  444. auto v1 = seam[1];
  445. auto v2 = seam[2];
  446. 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)){
  447. priorVertexIdx = Fcut(v1.f0, v1.k0);
  448. }
  449. else{
  450. priorVertexIdx = Fcut(v1.f1, v1.k1);
  451. 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));
  452. }
  453. }
  454. else{
  455. auto v1 = seam[1];
  456. priorVertexIdx = Fcut(v1.f0, v1.k0);
  457. }
  458. // Loop over each vertex of the seam
  459. for(auto it=seam.begin()+1; it != seam.end(); ++it){
  460. auto vertex = *it;
  461. // choose the correct side of the seam
  462. int f,k,ff,kk;
  463. if(priorVertexIdx == Fcut(vertex.f0, vertex.k0)){
  464. f = vertex.f0; ff = vertex.f1;
  465. k = vertex.k0; kk = vertex.k1;
  466. }
  467. else{
  468. f = vertex.f1; ff = vertex.f0;
  469. k = vertex.k1; kk = vertex.k0;
  470. assert(priorVertexIdx == Fcut(vertex.f1, vertex.k1));
  471. }
  472. int vtx0,vtx0p,vtx1,vtx1p;
  473. unsigned char MM;
  474. GetSeamInfo(f,ff,k,vtx0,vtx1,vtx0p,vtx1p,MM);
  475. Handle_SystemInfo.EdgeSeamInfo.push_back(SeamInfo(vtx0,vtx0p,MM,integerVar));
  476. if(it == seam.end() -1){
  477. Handle_SystemInfo.EdgeSeamInfo.push_back(SeamInfo(vtx1,vtx1p,MM,integerVar));
  478. }
  479. priorVertexIdx = vtx1;
  480. }
  481. // use the same integer for each seam
  482. integerVar++;
  483. }
  484. Handle_SystemInfo.num_integer_cuts = integerVar;
  485. #ifndef NDEBUG
  486. int totalNVerticesOnSeams = 0;
  487. for(auto seam : verticesPerSeam){
  488. totalNVerticesOnSeams += seam.size();
  489. }
  490. assert(Handle_SystemInfo.EdgeSeamInfo.size() == totalNVerticesOnSeams);
  491. #endif
  492. }
  493. template <typename DerivedV, typename DerivedF>
  494. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::SolvePoisson(Eigen::VectorXd Stiffness,
  495. double vector_field_scale,
  496. double grid_res,
  497. bool direct_round,
  498. int localIter,
  499. bool _integer_rounding,
  500. bool _singularity_rounding,
  501. std::vector<int> roundVertices,
  502. std::vector<std::vector<int> > hardFeatures)
  503. {
  504. Handle_Stiffness = Stiffness;
  505. //initialization of flags and data structures
  506. integer_rounding=_integer_rounding;
  507. ids_to_round.clear();
  508. clearUserConstraint();
  509. // copy the user constraints number
  510. for (size_t i = 0; i < hardFeatures.size(); ++i)
  511. {
  512. addSharpEdgeConstraint(hardFeatures[i][0],hardFeatures[i][1]);
  513. }
  514. ///Initializing Matrix
  515. int t0=clock();
  516. ///initialize the matrix ALLOCATING SPACE
  517. InitMatrix();
  518. if (DEBUGPRINT)
  519. printf("\n ALLOCATED THE MATRIX \n");
  520. ///build the laplacian system
  521. BuildLaplacianMatrix(vector_field_scale);
  522. // add seam constraints
  523. BuildSeamConstraintsExplicitTranslation();
  524. // add user defined constraints
  525. BuildUserDefinedConstraints();
  526. ////add the lagrange multiplier
  527. FixBlockedVertex();
  528. if (DEBUGPRINT)
  529. printf("\n BUILT THE MATRIX \n");
  530. if (integer_rounding)
  531. AddToRoundVertices(roundVertices);
  532. if (_singularity_rounding)
  533. AddSingularityRound();
  534. int t1=clock();
  535. if (DEBUGPRINT) printf("\n time:%d \n",t1-t0);
  536. if (DEBUGPRINT) printf("\n SOLVING \n");
  537. MixedIntegerSolve(grid_res,direct_round,localIter);
  538. int t2=clock();
  539. if (DEBUGPRINT) printf("\n time:%d \n",t2-t1);
  540. if (DEBUGPRINT) printf("\n ASSIGNING COORDS \n");
  541. MapCoords();
  542. int t3=clock();
  543. if (DEBUGPRINT) printf("\n time:%d \n",t3-t2);
  544. if (DEBUGPRINT) printf("\n FINISHED \n");
  545. }
  546. template <typename DerivedV, typename DerivedF>
  547. IGL_INLINE igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>
  548. ::PoissonSolver(const Eigen::PlainObjectBase<DerivedV> &_V,
  549. const Eigen::PlainObjectBase<DerivedF> &_F,
  550. const Eigen::PlainObjectBase<DerivedV> &_Vcut,
  551. const Eigen::PlainObjectBase<DerivedF> &_Fcut,
  552. const Eigen::PlainObjectBase<DerivedF> &_TT,
  553. const Eigen::PlainObjectBase<DerivedF> &_TTi,
  554. const Eigen::PlainObjectBase<DerivedV> &_PD1,
  555. const Eigen::PlainObjectBase<DerivedV> &_PD2,
  556. const Eigen::Matrix<int, Eigen::Dynamic, 1>&_Handle_Singular,
  557. const MeshSystemInfo &_Handle_SystemInfo
  558. ):
  559. V(_V),
  560. F(_F),
  561. Vcut(_Vcut),
  562. Fcut(_Fcut),
  563. TT(_TT),
  564. TTi(_TTi),
  565. PD1(_PD1),
  566. PD2(_PD2),
  567. Handle_Singular(_Handle_Singular),
  568. Handle_SystemInfo(_Handle_SystemInfo)
  569. {
  570. UV = Eigen::MatrixXd(V.rows(),2);
  571. WUV = Eigen::MatrixXd(F.rows(),6);
  572. UV_out = Eigen::MatrixXd(Vcut.rows(),2);
  573. igl::vertex_triangle_adjacency(V,F,VF,VFi);
  574. }
  575. ///START COMMON MATH FUNCTIONS
  576. ///return the complex encoding the rotation
  577. ///for a given missmatch interval
  578. template <typename DerivedV, typename DerivedF>
  579. IGL_INLINE std::complex<double> igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::GetRotationComplex(int interval)
  580. {
  581. assert((interval>=0)&&(interval<4));
  582. switch(interval)
  583. {
  584. case 0:return std::complex<double>(1,0);
  585. case 1:return std::complex<double>(0,1);
  586. case 2:return std::complex<double>(-1,0);
  587. default:return std::complex<double>(0,-1);
  588. }
  589. }
  590. ///END COMMON MATH FUNCTIONS
  591. ///START FIXING VERTICES
  592. ///set a given vertex as fixed
  593. template <typename DerivedV, typename DerivedF>
  594. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::AddFixedVertex(int v)
  595. {
  596. n_fixed_vars++;
  597. Hard_constraints.push_back(v);
  598. }
  599. ///find vertex to fix in case we're using
  600. ///a vector field NB: multiple components not handled
  601. template <typename DerivedV, typename DerivedF>
  602. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::FindFixedVertField()
  603. {
  604. Hard_constraints.clear();
  605. n_fixed_vars=0;
  606. //fix the first singularity
  607. for (unsigned int v=0;v<V.rows();v++)
  608. {
  609. if (Handle_Singular(v))
  610. {
  611. AddFixedVertex(v);
  612. UV.row(v) << 0,0;
  613. return;
  614. }
  615. }
  616. ///if anything fixed fix the first
  617. AddFixedVertex(0);
  618. UV.row(0) << 0,0;
  619. std::cerr << "No vertices to fix, I am fixing the first vertex to the origin!" << std::endl;
  620. }
  621. ///find hard constraint depending if using or not
  622. ///a vector field
  623. template <typename DerivedV, typename DerivedF>
  624. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::FindFixedVert()
  625. {
  626. Hard_constraints.clear();
  627. FindFixedVertField();
  628. }
  629. template <typename DerivedV, typename DerivedF>
  630. IGL_INLINE int igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::GetFirstVertexIndex(int v)
  631. {
  632. return Fcut(VF[v][0],VFi[v][0]);
  633. }
  634. ///fix the vertices which are flagged as fixed
  635. template <typename DerivedV, typename DerivedF>
  636. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::FixBlockedVertex()
  637. {
  638. int offset_row = num_cut_constraint*2;
  639. unsigned int constr_num = 0;
  640. for (unsigned int i=0;i<Hard_constraints.size();i++)
  641. {
  642. int v = Hard_constraints[i];
  643. ///get first index of the vertex that must blocked
  644. //int index=v->vertex_index[0];
  645. int index = GetFirstVertexIndex(v);
  646. ///multiply times 2 because of uv
  647. int indexvert = index*2;
  648. ///find the first free row to add the constraint
  649. int indexRow = (offset_row+constr_num*2);
  650. int indexCol = indexRow;
  651. ///add fixing constraint LHS
  652. Constraints.coeffRef(indexRow, indexvert) += 1;
  653. Constraints.coeffRef(indexRow+1,indexvert+1) += 1;
  654. ///add fixing constraint RHS
  655. constraints_rhs[indexCol] = UV(v,0);
  656. constraints_rhs[indexCol+1] = UV(v,1);
  657. constr_num++;
  658. }
  659. assert(constr_num==n_fixed_vars);
  660. }
  661. ///END FIXING VERTICES
  662. ///HANDLING SINGULARITY
  663. //set the singularity round to integer location
  664. template <typename DerivedV, typename DerivedF>
  665. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::AddSingularityRound()
  666. {
  667. for (unsigned int v=0;v<V.rows();v++)
  668. {
  669. if (Handle_Singular(v))
  670. {
  671. int index0=GetFirstVertexIndex(v);
  672. ids_to_round.push_back( index0*2 );
  673. ids_to_round.push_back((index0*2)+1);
  674. }
  675. }
  676. }
  677. template <typename DerivedV, typename DerivedF>
  678. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::AddToRoundVertices(std::vector<int> ids)
  679. {
  680. for (size_t i = 0; i < ids.size(); ++i)
  681. {
  682. if (ids[i] < 0 || ids[i] >= V.rows())
  683. std::cerr << "WARNING: Ignored round vertex constraint, vertex " << ids[i] << " does not exist in the mesh." << std::endl;
  684. int index0 = GetFirstVertexIndex(ids[i]);
  685. ids_to_round.push_back( index0*2 );
  686. ids_to_round.push_back((index0*2)+1);
  687. }
  688. }
  689. ///START GENERIC SYSTEM FUNCTIONS
  690. template <typename DerivedV, typename DerivedF>
  691. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::BuildLaplacianMatrix(double vfscale)
  692. {
  693. Eigen::VectorXi idx = Eigen::VectorXi::LinSpaced(Vcut.rows(), 0, 2*Vcut.rows()-2);
  694. Eigen::VectorXi idx2 = Eigen::VectorXi::LinSpaced(Vcut.rows(), 1, 2*Vcut.rows()-1);
  695. // get gradient matrix
  696. Eigen::SparseMatrix<double> G(Fcut.rows() * 3, Vcut.rows());
  697. igl::grad(Vcut, Fcut, G);
  698. // get triangle weights
  699. Eigen::VectorXd dblA(Fcut.rows());
  700. igl::doublearea(Vcut, Fcut, dblA);
  701. // compute intermediate result
  702. Eigen::SparseMatrix<double> G2;
  703. G2 = G.transpose() * dblA.replicate<3,1>().asDiagonal() * Handle_Stiffness.replicate<3,1>().asDiagonal();
  704. /// Compute LHS
  705. Eigen::SparseMatrix<double> Cotmatrix;
  706. Cotmatrix = 0.5 * G2 * G;
  707. igl::slice_into(Cotmatrix, idx, idx, Lhs);
  708. igl::slice_into(Cotmatrix, idx2, idx2, Lhs);
  709. /// Compute RHS
  710. // reshape nrosy vectors
  711. 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.
  712. 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.
  713. // multiply with weights
  714. Eigen::VectorXd rhs1 = G2 * u * 0.5 * vfscale;
  715. Eigen::VectorXd rhs2 = -G2 * v * 0.5 * vfscale;
  716. igl::slice_into(rhs1, idx, 1, rhs);
  717. igl::slice_into(rhs2, idx2, 1, rhs);
  718. }
  719. ///find different sized of the system
  720. template <typename DerivedV, typename DerivedF>
  721. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::FindSizes()
  722. {
  723. ///find the vertex that need to be fixed
  724. FindFixedVert();
  725. ///REAL PART
  726. n_vert_vars = Handle_SystemInfo.num_vert_variables;
  727. ///INTEGER PART
  728. ///the total number of integer variables
  729. n_integer_vars = Handle_SystemInfo.num_integer_cuts;
  730. ///CONSTRAINT PART
  731. num_cut_constraint = Handle_SystemInfo.EdgeSeamInfo.size();
  732. num_constraint_equations = num_cut_constraint * 2 + n_fixed_vars * 2 + num_userdefined_constraint;
  733. ///total variable of the system
  734. num_total_vars = (n_vert_vars+n_integer_vars) * 2;
  735. ///initialize matrix size
  736. if (DEBUGPRINT) printf("\n*** SYSTEM VARIABLES *** \n");
  737. if (DEBUGPRINT) printf("* NUM REAL VERTEX VARIABLES %d \n",n_vert_vars);
  738. if (DEBUGPRINT) printf("\n*** INTEGER VARIABLES *** \n");
  739. if (DEBUGPRINT) printf("* NUM INTEGER VARIABLES %d \n",(int)n_integer_vars);
  740. if (DEBUGPRINT) printf("\n*** CONSTRAINTS *** \n ");
  741. if (DEBUGPRINT) printf("* NUM FIXED CONSTRAINTS %d\n",n_fixed_vars);
  742. if (DEBUGPRINT) printf("* NUM CUTS CONSTRAINTS %d\n",num_cut_constraint);
  743. if (DEBUGPRINT) printf("* NUM USER DEFINED CONSTRAINTS %d\n",num_userdefined_constraint);
  744. if (DEBUGPRINT) printf("\n*** TOTAL SIZE *** \n");
  745. if (DEBUGPRINT) printf("* TOTAL VARIABLE SIZE (WITH INTEGER TRASL) %d \n",num_total_vars);
  746. if (DEBUGPRINT) printf("* TOTAL CONSTRAINTS %d \n",num_constraint_equations);
  747. }
  748. template <typename DerivedV, typename DerivedF>
  749. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::AllocateSystem()
  750. {
  751. Lhs.resize(n_vert_vars * 2, n_vert_vars * 2);
  752. Constraints.resize(num_constraint_equations, num_total_vars);
  753. rhs.resize(n_vert_vars * 2);
  754. constraints_rhs.resize(num_constraint_equations);
  755. printf("\n INITIALIZED SPARSE MATRIX OF %d x %d \n",n_vert_vars*2, n_vert_vars*2);
  756. printf("\n INITIALIZED SPARSE MATRIX OF %d x %d \n",num_constraint_equations, num_total_vars);
  757. printf("\n INITIALIZED VECTOR OF %d x 1 \n",n_vert_vars*2);
  758. printf("\n INITIALIZED VECTOR OF %d x 1 \n",num_constraint_equations);
  759. }
  760. ///intitialize the whole matrix
  761. template <typename DerivedV, typename DerivedF>
  762. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::InitMatrix()
  763. {
  764. FindSizes();
  765. AllocateSystem();
  766. }
  767. ///map back coordinates after that
  768. ///the system has been solved
  769. template <typename DerivedV, typename DerivedF>
  770. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::MapCoords()
  771. {
  772. ///map coords to faces
  773. for (unsigned int f=0;f<Fcut.rows();f++)
  774. {
  775. for (int k=0;k<3;k++)
  776. {
  777. //get the index of the variable in the system
  778. int indexUV = Fcut(f,k);
  779. ///then get U and V coords
  780. double U=X[indexUV*2];
  781. double V=X[indexUV*2+1];
  782. WUV(f,k*2 + 0) = U;
  783. WUV(f,k*2 + 1) = V;
  784. }
  785. }
  786. for(int i = 0; i < Vcut.rows(); i++){
  787. UV_out(i,0) = X[i*2];
  788. UV_out(i,1) = X[i*2+1];
  789. }
  790. }
  791. ///END GENERIC SYSTEM FUNCTIONS
  792. ///set the constraints for the inter-range cuts
  793. template <typename DerivedV, typename DerivedF>
  794. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::BuildSeamConstraintsExplicitTranslation()
  795. {
  796. ///current constraint row
  797. int constr_row = 0;
  798. for (unsigned int i=0; i<num_cut_constraint; i++)
  799. {
  800. unsigned char interval = Handle_SystemInfo.EdgeSeamInfo[i].MMatch;
  801. if (interval==1)
  802. interval=3;
  803. else
  804. if(interval==3)
  805. interval=1;
  806. int p0 = Handle_SystemInfo.EdgeSeamInfo[i].v0;
  807. int p0p = Handle_SystemInfo.EdgeSeamInfo[i].v0p;
  808. std::complex<double> rot = GetRotationComplex(interval);
  809. ///get the integer variable
  810. int integerVar = n_vert_vars + Handle_SystemInfo.EdgeSeamInfo[i].integerVar;
  811. if (integer_rounding)
  812. {
  813. ids_to_round.push_back(integerVar*2);
  814. ids_to_round.push_back(integerVar*2+1);
  815. }
  816. // cross boundary compatibility conditions
  817. Constraints.coeffRef(constr_row, 2*p0) += rot.real();
  818. Constraints.coeffRef(constr_row, 2*p0+1) += -rot.imag();
  819. Constraints.coeffRef(constr_row+1, 2*p0) += rot.imag();
  820. Constraints.coeffRef(constr_row+1, 2*p0+1) += rot.real();
  821. Constraints.coeffRef(constr_row, 2*p0p) += -1;
  822. Constraints.coeffRef(constr_row+1, 2*p0p+1) += -1;
  823. Constraints.coeffRef(constr_row, 2*integerVar) += 1;
  824. Constraints.coeffRef(constr_row+1, 2*integerVar+1) += 1;
  825. constraints_rhs[constr_row] = 0;
  826. constraints_rhs[constr_row+1] = 0;
  827. constr_row += 2;
  828. }
  829. }
  830. ///set the constraints for the inter-range cuts
  831. template <typename DerivedV, typename DerivedF>
  832. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::BuildUserDefinedConstraints()
  833. {
  834. /// the user defined constraints are at the end
  835. int offset_row = num_cut_constraint*2 + n_fixed_vars*2;
  836. ///current constraint row
  837. int constr_row = offset_row;
  838. assert(num_userdefined_constraint == userdefined_constraints.size());
  839. for (unsigned int i=0; i<num_userdefined_constraint; i++)
  840. {
  841. for (unsigned int j=0; j<userdefined_constraints[i].size()-1; ++j)
  842. {
  843. Constraints.coeffRef(constr_row, j) = userdefined_constraints[i][j];
  844. }
  845. constraints_rhs[constr_row] = userdefined_constraints[i][userdefined_constraints[i].size()-1];
  846. constr_row +=1;
  847. }
  848. }
  849. ///call of the mixed integer solver
  850. template <typename DerivedV, typename DerivedF>
  851. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::MixedIntegerSolve(double cone_grid_res,
  852. bool direct_round,
  853. int localIter)
  854. {
  855. X = std::vector<double>((n_vert_vars+n_integer_vars)*2);
  856. if (DEBUGPRINT)
  857. printf("\n ALLOCATED X \n");
  858. ///variables part
  859. int ScalarSize = n_vert_vars*2;
  860. int SizeMatrix = (n_vert_vars+n_integer_vars)*2;
  861. ///matrix A
  862. gmm::col_matrix< gmm::wsvector< double > > A(SizeMatrix,SizeMatrix); // lhs matrix variables
  863. ///constraints part
  864. int CsizeX = num_constraint_equations;
  865. int CsizeY = SizeMatrix+1;
  866. gmm::row_matrix< gmm::wsvector< double > > C(CsizeX,CsizeY); // constraints
  867. if (DEBUGPRINT)
  868. printf("\n ALLOCATED QMM STRUCTURES \n");
  869. std::vector<double> B(SizeMatrix,0); // rhs
  870. if (DEBUGPRINT)
  871. printf("\n ALLOCATED RHS STRUCTURES \n");
  872. //// copy LHS
  873. for (int k=0; k < Lhs.outerSize(); ++k){
  874. for (Eigen::SparseMatrix<double>::InnerIterator it(Lhs,k); it; ++it){
  875. int row = it.row();
  876. int col = it.col();
  877. A(row, col) += it.value();
  878. }
  879. }
  880. //// copy Constraints
  881. for (int k=0; k < Constraints.outerSize(); ++k){
  882. for (Eigen::SparseMatrix<double>::InnerIterator it(Constraints,k); it; ++it){
  883. int row = it.row();
  884. int col = it.col();
  885. C(row, col) += it.value();
  886. }
  887. }
  888. if (DEBUGPRINT)
  889. printf("\n SET %d INTEGER VALUES \n",n_integer_vars);
  890. ///add penalization term for integer variables
  891. double penalization = 0.000001;
  892. int offline_index = ScalarSize;
  893. for(unsigned int i = 0; i < (n_integer_vars)*2; ++i)
  894. {
  895. int index=offline_index+i;
  896. A(index,index)=penalization;
  897. }
  898. if (DEBUGPRINT)
  899. printf("\n SET RHS \n");
  900. // copy RHS
  901. for(int i = 0; i < (int)ScalarSize; ++i)
  902. {
  903. B[i] = rhs[i] * cone_grid_res;
  904. }
  905. // copy constraint RHS
  906. if (DEBUGPRINT)
  907. printf("\n SET %d CONSTRAINTS \n",num_constraint_equations);
  908. for(unsigned int i = 0; i < num_constraint_equations; ++i)
  909. {
  910. C(i, SizeMatrix) = -constraints_rhs[i] * cone_grid_res;
  911. }
  912. COMISO::ConstrainedSolver solver;
  913. solver.misolver().set_local_iters(localIter);
  914. solver.misolver().set_direct_rounding(direct_round);
  915. std::sort(ids_to_round.begin(),ids_to_round.end());
  916. std::vector<int>::iterator new_end=std::unique(ids_to_round.begin(),ids_to_round.end());
  917. int dist=distance(ids_to_round.begin(),new_end);
  918. ids_to_round.resize(dist);
  919. solver.solve( C, A, X, B, ids_to_round, 0.0, false, false);
  920. }
  921. template <typename DerivedV, typename DerivedF>
  922. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::clearUserConstraint()
  923. {
  924. num_userdefined_constraint = 0;
  925. userdefined_constraints.clear();
  926. }
  927. template <typename DerivedV, typename DerivedF>
  928. IGL_INLINE void igl::copyleft::comiso::PoissonSolver<DerivedV, DerivedF>::addSharpEdgeConstraint(int fid, int vid)
  929. {
  930. // prepare constraint
  931. std::vector<int> c(Handle_SystemInfo.num_vert_variables*2 + 1);
  932. for (size_t i = 0; i < c.size(); ++i)
  933. {
  934. c[i] = 0;
  935. }
  936. int v1 = Fcut(fid,vid);
  937. int v2 = Fcut(fid,(vid+1)%3);
  938. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> e = Vcut.row(v2) - Vcut.row(v1);
  939. e = e.normalized();
  940. double d1 = fabs(e.dot(PD1.row(fid).normalized()));
  941. double d2 = fabs(e.dot(PD2.row(fid).normalized()));
  942. int offset = 0;
  943. if (d1>d2)
  944. offset = 1;
  945. ids_to_round.push_back((v1 * 2) + offset);
  946. ids_to_round.push_back((v2 * 2) + offset);
  947. // add constraint
  948. c[(v1 * 2) + offset] = 1;
  949. c[(v2 * 2) + offset] = -1;
  950. // add to the user-defined constraints
  951. num_userdefined_constraint++;
  952. userdefined_constraints.push_back(c);
  953. }
  954. template <typename DerivedV, typename DerivedF, typename DerivedU>
  955. IGL_INLINE igl::copyleft::comiso::MIQ_class<DerivedV, DerivedF, DerivedU>::MIQ_class(const Eigen::PlainObjectBase<DerivedV> &V_,
  956. const Eigen::PlainObjectBase<DerivedF> &F_,
  957. const Eigen::PlainObjectBase<DerivedV> &PD1_combed,
  958. const Eigen::PlainObjectBase<DerivedV> &PD2_combed,
  959. const Eigen::Matrix<int, Eigen::Dynamic, 3> &Handle_MMatch,
  960. const Eigen::Matrix<int, Eigen::Dynamic, 1> &Handle_Singular,
  961. const Eigen::Matrix<int, Eigen::Dynamic, 3> &Handle_Seams,
  962. Eigen::PlainObjectBase<DerivedU> &UV,
  963. Eigen::PlainObjectBase<DerivedF> &FUV,
  964. double GradientSize,
  965. double Stiffness,
  966. bool DirectRound,
  967. int iter,
  968. int localIter,
  969. bool DoRound,
  970. bool SingularityRound,
  971. std::vector<int> roundVertices,
  972. std::vector<std::vector<int> > hardFeatures):
  973. V(V_),
  974. F(F_)
  975. {
  976. igl::cut_mesh(V, F, Handle_Seams, Vcut, Fcut);
  977. igl::local_basis(V,F,B1,B2,B3);
  978. igl::triangle_triangle_adjacency(F,TT,TTi);
  979. // Prepare indexing for the linear system
  980. VertexIndexing<DerivedV, DerivedF> VInd(V, F, Vcut, Fcut, TT, TTi, Handle_MMatch, Handle_Singular, Handle_Seams);
  981. VInd.InitSeamInfo();
  982. // Assemble the system and solve
  983. PoissonSolver<DerivedV, DerivedF> PSolver(V,
  984. F,
  985. Vcut,
  986. Fcut,
  987. TT,
  988. TTi,
  989. PD1_combed,
  990. PD2_combed,
  991. Handle_Singular,
  992. VInd.Handle_SystemInfo);
  993. Handle_Stiffness = Eigen::VectorXd::Constant(F.rows(),1);
  994. if (iter > 0) // do stiffening
  995. {
  996. for (int i=0;i<iter;i++)
  997. {
  998. PSolver.SolvePoisson(Handle_Stiffness, GradientSize,1.f,DirectRound,localIter,DoRound,SingularityRound,roundVertices,hardFeatures);
  999. int nflips=NumFlips(PSolver.WUV);
  1000. bool folded = updateStiffeningJacobianDistorsion(GradientSize,PSolver.WUV);
  1001. printf("ITERATION %d FLIPS %d \n",i,nflips);
  1002. if (!folded)break;
  1003. }
  1004. }
  1005. else
  1006. {
  1007. PSolver.SolvePoisson(Handle_Stiffness,GradientSize,1.f,DirectRound,localIter,DoRound,SingularityRound,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> &Handle_MMatch,
  1174. const Eigen::Matrix<int, Eigen::Dynamic, 1> &Handle_Singular,
  1175. const Eigen::Matrix<int, Eigen::Dynamic, 3> &Handle_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. Handle_MMatch,
  1194. Handle_Singular,
  1195. Handle_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. Eigen::PlainObjectBase<DerivedV> BIS1, BIS2;
  1228. igl::compute_frame_field_bisectors(V, F, PD1, PD2, BIS1, BIS2);
  1229. Eigen::PlainObjectBase<DerivedV> BIS1_combed, BIS2_combed;
  1230. igl::comb_cross_field(V, F, BIS1, BIS2, BIS1_combed, BIS2_combed);
  1231. Eigen::PlainObjectBase<DerivedF> Handle_MMatch;
  1232. igl::cross_field_missmatch(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. Eigen::PlainObjectBase<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 specialization
  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