miq.cpp 56 KB

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