miq.cpp 57 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210121112121213121412151216121712181219122012211222122312241225122612271228122912301231123212331234123512361237123812391240124112421243124412451246124712481249125012511252125312541255125612571258125912601261126212631264126512661267126812691270127112721273127412751276127712781279128012811282128312841285128612871288128912901291129212931294129512961297129812991300130113021303130413051306130713081309131013111312131313141315131613171318131913201321132213231324132513261327132813291330133113321333133413351336133713381339134013411342134313441345134613471348134913501351135213531354135513561357135813591360136113621363136413651366136713681369137013711372137313741375137613771378137913801381138213831384138513861387138813891390139113921393139413951396139713981399140014011402140314041405140614071408140914101411141214131414141514161417141814191420142114221423142414251426142714281429143014311432143314341435143614371438143914401441144214431444144514461447144814491450145114521453145414551456145714581459146014611462146314641465146614671468146914701471147214731474147514761477147814791480148114821483148414851486148714881489149014911492149314941495149614971498149915001501150215031504150515061507150815091510151115121513151415151516151715181519152015211522152315241525152615271528152915301531153215331534153515361537153815391540154115421543154415451546154715481549155015511552155315541555155615571558155915601561156215631564156515661567156815691570157115721573157415751576157715781579158015811582158315841585158615871588158915901591159215931594159515961597159815991600160116021603160416051606160716081609161016111612161316141615161616171618161916201621162216231624162516261627162816291630163116321633163416351636163716381639164016411642164316441645164616471648164916501651165216531654165516561657165816591660166116621663166416651666166716681669
  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. std::vector<std::vector<int> >lEdgeSeamInfo; //tmp
  394. // for every vertex, keep track of their adjacent vertices on seams.
  395. std::vector<std::list<int> > VVSeam(V.rows());
  396. Eigen::MatrixXi EV, FE, EF;
  397. igl::edge_topology(V, F, EV, FE, EF);
  398. for (unsigned int e=0;e<EF.rows();e++)
  399. {
  400. int f0 = EF(e,0);
  401. int f1 = EF(e,1);
  402. if (f1 == -1)
  403. continue;
  404. int k=0;
  405. while(k<3)
  406. {
  407. if(FE(f0,k) == e)
  408. break;
  409. k++;
  410. }
  411. bool seam = Handle_Seams(f0,k);
  412. if (seam)
  413. {
  414. int v0 = F(f0, k);
  415. int v1 = F(f0, (k+1)%3);
  416. VVSeam[v0].push_back(v1);
  417. VVSeam[v1].push_back(v0);
  418. }
  419. }
  420. // Find start vertices
  421. std::vector<int> startVertices;
  422. std::vector<bool> isStartVertex(V.rows());
  423. for (unsigned int i=0;i<V.rows();i++)
  424. {
  425. isStartVertex[i] = false;
  426. if (VVSeam[i].size() > 0 && VVSeam[i].size() != 2)
  427. {
  428. startVertices.push_back(i);
  429. isStartVertex[i] = true;
  430. }
  431. }
  432. for (unsigned int i=0;i<startVertices.size();i++)
  433. {
  434. auto startVertex = &VVSeam[startVertices[i]];
  435. for (unsigned int j=0;j<startVertex->size();j++)
  436. {
  437. auto currentVertex = startVertex;
  438. int currentVertexIndex = startVertices[i];
  439. std::vector<int> thisSeam;
  440. thisSeam.push_back(currentVertexIndex);
  441. // walk along the seam
  442. int nextVertexIndex = currentVertex->front();
  443. currentVertex->pop_front();
  444. int prevVertexIndex;
  445. while (true)
  446. {
  447. // update indices (move to the next vertex)
  448. prevVertexIndex = currentVertexIndex;
  449. currentVertexIndex = nextVertexIndex;
  450. currentVertex = &VVSeam[nextVertexIndex];
  451. // add current vertex to this seam
  452. thisSeam.push_back(currentVertexIndex);
  453. // remove the previous vertex
  454. auto it = std::find(currentVertex->begin(), currentVertex->end(), prevVertexIndex);
  455. currentVertex->erase(it);
  456. if (currentVertex->size() == 1 && !isStartVertex[currentVertexIndex])
  457. {
  458. nextVertexIndex = currentVertex->front();
  459. currentVertex->pop_front();
  460. }
  461. else
  462. break;
  463. }
  464. lEdgeSeamInfo.push_back(thisSeam);
  465. }
  466. }
  467. for(auto elem : lEdgeSeamInfo){
  468. for(auto elem2 : elem){
  469. std::cout << elem2 << "\t";
  470. }
  471. std::cout << std::endl;
  472. }
  473. Handle_SystemInfo.EdgeSeamInfo.clear();
  474. std::set<int> hasConstraint;
  475. int integerVar = 0;
  476. for (unsigned int f0=0;f0<F.rows();f0++)
  477. {
  478. for (int k=0;k<3;k++)
  479. {
  480. int f1 = TT(f0,k);
  481. if (f1 == -1)
  482. continue;
  483. bool seam = Handle_Seams(f0,k);
  484. auto search = hasConstraint.find(3*f0 + k);
  485. if (seam && search == hasConstraint.end())
  486. {
  487. int v0,v0p,v1,v1p;
  488. unsigned char MM;
  489. GetSeamInfo(f0,f1,k,v0,v1,v0p,v1p,MM);
  490. Handle_SystemInfo.EdgeSeamInfo.push_back(SeamInfo(v0,v1,v0p,v1p,MM,integerVar));
  491. // mark this face-edge pair and face-edge pair across seam as constrained
  492. hasConstraint.insert(3*f0 + k);
  493. hasConstraint.insert(3*f1 + TTi(f0,k));
  494. integerVar++;
  495. }
  496. }
  497. }
  498. assert(integerVar == Handle_SystemInfo.num_integer_cuts);
  499. }
  500. template <typename DerivedV, typename DerivedF>
  501. IGL_INLINE void igl::comiso::PoissonSolver<DerivedV, DerivedF>::SolvePoisson(Eigen::VectorXd Stiffness,
  502. double vector_field_scale,
  503. double grid_res,
  504. bool direct_round,
  505. int localIter,
  506. bool _integer_rounding,
  507. bool _singularity_rounding,
  508. std::vector<int> roundVertices,
  509. std::vector<std::vector<int> > hardFeatures)
  510. {
  511. Handle_Stiffness = Stiffness;
  512. //initialization of flags and data structures
  513. integer_rounding=_integer_rounding;
  514. ids_to_round.clear();
  515. clearUserConstraint();
  516. // copy the user constraints number
  517. for (size_t i = 0; i < hardFeatures.size(); ++i)
  518. {
  519. addSharpEdgeConstraint(hardFeatures[i][0],hardFeatures[i][1]);
  520. }
  521. ///Initializing Matrix
  522. int t0=clock();
  523. ///initialize the matrix ALLOCATING SPACE
  524. InitMatrix();
  525. if (DEBUGPRINT)
  526. printf("\n ALLOCATED THE MATRIX \n");
  527. ///build the laplacian system
  528. BuildLaplacianMatrix(vector_field_scale);
  529. // add seam constraints
  530. BuildSeamConstraintsExplicitTranslation();
  531. // add user defined constraints
  532. BuildUserDefinedConstraints();
  533. ////add the lagrange multiplier
  534. FixBlockedVertex();
  535. if (DEBUGPRINT)
  536. printf("\n BUILT THE MATRIX \n");
  537. if (integer_rounding)
  538. AddToRoundVertices(roundVertices);
  539. if (_singularity_rounding)
  540. AddSingularityRound();
  541. int t1=clock();
  542. if (DEBUGPRINT) printf("\n time:%d \n",t1-t0);
  543. if (DEBUGPRINT) printf("\n SOLVING \n");
  544. MixedIntegerSolve(grid_res,direct_round,localIter);
  545. int t2=clock();
  546. if (DEBUGPRINT) printf("\n time:%d \n",t2-t1);
  547. if (DEBUGPRINT) printf("\n ASSIGNING COORDS \n");
  548. MapCoords();
  549. int t3=clock();
  550. if (DEBUGPRINT) printf("\n time:%d \n",t3-t2);
  551. if (DEBUGPRINT) printf("\n FINISHED \n");
  552. }
  553. template <typename DerivedV, typename DerivedF>
  554. IGL_INLINE igl::comiso::PoissonSolver<DerivedV, DerivedF>
  555. ::PoissonSolver(const Eigen::PlainObjectBase<DerivedV> &_V,
  556. const Eigen::PlainObjectBase<DerivedF> &_F,
  557. const Eigen::PlainObjectBase<DerivedV> &_Vcut,
  558. const Eigen::PlainObjectBase<DerivedF> &_Fcut,
  559. const Eigen::PlainObjectBase<DerivedF> &_TT,
  560. const Eigen::PlainObjectBase<DerivedF> &_TTi,
  561. const Eigen::PlainObjectBase<DerivedV> &_PD1,
  562. const Eigen::PlainObjectBase<DerivedV> &_PD2,
  563. const Eigen::Matrix<int, Eigen::Dynamic, 1>&_Handle_Singular,
  564. const MeshSystemInfo &_Handle_SystemInfo //todo: const?
  565. ):
  566. V(_V),
  567. F(_F),
  568. Vcut(_Vcut),
  569. Fcut(_Fcut),
  570. TT(_TT),
  571. TTi(_TTi),
  572. PD1(_PD1),
  573. PD2(_PD2),
  574. Handle_Singular(_Handle_Singular),
  575. Handle_SystemInfo(_Handle_SystemInfo)
  576. {
  577. UV = Eigen::MatrixXd(V.rows(),2);
  578. WUV = Eigen::MatrixXd(F.rows(),6);
  579. UV_out = Eigen::MatrixXd(Vcut.rows(),2);
  580. igl::doublearea(V,F,doublearea);
  581. igl::per_face_normals(V,F,N);
  582. igl::vertex_triangle_adjacency(V,F,VF,VFi);
  583. }
  584. ///START COMMON MATH FUNCTIONS
  585. ///return the complex encoding the rotation
  586. ///for a given missmatch interval
  587. template <typename DerivedV, typename DerivedF>
  588. IGL_INLINE std::complex<double> igl::comiso::PoissonSolver<DerivedV, DerivedF>::GetRotationComplex(int interval)
  589. {
  590. assert((interval>=0)&&(interval<4));
  591. switch(interval)
  592. {
  593. case 0:return std::complex<double>(1,0);
  594. case 1:return std::complex<double>(0,1);
  595. case 2:return std::complex<double>(-1,0);
  596. default:return std::complex<double>(0,-1);
  597. }
  598. }
  599. ///END COMMON MATH FUNCTIONS
  600. ///START FIXING VERTICES
  601. ///set a given vertex as fixed
  602. template <typename DerivedV, typename DerivedF>
  603. IGL_INLINE void igl::comiso::PoissonSolver<DerivedV, DerivedF>::AddFixedVertex(int v)
  604. {
  605. n_fixed_vars++;
  606. Hard_constraints.push_back(v);
  607. }
  608. ///find vertex to fix in case we're using
  609. ///a vector field NB: multiple components not handled
  610. template <typename DerivedV, typename DerivedF>
  611. IGL_INLINE void igl::comiso::PoissonSolver<DerivedV, DerivedF>::FindFixedVertField()
  612. {
  613. Hard_constraints.clear();
  614. n_fixed_vars=0;
  615. //fix the first singularity
  616. for (unsigned int v=0;v<V.rows();v++)
  617. {
  618. if (Handle_Singular(v))
  619. {
  620. AddFixedVertex(v);
  621. UV.row(v) << 0,0;
  622. return;
  623. }
  624. }
  625. ///if anything fixed fix the first
  626. AddFixedVertex(0); // TODO HERE IT ISSSSSS
  627. UV.row(0) << 0,0;
  628. std::cerr << "No vertices to fix, I am fixing the first vertex to the origin!" << std::endl;
  629. }
  630. ///find hard constraint depending if using or not
  631. ///a vector field
  632. template <typename DerivedV, typename DerivedF>
  633. IGL_INLINE void igl::comiso::PoissonSolver<DerivedV, DerivedF>::FindFixedVert()
  634. {
  635. Hard_constraints.clear();
  636. FindFixedVertField();
  637. }
  638. template <typename DerivedV, typename DerivedF>
  639. IGL_INLINE int igl::comiso::PoissonSolver<DerivedV, DerivedF>::GetFirstVertexIndex(int v)
  640. {
  641. return Fcut(VF[v][0],VFi[v][0]);
  642. }
  643. ///fix the vertices which are flagged as fixed
  644. template <typename DerivedV, typename DerivedF>
  645. IGL_INLINE void igl::comiso::PoissonSolver<DerivedV, DerivedF>::FixBlockedVertex()
  646. {
  647. int offset_row = num_cut_constraint*2;
  648. unsigned int constr_num = 0;
  649. for (unsigned int i=0;i<Hard_constraints.size();i++)
  650. {
  651. int v = Hard_constraints[i];
  652. ///get first index of the vertex that must blocked
  653. //int index=v->vertex_index[0];
  654. int index = GetFirstVertexIndex(v);
  655. ///multiply times 2 because of uv
  656. int indexvert = index*2;
  657. ///find the first free row to add the constraint
  658. int indexRow = (offset_row+constr_num*2);
  659. int indexCol = indexRow;
  660. ///add fixing constraint LHS
  661. Constraints.coeffRef(indexRow, indexvert) += 1;
  662. Constraints.coeffRef(indexRow+1,indexvert+1) += 1;
  663. ///add fixing constraint RHS
  664. constraints_rhs[indexCol] = UV(v,0);
  665. constraints_rhs[indexCol+1] = UV(v,1);
  666. constr_num++;
  667. }
  668. assert(constr_num==n_fixed_vars);
  669. }
  670. ///END FIXING VERTICES
  671. ///HANDLING SINGULARITY
  672. //set the singularity round to integer location
  673. template <typename DerivedV, typename DerivedF>
  674. IGL_INLINE void igl::comiso::PoissonSolver<DerivedV, DerivedF>::AddSingularityRound()
  675. {
  676. for (unsigned int v=0;v<V.rows();v++)
  677. {
  678. if (Handle_Singular(v))
  679. {
  680. int index0=GetFirstVertexIndex(v);
  681. ids_to_round.push_back( index0*2 );
  682. ids_to_round.push_back((index0*2)+1);
  683. }
  684. }
  685. }
  686. template <typename DerivedV, typename DerivedF>
  687. IGL_INLINE void igl::comiso::PoissonSolver<DerivedV, DerivedF>::AddToRoundVertices(std::vector<int> ids)
  688. {
  689. for (size_t i = 0; i < ids.size(); ++i)
  690. {
  691. if (ids[i] < 0 || ids[i] >= V.rows())
  692. std::cerr << "WARNING: Ignored round vertex constraint, vertex " << ids[i] << " does not exist in the mesh." << std::endl;
  693. int index0 = GetFirstVertexIndex(ids[i]);
  694. ids_to_round.push_back( index0*2 );
  695. ids_to_round.push_back((index0*2)+1);
  696. }
  697. }
  698. ///START GENERIC SYSTEM FUNCTIONS
  699. //build the laplacian matrix cyclyng over all rangemaps
  700. //and over all faces
  701. template <typename DerivedV, typename DerivedF>
  702. IGL_INLINE void igl::comiso::PoissonSolver<DerivedV, DerivedF>::BuildLaplacianMatrix(double vfscale)
  703. {
  704. Eigen::VectorXi idx = Eigen::VectorXi::LinSpaced(Vcut.rows(), 0, 2*Vcut.rows()-2);
  705. Eigen::VectorXi idx2 = Eigen::VectorXi::LinSpaced(Vcut.rows(), 1, 2*Vcut.rows()-1);
  706. // get gradient matrix
  707. Eigen::SparseMatrix<double> G(Fcut.rows() * 3, Vcut.rows());
  708. igl::grad(Vcut, Fcut, G);
  709. // get triangle weights
  710. Eigen::VectorXd dblA(Fcut.rows());
  711. igl::doublearea(Vcut, Fcut, dblA);
  712. // compute intermediate result
  713. Eigen::SparseMatrix<double> G2;
  714. G2 = G.transpose() * dblA.replicate<3,1>().asDiagonal() * Handle_Stiffness.replicate<3,1>().asDiagonal();
  715. /// Compute LHS
  716. Eigen::SparseMatrix<double> Cotmatrix;
  717. Cotmatrix = 0.5 * G2 * G;
  718. igl::slice_into(Cotmatrix, idx, idx, Lhs);
  719. igl::slice_into(Cotmatrix, idx2, idx2, Lhs);
  720. /// Compute RHS
  721. // reshape nrosy vectors
  722. 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.
  723. 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.
  724. // multiply with weights
  725. Eigen::VectorXd rhs1 = G2 * u * 0.5 * vfscale;
  726. Eigen::VectorXd rhs2 = -G2 * v * 0.5 * vfscale;
  727. igl::slice_into(rhs1, idx, 1, rhs);
  728. igl::slice_into(rhs2, idx2, 1, rhs);
  729. }
  730. ///find different sized of the system
  731. template <typename DerivedV, typename DerivedF>
  732. IGL_INLINE void igl::comiso::PoissonSolver<DerivedV, DerivedF>::FindSizes()
  733. {
  734. ///find the vertex that need to be fixed
  735. FindFixedVert();
  736. ///REAL PART
  737. n_vert_vars = Handle_SystemInfo.num_vert_variables;
  738. ///INTEGER PART
  739. ///the total number of integer variables
  740. n_integer_vars = Handle_SystemInfo.num_integer_cuts;
  741. ///CONSTRAINT PART
  742. num_cut_constraint = Handle_SystemInfo.EdgeSeamInfo.size()*2;
  743. num_constraint_equations = num_cut_constraint * 2 + n_fixed_vars * 2 + num_userdefined_constraint;
  744. ///total variable of the system
  745. num_total_vars = (n_vert_vars+n_integer_vars) * 2;
  746. ///initialize matrix size
  747. system_size = num_total_vars + num_constraint_equations;
  748. if (DEBUGPRINT) printf("\n*** SYSTEM VARIABLES *** \n");
  749. if (DEBUGPRINT) printf("* NUM REAL VERTEX VARIABLES %d \n",n_vert_vars);
  750. if (DEBUGPRINT) printf("\n*** SINGULARITY *** \n ");
  751. if (DEBUGPRINT) printf("* NUM SINGULARITY %d\n",(int)ids_to_round.size()/2);
  752. if (DEBUGPRINT) printf("\n*** INTEGER VARIABLES *** \n");
  753. if (DEBUGPRINT) printf("* NUM INTEGER VARIABLES %d \n",(int)n_integer_vars);
  754. if (DEBUGPRINT) printf("\n*** CONSTRAINTS *** \n ");
  755. if (DEBUGPRINT) printf("* NUM FIXED CONSTRAINTS %d\n",n_fixed_vars);
  756. if (DEBUGPRINT) printf("* NUM CUTS CONSTRAINTS %d\n",num_cut_constraint);
  757. if (DEBUGPRINT) printf("* NUM USER DEFINED CONSTRAINTS %d\n",num_userdefined_constraint);
  758. if (DEBUGPRINT) printf("\n*** TOTAL SIZE *** \n");
  759. if (DEBUGPRINT) printf("* TOTAL VARIABLE SIZE (WITH INTEGER TRASL) %d \n",num_total_vars);
  760. if (DEBUGPRINT) printf("* TOTAL CONSTRAINTS %d \n",num_constraint_equations);
  761. if (DEBUGPRINT) printf("* MATRIX SIZE %d \n",system_size);
  762. }
  763. template <typename DerivedV, typename DerivedF>
  764. IGL_INLINE void igl::comiso::PoissonSolver<DerivedV, DerivedF>::AllocateSystem()
  765. {
  766. Lhs.resize(n_vert_vars * 2, n_vert_vars * 2);
  767. Constraints.resize(num_constraint_equations, system_size);
  768. rhs.resize(system_size);
  769. constraints_rhs.resize(num_constraint_equations);
  770. printf("\n INITIALIZED SPARSE MATRIX OF %d x %d \n",system_size, system_size);
  771. printf("\n INITIALIZED SPARSE MATRIX OF %d x %d \n",num_constraint_equations, system_size);
  772. printf("\n INITIALIZED VECTOR OF %d x 1 \n",system_size);
  773. printf("\n INITIALIZED VECTOR OF %d x 1 \n",num_constraint_equations);
  774. }
  775. ///intitialize the whole matrix
  776. template <typename DerivedV, typename DerivedF>
  777. IGL_INLINE void igl::comiso::PoissonSolver<DerivedV, DerivedF>::InitMatrix()
  778. {
  779. FindSizes();
  780. AllocateSystem();
  781. }
  782. ///map back coordinates after that
  783. ///the system has been solved
  784. template <typename DerivedV, typename DerivedF>
  785. IGL_INLINE void igl::comiso::PoissonSolver<DerivedV, DerivedF>::MapCoords()
  786. {
  787. ///map coords to faces
  788. for (unsigned int f=0;f<Fcut.rows();f++)
  789. {
  790. for (int k=0;k<3;k++)
  791. {
  792. //get the index of the variable in the system
  793. int indexUV = Fcut(f,k);
  794. ///then get U and V coords
  795. double U=X[indexUV*2];
  796. double V=X[indexUV*2+1];
  797. WUV(f,k*2 + 0) = U;
  798. WUV(f,k*2 + 1) = V;
  799. }
  800. }
  801. for(int i = 0; i < Vcut.rows(); i++){
  802. UV_out(i,0) = X[i*2];
  803. UV_out(i,1) = X[i*2+1];
  804. }
  805. #if 0
  806. ///initialize the vector of integer variables to return their values
  807. Handle_SystemInfo.IntegerValues.resize(n_integer_vars*2);
  808. int baseIndex = (n_vert_vars)*2;
  809. int endIndex = baseIndex+n_integer_vars*2;
  810. int index = 0;
  811. for (int i=baseIndex; i<endIndex; i++)
  812. {
  813. ///assert that the value is an integer value
  814. double value=X[i];
  815. double diff = value-(int)floor(value+0.5);
  816. assert(diff<0.00000001);
  817. Handle_SystemInfo.IntegerValues[index] = value;
  818. index++;
  819. }
  820. #endif
  821. }
  822. ///END GENERIC SYSTEM FUNCTIONS
  823. ///set the constraints for the inter-range cuts
  824. template <typename DerivedV, typename DerivedF>
  825. IGL_INLINE void igl::comiso::PoissonSolver<DerivedV, DerivedF>::BuildSeamConstraintsExplicitTranslation()
  826. {
  827. ///current constraint row
  828. int constr_row = 0;
  829. for (unsigned int i=0; i<num_cut_constraint/2; i++)
  830. {
  831. unsigned char interval = Handle_SystemInfo.EdgeSeamInfo[i].MMatch;
  832. if (interval==1)
  833. interval=3;
  834. else
  835. if(interval==3)
  836. interval=1;
  837. int p0 = Handle_SystemInfo.EdgeSeamInfo[i].v0;
  838. int p1 = Handle_SystemInfo.EdgeSeamInfo[i].v1;
  839. int p0p = Handle_SystemInfo.EdgeSeamInfo[i].v0p;
  840. int p1p = Handle_SystemInfo.EdgeSeamInfo[i].v1p;
  841. std::complex<double> rot = GetRotationComplex(interval);
  842. ///get the integer variable
  843. int integerVar = n_vert_vars + Handle_SystemInfo.EdgeSeamInfo[i].integerVar;
  844. if (integer_rounding)
  845. {
  846. ids_to_round.push_back(integerVar*2);
  847. ids_to_round.push_back(integerVar*2+1);
  848. }
  849. // TODO: exploit fact that rotations have either zeros on diagonal (real) or off-diagonal (imag). don't explicitly store the zeros.
  850. // cross boundary compatibility conditions
  851. // constraints for start vertex of edge
  852. Constraints.coeffRef(constr_row, 2*p0) += rot.real();
  853. Constraints.coeffRef(constr_row, 2*p0+1) += -rot.imag();
  854. Constraints.coeffRef(constr_row+1, 2*p0) += rot.imag();
  855. Constraints.coeffRef(constr_row+1, 2*p0+1) += rot.real();
  856. Constraints.coeffRef(constr_row, 2*p0p) += -1;
  857. Constraints.coeffRef(constr_row+1, 2*p0p+1) += -1;
  858. Constraints.coeffRef(constr_row, 2*integerVar) += 1;
  859. Constraints.coeffRef(constr_row+1, 2*integerVar+1) += 1;
  860. constraints_rhs[constr_row] = 0;
  861. constraints_rhs[constr_row+1] = 0;
  862. constr_row += 2;
  863. // constraints for end vertex of edge
  864. Constraints.coeffRef(constr_row, 2*p1) += rot.real();
  865. Constraints.coeffRef(constr_row, 2*p1+1) += -rot.imag();
  866. Constraints.coeffRef(constr_row+1, 2*p1) += rot.imag();
  867. Constraints.coeffRef(constr_row+1, 2*p1+1) += rot.real();
  868. Constraints.coeffRef(constr_row, 2*p1p) += -1;
  869. Constraints.coeffRef(constr_row+1, 2*p1p+1) += -1;
  870. Constraints.coeffRef(constr_row, 2*integerVar) += 1;
  871. Constraints.coeffRef(constr_row+1, 2*integerVar+1) += 1;
  872. constraints_rhs[constr_row] = 0;
  873. constraints_rhs[constr_row+1] = 0;
  874. constr_row += 2;
  875. }
  876. }
  877. ///set the constraints for the inter-range cuts
  878. template <typename DerivedV, typename DerivedF>
  879. IGL_INLINE void igl::comiso::PoissonSolver<DerivedV, DerivedF>::BuildUserDefinedConstraints()
  880. {
  881. /// the user defined constraints are at the end
  882. int offset_row = num_cut_constraint*2 + n_fixed_vars*2;
  883. ///current constraint row
  884. int constr_row = offset_row;
  885. assert(num_userdefined_constraint == userdefined_constraints.size());
  886. for (unsigned int i=0; i<num_userdefined_constraint; i++)
  887. {
  888. for (unsigned int j=0; j<userdefined_constraints[i].size()-1; ++j)
  889. {
  890. Constraints.coeffRef(constr_row, j) = userdefined_constraints[i][j];
  891. }
  892. constraints_rhs[constr_row] = userdefined_constraints[i][userdefined_constraints[i].size()-1];
  893. constr_row +=1;
  894. }
  895. }
  896. ///call of the mixed integer solver
  897. template <typename DerivedV, typename DerivedF>
  898. IGL_INLINE void igl::comiso::PoissonSolver<DerivedV, DerivedF>::MixedIntegerSolve(double cone_grid_res,
  899. bool direct_round,
  900. int localIter)
  901. {
  902. X = std::vector<double>((n_vert_vars+n_integer_vars)*2);
  903. ///variables part
  904. int ScalarSize = n_vert_vars*2;
  905. int SizeMatrix = (n_vert_vars+n_integer_vars)*2;
  906. if (DEBUGPRINT)
  907. printf("\n ALLOCATED X \n");
  908. ///matrix A
  909. gmm::col_matrix< gmm::wsvector< double > > A(SizeMatrix,SizeMatrix); // lhs matrix variables +
  910. ///constraints part
  911. int CsizeX = num_constraint_equations;
  912. int CsizeY = SizeMatrix+1;
  913. gmm::row_matrix< gmm::wsvector< double > > C(CsizeX,CsizeY); // constraints
  914. if (DEBUGPRINT)
  915. printf("\n ALLOCATED QMM STRUCTURES \n");
  916. std::vector<double> B(SizeMatrix,0); // rhs
  917. if (DEBUGPRINT)
  918. printf("\n ALLOCATED RHS STRUCTURES \n");
  919. //// copy LHS
  920. for (int k=0; k < Lhs.outerSize(); ++k){
  921. for (Eigen::SparseMatrix<double>::InnerIterator it(Lhs,k); it; ++it){
  922. int row = it.row();
  923. int col = it.col();
  924. A(row, col) += it.value();
  925. }
  926. }
  927. //// copy Constraints
  928. for (int k=0; k < Constraints.outerSize(); ++k){
  929. for (Eigen::SparseMatrix<double>::InnerIterator it(Constraints,k); it; ++it){
  930. int row = it.row();
  931. int col = it.col();
  932. C(row, col) += it.value();
  933. }
  934. }
  935. if (DEBUGPRINT)
  936. printf("\n SET %d INTEGER VALUES \n",n_integer_vars);
  937. ///add penalization term for integer variables
  938. double penalization = 0.000001;
  939. int offline_index = ScalarSize;
  940. for(unsigned int i = 0; i < (n_integer_vars)*2; ++i)
  941. {
  942. int index=offline_index+i;
  943. A(index,index)=penalization;
  944. }
  945. if (DEBUGPRINT)
  946. printf("\n SET RHS \n");
  947. // copy RHS
  948. for(int i = 0; i < (int)ScalarSize; ++i)
  949. {
  950. B[i] = rhs[i] * cone_grid_res;
  951. }
  952. // copy constraint RHS
  953. if (DEBUGPRINT)
  954. printf("\n SET %d CONSTRAINTS \n",num_constraint_equations);
  955. for(unsigned int i = 0; i < num_constraint_equations; ++i)
  956. {
  957. C(i, SizeMatrix) = -constraints_rhs[i] * cone_grid_res;
  958. }
  959. ///copy values back into S
  960. COMISO::ConstrainedSolver solver;
  961. solver.misolver().set_local_iters(localIter);
  962. solver.misolver().set_direct_rounding(direct_round);
  963. std::sort(ids_to_round.begin(),ids_to_round.end());
  964. std::vector<int>::iterator new_end=std::unique(ids_to_round.begin(),ids_to_round.end());
  965. int dist=distance(ids_to_round.begin(),new_end);
  966. ids_to_round.resize(dist);
  967. solver.solve( C, A, X, B, ids_to_round, 0.0, false, false);
  968. ////DEBUG OUTPUT
  969. std::ofstream idsout("ids.txt");
  970. for(auto elem : ids_to_round){
  971. idsout << elem << std::endl;
  972. }
  973. idsout.close();
  974. std::ofstream consout("Cmat.txt");
  975. Eigen::SparseMatrix<double, Eigen::RowMajor> Cmat = Constraints;
  976. for (int k=0; k < Cmat.outerSize(); ++k){
  977. for (Eigen::SparseMatrix<double, Eigen::RowMajor>::InnerIterator it(Cmat,k); it; ++it){
  978. int row = it.row();
  979. int col = it.col();
  980. consout << "(" << row << ", " << col << ")" << "\t" << it.value() << std::endl;
  981. }
  982. }
  983. consout.close();
  984. std::ofstream rhsCout("rhsC.txt");
  985. rhsCout << rhs;
  986. rhsCout.close();
  987. std::ofstream xout("Xout.txt");
  988. for(auto it = X.begin(); it != X.end(); it+=2){
  989. xout << *it << "\t" << *(it+1) << std::endl;
  990. }
  991. xout.close();
  992. }
  993. template <typename DerivedV, typename DerivedF>
  994. IGL_INLINE void igl::comiso::PoissonSolver<DerivedV, DerivedF>::clearUserConstraint()
  995. {
  996. num_userdefined_constraint = 0;
  997. userdefined_constraints.clear();
  998. }
  999. template <typename DerivedV, typename DerivedF>
  1000. IGL_INLINE void igl::comiso::PoissonSolver<DerivedV, DerivedF>::addSharpEdgeConstraint(int fid, int vid)
  1001. {
  1002. // prepare constraint
  1003. std::vector<int> c(Handle_SystemInfo.num_vert_variables*2 + 1);
  1004. for (size_t i = 0; i < c.size(); ++i)
  1005. {
  1006. c[i] = 0;
  1007. }
  1008. int v1 = Fcut(fid,vid);
  1009. int v2 = Fcut(fid,(vid+1)%3);
  1010. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> e = Vcut.row(v2) - Vcut.row(v1);
  1011. e = e.normalized();
  1012. double d1 = fabs(e.dot(PD1.row(fid).normalized()));
  1013. double d2 = fabs(e.dot(PD2.row(fid).normalized()));
  1014. int offset = 0;
  1015. if (d1>d2)
  1016. offset = 1;
  1017. ids_to_round.push_back((v1 * 2) + offset);
  1018. ids_to_round.push_back((v2 * 2) + offset);
  1019. // add constraint
  1020. c[(v1 * 2) + offset] = 1;
  1021. c[(v2 * 2) + offset] = -1;
  1022. // add to the user-defined constraints
  1023. num_userdefined_constraint++;
  1024. userdefined_constraints.push_back(c);
  1025. }
  1026. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1027. IGL_INLINE igl::comiso::MIQ_class<DerivedV, DerivedF, DerivedU>::MIQ_class(const Eigen::PlainObjectBase<DerivedV> &V_,
  1028. const Eigen::PlainObjectBase<DerivedF> &F_,
  1029. const Eigen::PlainObjectBase<DerivedV> &PD1_combed,
  1030. const Eigen::PlainObjectBase<DerivedV> &PD2_combed,
  1031. // const Eigen::PlainObjectBase<DerivedV> &BIS1_combed,
  1032. // const Eigen::PlainObjectBase<DerivedV> &BIS2_combed,
  1033. const Eigen::Matrix<int, Eigen::Dynamic, 3> &Handle_MMatch,
  1034. const Eigen::Matrix<int, Eigen::Dynamic, 1> &Handle_Singular,
  1035. // const Eigen::Matrix<int, Eigen::Dynamic, 1> &Handle_SingularDegree,
  1036. const Eigen::Matrix<int, Eigen::Dynamic, 3> &Handle_Seams,
  1037. Eigen::PlainObjectBase<DerivedU> &UV,
  1038. Eigen::PlainObjectBase<DerivedF> &FUV,
  1039. double GradientSize,
  1040. double Stiffness,
  1041. bool DirectRound,
  1042. int iter,
  1043. int localIter,
  1044. bool DoRound,
  1045. bool SingularityRound,
  1046. std::vector<int> roundVertices,
  1047. std::vector<std::vector<int> > hardFeatures):
  1048. V(V_),
  1049. F(F_)
  1050. {
  1051. igl::cut_mesh(V, F, Handle_Seams, Vcut, Fcut);
  1052. igl::local_basis(V,F,B1,B2,B3);
  1053. igl::triangle_triangle_adjacency(V,F,TT,TTi);
  1054. // Prepare indexing for the linear system
  1055. VertexIndexing<DerivedV, DerivedF> VInd(V, F, Vcut, Fcut, TT, TTi, /*BIS1_combed, BIS2_combed,*/ Handle_MMatch, /*Handle_Singular, Handle_SingularDegree,*/ Handle_Seams);
  1056. VInd.InitFaceIntegerVal();
  1057. VInd.InitSeamInfo();
  1058. // Eigen::PlainObjectBase<DerivedV> PD1_combed_for_poisson, PD2_combed_for_poisson;
  1059. // // Rotate by 90 degrees CCW
  1060. // PD1_combed_for_poisson.setZero(PD1_combed.rows(),3);
  1061. // PD2_combed_for_poisson.setZero(PD2_combed.rows(),3);
  1062. // for (unsigned i=0; i<PD1_combed.rows();++i)
  1063. // {
  1064. // double n1 = PD1_combed.row(i).norm();
  1065. // double n2 = PD2_combed.row(i).norm();
  1066. //
  1067. // double a1 = atan2(B2.row(i).dot(PD1_combed.row(i)),B1.row(i).dot(PD1_combed.row(i)));
  1068. // double a2 = atan2(B2.row(i).dot(PD2_combed.row(i)),B1.row(i).dot(PD2_combed.row(i)));
  1069. //
  1070. // // a1 += M_PI/2;
  1071. // // a2 += M_PI/2;
  1072. //
  1073. //
  1074. // PD1_combed_for_poisson.row(i) = cos(a1) * B1.row(i) + sin(a1) * B2.row(i);
  1075. // PD2_combed_for_poisson.row(i) = cos(a2) * B1.row(i) + sin(a2) * B2.row(i);
  1076. //
  1077. // PD1_combed_for_poisson.row(i) = PD1_combed_for_poisson.row(i).normalized() * n1;
  1078. // PD2_combed_for_poisson.row(i) = PD2_combed_for_poisson.row(i).normalized() * n2;
  1079. // }
  1080. // Assemble the system and solve
  1081. PoissonSolver<DerivedV, DerivedF> PSolver(V,
  1082. F,
  1083. Vcut,
  1084. Fcut,
  1085. TT,
  1086. TTi,
  1087. PD1_combed,
  1088. PD2_combed,
  1089. /*VInd.Handle_Singular*/Handle_Singular,
  1090. VInd.Handle_SystemInfo);
  1091. Handle_Stiffness = Eigen::VectorXd::Constant(F.rows(),1);
  1092. if (iter > 0) // do stiffening
  1093. {
  1094. for (int i=0;i<iter;i++)
  1095. {
  1096. PSolver.SolvePoisson(Handle_Stiffness, GradientSize,1.f,DirectRound,localIter,DoRound,SingularityRound,roundVertices,hardFeatures);
  1097. int nflips=NumFlips(PSolver.WUV);
  1098. bool folded = updateStiffeningJacobianDistorsion(GradientSize,PSolver.WUV);
  1099. printf("ITERATION %d FLIPS %d \n",i,nflips);
  1100. if (!folded)break;
  1101. }
  1102. }
  1103. else
  1104. {
  1105. PSolver.SolvePoisson(Handle_Stiffness,GradientSize,1.f,DirectRound,localIter,DoRound,SingularityRound,roundVertices,hardFeatures);
  1106. }
  1107. int nflips=NumFlips(PSolver.WUV);
  1108. printf("**** END OPTIMIZING #FLIPS %d ****\n",nflips);
  1109. UV_out = PSolver.UV_out;
  1110. FUV_out = PSolver.Fcut;
  1111. fflush(stdout);
  1112. }
  1113. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1114. IGL_INLINE void igl::comiso::MIQ_class<DerivedV, DerivedF, DerivedU>::extractUV(Eigen::PlainObjectBase<DerivedU> &UV_out,
  1115. Eigen::PlainObjectBase<DerivedF> &FUV_out)
  1116. {
  1117. UV_out = this->UV_out;
  1118. FUV_out = this->FUV_out;
  1119. }
  1120. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1121. IGL_INLINE int igl::comiso::MIQ_class<DerivedV, DerivedF, DerivedU>::NumFlips(const Eigen::MatrixXd& WUV)
  1122. {
  1123. int numFl=0;
  1124. for (unsigned int i=0;i<F.rows();i++)
  1125. {
  1126. if (IsFlipped(i, WUV))
  1127. numFl++;
  1128. }
  1129. return numFl;
  1130. }
  1131. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1132. IGL_INLINE double igl::comiso::MIQ_class<DerivedV, DerivedF, DerivedU>::Distortion(int f, double h, const Eigen::MatrixXd& WUV)
  1133. {
  1134. assert(h > 0);
  1135. Eigen::Vector2d uv0,uv1,uv2;
  1136. uv0 << WUV(f,0), WUV(f,1);
  1137. uv1 << WUV(f,2), WUV(f,3);
  1138. uv2 << WUV(f,4), WUV(f,5);
  1139. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> p0 = V.row(Fcut(f,0));
  1140. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> p1 = V.row(Fcut(f,1));
  1141. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> p2 = V.row(Fcut(f,2));
  1142. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> norm = (p1 - p0).cross(p2 - p0);
  1143. double area2 = norm.norm();
  1144. double area2_inv = 1.0 / area2;
  1145. norm *= area2_inv;
  1146. if (area2 > 0)
  1147. {
  1148. // Singular values of the Jacobian
  1149. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> neg_t0 = norm.cross(p2 - p1);
  1150. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> neg_t1 = norm.cross(p0 - p2);
  1151. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> neg_t2 = norm.cross(p1 - p0);
  1152. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> diffu = (neg_t0 * uv0(0) +neg_t1 *uv1(0) + neg_t2 * uv2(0) )*area2_inv;
  1153. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> diffv = (neg_t0 * uv0(1) + neg_t1*uv1(1) + neg_t2*uv2(1) )*area2_inv;
  1154. // first fundamental form
  1155. double I00 = diffu.dot(diffu); // guaranteed non-neg
  1156. double I01 = diffu.dot(diffv); // I01 = I10
  1157. double I11 = diffv.dot(diffv); // guaranteed non-neg
  1158. // eigenvalues of a 2x2 matrix
  1159. // [a00 a01]
  1160. // [a10 a11]
  1161. // 1/2 * [ (a00 + a11) +/- sqrt((a00 - a11)^2 + 4 a01 a10) ]
  1162. double trI = I00 + I11; // guaranteed non-neg
  1163. double diffDiag = I00 - I11; // guaranteed non-neg
  1164. double sqrtDet = sqrt(std::max(0.0, diffDiag*diffDiag +
  1165. 4 * I01 * I01)); // guaranteed non-neg
  1166. double sig1 = 0.5 * (trI + sqrtDet); // higher singular value
  1167. double sig2 = 0.5 * (trI - sqrtDet); // lower singular value
  1168. // Avoid sig2 < 0 due to numerical error
  1169. if (fabs(sig2) < 1.0e-8)
  1170. sig2 = 0;
  1171. assert(sig1 >= 0);
  1172. assert(sig2 >= 0);
  1173. if (sig2 < 0) {
  1174. printf("Distortion will be NaN! sig1^2 is negative (%lg)\n",
  1175. sig2);
  1176. }
  1177. // The singular values of the Jacobian are the sqrts of the
  1178. // eigenvalues of the first fundamental form.
  1179. sig1 = sqrt(sig1);
  1180. sig2 = sqrt(sig2);
  1181. // distortion
  1182. double tao = IsFlipped(f,WUV) ? -1 : 1;
  1183. double factor = tao / h;
  1184. double lam = fabs(factor * sig1 - 1) + fabs(factor * sig2 - 1);
  1185. return lam;
  1186. }
  1187. else {
  1188. return 10; // something "large"
  1189. }
  1190. }
  1191. ////////////////////////////////////////////////////////////////////////////
  1192. // Approximate the distortion laplacian using a uniform laplacian on
  1193. // the dual mesh:
  1194. // ___________
  1195. // \-1 / \-1 /
  1196. // \ / 3 \ /
  1197. // \-----/
  1198. // \-1 /
  1199. // \ /
  1200. //
  1201. // @param[in] f facet on which to compute distortion laplacian
  1202. // @param[in] h scaling factor applied to cross field
  1203. // @return distortion laplacian for f
  1204. ///////////////////////////////////////////////////////////////////////////
  1205. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1206. IGL_INLINE double igl::comiso::MIQ_class<DerivedV, DerivedF, DerivedU>::LaplaceDistortion(const int f, double h, const Eigen::MatrixXd& WUV)
  1207. {
  1208. double mydist = Distortion(f, h, WUV);
  1209. double lapl=0;
  1210. for (int i=0;i<3;i++)
  1211. {
  1212. if (TT(f,i) != -1)
  1213. lapl += (mydist - Distortion(TT(f,i), h, WUV));
  1214. }
  1215. return lapl;
  1216. }
  1217. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1218. IGL_INLINE bool igl::comiso::MIQ_class<DerivedV, DerivedF, DerivedU>::updateStiffeningJacobianDistorsion(double grad_size, const Eigen::MatrixXd& WUV)
  1219. {
  1220. bool flipped = NumFlips(WUV)>0;
  1221. if (!flipped)
  1222. return false;
  1223. double maxL=0;
  1224. double maxD=0;
  1225. if (flipped)
  1226. {
  1227. const double c = 1.0;
  1228. const double d = 5.0;
  1229. for (unsigned int i = 0; i < Fcut.rows(); ++i)
  1230. {
  1231. double dist=Distortion(i,grad_size,WUV);
  1232. if (dist > maxD)
  1233. maxD=dist;
  1234. double absLap=fabs(LaplaceDistortion(i, grad_size,WUV));
  1235. if (absLap > maxL)
  1236. maxL = absLap;
  1237. double stiffDelta = std::min(c * absLap, d);
  1238. Handle_Stiffness[i]+=stiffDelta;
  1239. }
  1240. }
  1241. printf("Maximum Distorsion %4.4f \n",maxD);
  1242. printf("Maximum Laplacian %4.4f \n",maxL);
  1243. return flipped;
  1244. }
  1245. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1246. IGL_INLINE bool igl::comiso::MIQ_class<DerivedV, DerivedF, DerivedU>::IsFlipped(const Eigen::Vector2d &uv0,
  1247. const Eigen::Vector2d &uv1,
  1248. const Eigen::Vector2d &uv2)
  1249. {
  1250. Eigen::Vector2d e0 = (uv1-uv0);
  1251. Eigen::Vector2d e1 = (uv2-uv0);
  1252. double Area = e0(0)*e1(1) - e0(1)*e1(0);
  1253. return (Area<=0);
  1254. }
  1255. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1256. IGL_INLINE bool igl::comiso::MIQ_class<DerivedV, DerivedF, DerivedU>::IsFlipped(
  1257. const int i, const Eigen::MatrixXd& WUV)
  1258. {
  1259. Eigen::Vector2d uv0,uv1,uv2;
  1260. uv0 << WUV(i,0), WUV(i,1);
  1261. uv1 << WUV(i,2), WUV(i,3);
  1262. uv2 << WUV(i,4), WUV(i,5);
  1263. return (IsFlipped(uv0,uv1,uv2));
  1264. }
  1265. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1266. IGL_INLINE void igl::comiso::miq(
  1267. const Eigen::PlainObjectBase<DerivedV> &V,
  1268. const Eigen::PlainObjectBase<DerivedF> &F,
  1269. const Eigen::PlainObjectBase<DerivedV> &PD1_combed,
  1270. const Eigen::PlainObjectBase<DerivedV> &PD2_combed,
  1271. // const Eigen::PlainObjectBase<DerivedV> &BIS1_combed,
  1272. // const Eigen::PlainObjectBase<DerivedV> &BIS2_combed,
  1273. const Eigen::Matrix<int, Eigen::Dynamic, 3> &Handle_MMatch,
  1274. const Eigen::Matrix<int, Eigen::Dynamic, 1> &Handle_Singular,
  1275. // const Eigen::Matrix<int, Eigen::Dynamic, 1> &Handle_SingularDegree,
  1276. const Eigen::Matrix<int, Eigen::Dynamic, 3> &Handle_Seams,
  1277. Eigen::PlainObjectBase<DerivedU> &UV,
  1278. Eigen::PlainObjectBase<DerivedF> &FUV,
  1279. double GradientSize,
  1280. double Stiffness,
  1281. bool DirectRound,
  1282. int iter,
  1283. int localIter,
  1284. bool DoRound,
  1285. bool SingularityRound,
  1286. std::vector<int> roundVertices,
  1287. std::vector<std::vector<int> > hardFeatures)
  1288. {
  1289. GradientSize = GradientSize/(V.colwise().maxCoeff()-V.colwise().minCoeff()).norm();
  1290. igl::comiso::MIQ_class<DerivedV, DerivedF, DerivedU> miq(V,
  1291. F,
  1292. PD1_combed,
  1293. PD2_combed,
  1294. // BIS1_combed,
  1295. // BIS2_combed,
  1296. Handle_MMatch,
  1297. Handle_Singular,
  1298. // Handle_SingularDegree,
  1299. Handle_Seams,
  1300. UV,
  1301. FUV,
  1302. GradientSize,
  1303. Stiffness,
  1304. DirectRound,
  1305. iter,
  1306. localIter,
  1307. DoRound,
  1308. SingularityRound,
  1309. roundVertices,
  1310. hardFeatures);
  1311. miq.extractUV(UV,FUV);
  1312. }
  1313. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1314. IGL_INLINE void igl::comiso::miq(
  1315. const Eigen::PlainObjectBase<DerivedV> &V,
  1316. const Eigen::PlainObjectBase<DerivedF> &F,
  1317. const Eigen::PlainObjectBase<DerivedV> &PD1,
  1318. const Eigen::PlainObjectBase<DerivedV> &PD2,
  1319. Eigen::PlainObjectBase<DerivedU> &UV,
  1320. Eigen::PlainObjectBase<DerivedF> &FUV,
  1321. double GradientSize,
  1322. double Stiffness,
  1323. bool DirectRound,
  1324. int iter,
  1325. int localIter,
  1326. bool DoRound,
  1327. bool SingularityRound,
  1328. std::vector<int> roundVertices,
  1329. std::vector<std::vector<int> > hardFeatures)
  1330. {
  1331. // Eigen::MatrixXd PD2i = PD2;
  1332. // if (PD2i.size() == 0)
  1333. // {
  1334. // Eigen::MatrixXd B1, B2, B3;
  1335. // igl::local_basis(V,F,B1,B2,B3);
  1336. // PD2i = igl::rotate_vectors(V,Eigen::VectorXd::Constant(1,M_PI/2),B1,B2);
  1337. // }
  1338. Eigen::PlainObjectBase<DerivedV> BIS1, BIS2;
  1339. igl::compute_frame_field_bisectors(V, F, PD1, PD2, BIS1, BIS2);
  1340. Eigen::PlainObjectBase<DerivedV> BIS1_combed, BIS2_combed;
  1341. igl::comb_cross_field(V, F, BIS1, BIS2, BIS1_combed, BIS2_combed);
  1342. Eigen::PlainObjectBase<DerivedF> Handle_MMatch;
  1343. igl::cross_field_missmatch(V, F, BIS1_combed, BIS2_combed, true, Handle_MMatch);
  1344. Eigen::Matrix<int, Eigen::Dynamic, 1> isSingularity, singularityIndex;
  1345. igl::find_cross_field_singularities(V, F, Handle_MMatch, isSingularity, singularityIndex);
  1346. Eigen::Matrix<int, Eigen::Dynamic, 3> Handle_Seams;
  1347. igl::cut_mesh_from_singularities(V, F, Handle_MMatch, Handle_Seams);
  1348. Eigen::PlainObjectBase<DerivedV> PD1_combed, PD2_combed;
  1349. igl::comb_frame_field(V, F, PD1, PD2, BIS1_combed, BIS2_combed, PD1_combed, PD2_combed);
  1350. igl::comiso::miq(V,
  1351. F,
  1352. PD1_combed,
  1353. PD2_combed,
  1354. // BIS1_combed,
  1355. // BIS2_combed,
  1356. Handle_MMatch,
  1357. isSingularity,
  1358. // singularityIndex,
  1359. Handle_Seams,
  1360. UV,
  1361. FUV,
  1362. GradientSize,
  1363. Stiffness,
  1364. DirectRound,
  1365. iter,
  1366. localIter,
  1367. DoRound,
  1368. SingularityRound,
  1369. roundVertices,
  1370. hardFeatures);
  1371. }
  1372. #ifdef IGL_STATIC_LIBRARY
  1373. // Explicit template specialization
  1374. 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> > > >);
  1375. 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> > > >);
  1376. 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> > > >);
  1377. #endif