miq.cpp 59 KB

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