miq.cpp 56 KB

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