miq.cpp 55 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365136613671368136913701371137213731374137513761377137813791380138113821383138413851386138713881389139013911392139313941395139613971398139914001401140214031404140514061407140814091410141114121413141414151416141714181419142014211422142314241425142614271428142914301431143214331434143514361437143814391440144114421443144414451446144714481449145014511452145314541455145614571458145914601461146214631464146514661467146814691470147114721473147414751476147714781479148014811482148314841485148614871488148914901491149214931494149514961497149814991500150115021503150415051506150715081509151015111512151315141515151615171518151915201521152215231524152515261527152815291530153115321533153415351536153715381539154015411542154315441545154615471548154915501551155215531554155515561557155815591560156115621563156415651566156715681569157015711572157315741575157615771578157915801581158215831584158515861587
  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. #include <igl/cut_mesh.h>
  12. // includes for VertexIndexing
  13. #include <igl/HalfEdgeIterator.h>
  14. #include <igl/is_border_vertex.h>
  15. #include <igl/vertex_triangle_adjacency.h>
  16. // includes for PoissonSolver
  17. #include <igl/slice_into.h>
  18. #include <igl/grad.h>
  19. #include <igl/cotmatrix.h>
  20. #include <igl/doublearea.h>
  21. #include <gmm/gmm.h>
  22. #include <CoMISo/Solver/ConstrainedSolver.hh>
  23. #include <CoMISo/Solver/MISolver.hh>
  24. #include <CoMISo/Solver/GMM_Tools.hh>
  25. //
  26. #include <igl/cross_field_missmatch.h>
  27. #include <igl/comb_frame_field.h>
  28. #include <igl/comb_cross_field.h>
  29. #include <igl/cut_mesh_from_singularities.h>
  30. #include <igl/find_cross_field_singularities.h>
  31. #include <igl/compute_frame_field_bisectors.h>
  32. #include <igl/rotate_vectors.h>
  33. #ifndef NDEBUG
  34. #include <fstream>
  35. #endif
  36. #include <iostream>
  37. #include <igl/matlab_format.h>
  38. 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;
  46. int integerVar;
  47. unsigned char MMatch;
  48. IGL_INLINE SeamInfo(int _v0,
  49. int _v0p,
  50. int _MMatch,
  51. int _integerVar);
  52. IGL_INLINE SeamInfo(const SeamInfo &S1);
  53. };
  54. struct MeshSystemInfo
  55. {
  56. ////number of vertices variables
  57. int num_vert_variables;
  58. ///num of integer for cuts
  59. int num_integer_cuts;
  60. ///this are used for drawing purposes
  61. std::vector<SeamInfo> EdgeSeamInfo;
  62. #if 0
  63. ///this are values of integer variables after optimization
  64. std::vector<int> IntegerValues;
  65. #endif
  66. };
  67. template <typename DerivedV, typename DerivedF>
  68. class VertexIndexing
  69. {
  70. public:
  71. // Input:
  72. const Eigen::PlainObjectBase<DerivedV> &V;
  73. const Eigen::PlainObjectBase<DerivedF> &F;
  74. const Eigen::PlainObjectBase<DerivedV> &Vcut;
  75. const Eigen::PlainObjectBase<DerivedF> &Fcut;
  76. const Eigen::PlainObjectBase<DerivedF> &TT;
  77. const Eigen::PlainObjectBase<DerivedF> &TTi;
  78. const Eigen::Matrix<int, Eigen::Dynamic, 3> &Handle_MMatch;
  79. const Eigen::Matrix<int, Eigen::Dynamic, 1> &Handle_Singular; // bool
  80. const Eigen::Matrix<int, Eigen::Dynamic, 3> &Handle_Seams; // 3 bool
  81. ///this handle for mesh TODO: move with the other global variables
  82. MeshSystemInfo Handle_SystemInfo;
  83. IGL_INLINE VertexIndexing(const Eigen::PlainObjectBase<DerivedV> &_V,
  84. const Eigen::PlainObjectBase<DerivedF> &_F,
  85. const Eigen::PlainObjectBase<DerivedV> &_Vcut,
  86. const Eigen::PlainObjectBase<DerivedF> &_Fcut,
  87. const Eigen::PlainObjectBase<DerivedF> &_TT,
  88. const Eigen::PlainObjectBase<DerivedF> &_TTi,
  89. const Eigen::Matrix<int, Eigen::Dynamic, 3> &_Handle_MMatch,
  90. const Eigen::Matrix<int, Eigen::Dynamic, 1> &_Handle_Singular,
  91. const Eigen::Matrix<int, Eigen::Dynamic, 3> &_Handle_Seams
  92. );
  93. // provide information about every vertex per seam
  94. IGL_INLINE void InitSeamInfo();
  95. private:
  96. struct VertexInfo{
  97. int v; // vertex index (according to V)
  98. int f0, k0; // face and local edge information of the edge that connects this vertex to the previous vertex (previous in the vector)
  99. int f1, k1; // face and local edge information of the other face corresponding to the same edge
  100. VertexInfo(int _v, int _f0, int _k0, int _f1, int _k1) :
  101. v(_v), f0(_f0), k0(_k0), f1(_f1), k1(_k1){}
  102. bool operator==(VertexInfo const& other){
  103. return other.v == v;
  104. }
  105. };
  106. IGL_INLINE void GetSeamInfo(const int f0,
  107. const int f1,
  108. const int indexE,
  109. int &v0,int &v1,
  110. int &v0p,int &v1p,
  111. unsigned char &_MMatch);
  112. IGL_INLINE std::vector<std::vector<VertexInfo> > GetVerticesPerSeam();
  113. };
  114. template <typename DerivedV, typename DerivedF>
  115. class PoissonSolver
  116. {
  117. public:
  118. IGL_INLINE void SolvePoisson(Eigen::VectorXd Stiffness,
  119. double vector_field_scale=0.1f,
  120. double grid_res=1.f,
  121. bool direct_round=true,
  122. int localIter=0,
  123. bool _integer_rounding=true,
  124. bool _singularity_rounding=true,
  125. std::vector<int> roundVertices = std::vector<int>(),
  126. std::vector<std::vector<int> > hardFeatures = std::vector<std::vector<int> >());
  127. IGL_INLINE PoissonSolver(const Eigen::PlainObjectBase<DerivedV> &_V,
  128. const Eigen::PlainObjectBase<DerivedF> &_F,
  129. const Eigen::PlainObjectBase<DerivedV> &_Vcut,
  130. const Eigen::PlainObjectBase<DerivedF> &_Fcut,
  131. const Eigen::PlainObjectBase<DerivedF> &_TT,
  132. const Eigen::PlainObjectBase<DerivedF> &_TTi,
  133. const Eigen::PlainObjectBase<DerivedV> &_PD1,
  134. const Eigen::PlainObjectBase<DerivedV> &_PD2,
  135. const Eigen::Matrix<int, Eigen::Dynamic, 1>&_Handle_Singular,
  136. const MeshSystemInfo &_Handle_SystemInfo
  137. );
  138. const Eigen::PlainObjectBase<DerivedV> &V;
  139. const Eigen::PlainObjectBase<DerivedF> &F;
  140. const Eigen::PlainObjectBase<DerivedV> &Vcut;
  141. const Eigen::PlainObjectBase<DerivedF> &Fcut;
  142. const Eigen::PlainObjectBase<DerivedF> &TT;
  143. const Eigen::PlainObjectBase<DerivedF> &TTi;
  144. const Eigen::PlainObjectBase<DerivedV> &PD1;
  145. const Eigen::PlainObjectBase<DerivedV> &PD2;
  146. const Eigen::Matrix<int, Eigen::Dynamic, 1> &Handle_Singular; // bool
  147. const MeshSystemInfo &Handle_SystemInfo;
  148. // Internal:
  149. Eigen::VectorXd Handle_Stiffness;
  150. std::vector<std::vector<int> > VF;
  151. std::vector<std::vector<int> > VFi;
  152. Eigen::MatrixXd UV; // this is probably useless
  153. // Output:
  154. // per wedge UV coordinates, 6 coordinates (1 face) per row
  155. Eigen::MatrixXd WUV;
  156. // per vertex UV coordinates, Vcut.rows() x 2
  157. Eigen::MatrixXd UV_out;
  158. // Matrices
  159. Eigen::SparseMatrix<double> Lhs;
  160. Eigen::SparseMatrix<double> Constraints;
  161. Eigen::VectorXd rhs;
  162. Eigen::VectorXd constraints_rhs;
  163. ///vector of unknowns
  164. std::vector< double > X;
  165. ////REAL PART
  166. ///number of fixed vertex
  167. unsigned int n_fixed_vars;
  168. ///the number of REAL variables for vertices
  169. unsigned int n_vert_vars;
  170. ///total number of variables of the system,
  171. ///do not consider constraints, but consider integer vars
  172. unsigned int num_total_vars;
  173. //////INTEGER PART
  174. ///the total number of integer variables
  175. unsigned int n_integer_vars;
  176. ///CONSTRAINT PART
  177. ///number of cuts constraints
  178. unsigned int num_cut_constraint;
  179. // number of user-defined constraints
  180. unsigned int num_userdefined_constraint;
  181. ///total number of constraints equations
  182. unsigned int num_constraint_equations;
  183. ///vector of blocked vertices
  184. std::vector<int> Hard_constraints;
  185. ///vector of indexes to round
  186. std::vector<int> ids_to_round;
  187. ///vector of indexes to round
  188. std::vector<std::vector<int > > userdefined_constraints;
  189. ///boolean that is true if rounding to integer is needed
  190. bool integer_rounding;
  191. ///START COMMON MATH FUNCTIONS
  192. ///return the complex encoding the rotation
  193. ///for a given missmatch interval
  194. IGL_INLINE std::complex<double> GetRotationComplex(int interval);
  195. ///END COMMON MATH FUNCTIONS
  196. ///START FIXING VERTICES
  197. ///set a given vertex as fixed
  198. IGL_INLINE void AddFixedVertex(int v);
  199. ///find vertex to fix in case we're using
  200. ///a vector field NB: multiple components not handled
  201. IGL_INLINE void FindFixedVertField();
  202. ///find hard constraint depending if using or not
  203. ///a vector field
  204. IGL_INLINE void FindFixedVert();
  205. IGL_INLINE int GetFirstVertexIndex(int v);
  206. ///fix the vertices which are flagged as fixed
  207. IGL_INLINE void FixBlockedVertex();
  208. ///END FIXING VERTICES
  209. ///HANDLING SINGULARITY
  210. //set the singularity round to integer location
  211. IGL_INLINE void AddSingularityRound();
  212. IGL_INLINE void AddToRoundVertices(std::vector<int> ids);
  213. ///START GENERIC SYSTEM FUNCTIONS
  214. //build the laplacian matrix cyclyng over all rangemaps
  215. //and over all faces
  216. IGL_INLINE void BuildLaplacianMatrix(double vfscale=1);
  217. ///find different sized of the system
  218. IGL_INLINE void FindSizes();
  219. IGL_INLINE void AllocateSystem();
  220. ///intitialize the whole matrix
  221. IGL_INLINE void InitMatrix();
  222. ///map back coordinates after that
  223. ///the system has been solved
  224. IGL_INLINE void MapCoords();
  225. ///END GENERIC SYSTEM FUNCTIONS
  226. ///set the constraints for the inter-range cuts
  227. IGL_INLINE void BuildSeamConstraintsExplicitTranslation();
  228. ///set the constraints for the inter-range cuts
  229. IGL_INLINE void BuildUserDefinedConstraints();
  230. ///call of the mixed integer solver
  231. IGL_INLINE void MixedIntegerSolve(double cone_grid_res=1,
  232. bool direct_round=true,
  233. int localIter=0);
  234. IGL_INLINE void clearUserConstraint();
  235. IGL_INLINE void addSharpEdgeConstraint(int fid, int vid);
  236. };
  237. template <typename DerivedV, typename DerivedF, typename DerivedU>
  238. class MIQ_class
  239. {
  240. private:
  241. const Eigen::PlainObjectBase<DerivedV> &V;
  242. const Eigen::PlainObjectBase<DerivedF> &F;
  243. Eigen::PlainObjectBase<DerivedV> Vcut;
  244. Eigen::PlainObjectBase<DerivedF> Fcut;
  245. Eigen::MatrixXd UV_out;
  246. Eigen::PlainObjectBase<DerivedF> FUV_out;
  247. // internal
  248. Eigen::PlainObjectBase<DerivedF> TT;
  249. Eigen::PlainObjectBase<DerivedF> TTi;
  250. // Stiffness per face
  251. Eigen::VectorXd Handle_Stiffness;
  252. Eigen::PlainObjectBase<DerivedV> B1, B2, B3;
  253. public:
  254. IGL_INLINE MIQ_class(const Eigen::PlainObjectBase<DerivedV> &V_,
  255. const Eigen::PlainObjectBase<DerivedF> &F_,
  256. const Eigen::PlainObjectBase<DerivedV> &PD1_combed,
  257. const Eigen::PlainObjectBase<DerivedV> &PD2_combed,
  258. const Eigen::Matrix<int, Eigen::Dynamic, 3> &Handle_MMatch,
  259. const Eigen::Matrix<int, Eigen::Dynamic, 1> &Handle_Singular,
  260. const Eigen::Matrix<int, Eigen::Dynamic, 3> &Handle_Seams,
  261. Eigen::PlainObjectBase<DerivedU> &UV,
  262. Eigen::PlainObjectBase<DerivedF> &FUV,
  263. double GradientSize = 30.0,
  264. double Stiffness = 5.0,
  265. bool DirectRound = false,
  266. int iter = 5,
  267. int localIter = 5,
  268. bool DoRound = true,
  269. bool SingularityRound=true,
  270. std::vector<int> roundVertices = std::vector<int>(),
  271. std::vector<std::vector<int> > hardFeatures = std::vector<std::vector<int> >());
  272. IGL_INLINE void extractUV(Eigen::PlainObjectBase<DerivedU> &UV_out,
  273. Eigen::PlainObjectBase<DerivedF> &FUV_out);
  274. private:
  275. IGL_INLINE int NumFlips(const Eigen::MatrixXd& WUV);
  276. IGL_INLINE double Distortion(int f, double h, const Eigen::MatrixXd& WUV);
  277. IGL_INLINE double LaplaceDistortion(const int f, double h, const Eigen::MatrixXd& WUV);
  278. IGL_INLINE bool updateStiffeningJacobianDistorsion(double grad_size, const Eigen::MatrixXd& WUV);
  279. IGL_INLINE bool IsFlipped(const Eigen::Vector2d &uv0,
  280. const Eigen::Vector2d &uv1,
  281. const Eigen::Vector2d &uv2);
  282. IGL_INLINE bool IsFlipped(const int i, const Eigen::MatrixXd& WUV);
  283. };
  284. };
  285. }
  286. IGL_INLINE igl::comiso::SeamInfo::SeamInfo(int _v0,
  287. int _v0p,
  288. int _MMatch,
  289. int _integerVar)
  290. {
  291. v0=_v0;
  292. v0p=_v0p;
  293. integerVar=_integerVar;
  294. MMatch=_MMatch;
  295. }
  296. IGL_INLINE igl::comiso::SeamInfo::SeamInfo(const SeamInfo &S1)
  297. {
  298. v0=S1.v0;
  299. v0p=S1.v0p;
  300. integerVar=S1.integerVar;
  301. MMatch=S1.MMatch;
  302. }
  303. template <typename DerivedV, typename DerivedF>
  304. IGL_INLINE igl::comiso::VertexIndexing<DerivedV, DerivedF>::VertexIndexing(const Eigen::PlainObjectBase<DerivedV> &_V,
  305. const Eigen::PlainObjectBase<DerivedF> &_F,
  306. const Eigen::PlainObjectBase<DerivedV> &_Vcut,
  307. const Eigen::PlainObjectBase<DerivedF> &_Fcut,
  308. const Eigen::PlainObjectBase<DerivedF> &_TT,
  309. const Eigen::PlainObjectBase<DerivedF> &_TTi,
  310. const Eigen::Matrix<int, Eigen::Dynamic, 3> &_Handle_MMatch,
  311. const Eigen::Matrix<int, Eigen::Dynamic, 1> &_Handle_Singular,
  312. const Eigen::Matrix<int, Eigen::Dynamic, 3> &_Handle_Seams
  313. ):
  314. V(_V),
  315. F(_F),
  316. Vcut(_Vcut),
  317. Fcut(_Fcut),
  318. TT(_TT),
  319. TTi(_TTi),
  320. Handle_MMatch(_Handle_MMatch),
  321. Handle_Singular(_Handle_Singular),
  322. Handle_Seams(_Handle_Seams)
  323. {
  324. #ifdef DEBUG_PRINT
  325. cerr<<igl::matlab_format(Handle_Seams,"Handle_Seams");
  326. #endif
  327. Handle_SystemInfo.num_vert_variables=Vcut.rows();
  328. Handle_SystemInfo.num_integer_cuts=0;
  329. }
  330. template <typename DerivedV, typename DerivedF>
  331. IGL_INLINE void igl::comiso::VertexIndexing<DerivedV, DerivedF>::GetSeamInfo(const int f0,
  332. const int f1,
  333. const int indexE,
  334. int &v0,int &v1,
  335. int &v0p,int &v1p,
  336. unsigned char &_MMatch)
  337. {
  338. int edgef0 = indexE;
  339. v0 = Fcut(f0,edgef0);
  340. v1 = Fcut(f0,(edgef0+1)%3);
  341. ////get the index on opposite side
  342. assert(TT(f0,edgef0) == f1);
  343. int edgef1 = TTi(f0,edgef0);
  344. v1p = Fcut(f1,edgef1);
  345. v0p = Fcut(f1,(edgef1+1)%3);
  346. _MMatch = Handle_MMatch(f0,edgef0);
  347. assert(F(f0,edgef0) == F(f1,((edgef1+1)%3)));
  348. assert(F(f0,((edgef0+1)%3)) == F(f1,edgef1));
  349. }
  350. template <typename DerivedV, typename DerivedF>
  351. IGL_INLINE std::vector<std::vector<typename igl::comiso::VertexIndexing<DerivedV, DerivedF>::VertexInfo> > igl::comiso::VertexIndexing<DerivedV, DerivedF>::GetVerticesPerSeam()
  352. {
  353. // Return value
  354. std::vector<std::vector<VertexInfo> >verticesPerSeam;
  355. // for every vertex, keep track of their adjacent vertices on seams.
  356. // regular vertices have two neighbors on a seam, start- and endvertices may have any other numbers of neighbors (e.g. 1 or 3)
  357. std::vector<std::list<VertexInfo> > VVSeam(V.rows());
  358. Eigen::MatrixXi F_hit = Eigen::MatrixXi::Zero(F.rows(), 3);
  359. for (unsigned int f=0; f<F.rows();f++)
  360. {
  361. int f0 = f;
  362. for(int k0=0; k0<3; k0++){
  363. int f1 = TT(f0,k0);
  364. if(f1 == -1)
  365. continue;
  366. bool seam = Handle_Seams(f0,k0);
  367. if (seam && F_hit(f0,k0) == 0)
  368. {
  369. int v0 = F(f0, k0);
  370. int v1 = F(f0, (k0+1)%3);
  371. int k1 = TTi(f0,k0);
  372. VVSeam[v0].push_back(VertexInfo(v1, f0, k0, f1, k1));
  373. VVSeam[v1].push_back(VertexInfo(v0, f0, k0, f1, k1));
  374. F_hit(f0, k0) = 1;
  375. F_hit(f1, k1) = 1;
  376. }
  377. }
  378. }
  379. // Find start vertices, i.e. vertices that start or end a seam branch
  380. std::vector<int> startVertexIndices;
  381. std::vector<bool> isStartVertex(V.rows());
  382. for (unsigned int i=0;i<V.rows();i++)
  383. {
  384. isStartVertex[i] = false;
  385. // vertices with two neighbors are regular vertices, unless the vertex is a singularity, in which case it qualifies as a start vertex
  386. if (VVSeam[i].size() > 0 && VVSeam[i].size() != 2 || Handle_Singular(i) == true)
  387. {
  388. startVertexIndices.push_back(i);
  389. isStartVertex[i] = true;
  390. }
  391. }
  392. // For each startVertex, walk along its seam
  393. for (unsigned int i=0;i<startVertexIndices.size();i++)
  394. {
  395. auto startVertexNeighbors = &VVSeam[startVertexIndices[i]];
  396. const int neighborSize = startVertexNeighbors->size();
  397. // explore every seam to which this vertex is a start vertex
  398. // note: a vertex can never be a start vertex and a regular vertex simultaneously
  399. for (unsigned int j=0;j<neighborSize;j++)
  400. {
  401. std::vector<VertexInfo> thisSeam; // temporary container
  402. // Create vertexInfo struct for start vertex
  403. auto startVertex = VertexInfo(startVertexIndices[i], -1, -1, -1, -1);// -1 values are arbitrary (will never be used)
  404. auto currentVertex = startVertex;
  405. // Add start vertex to the seam
  406. thisSeam.push_back(currentVertex);
  407. // advance on the seam
  408. auto currentVertexNeighbors = startVertexNeighbors;
  409. auto nextVertex = currentVertexNeighbors->front();
  410. currentVertexNeighbors->pop_front();
  411. auto prevVertex = startVertex; // bogus initialization to get the type
  412. while (true)
  413. {
  414. // move to the next vertex
  415. prevVertex = currentVertex;
  416. currentVertex = nextVertex;
  417. currentVertexNeighbors = &VVSeam[nextVertex.v];
  418. // add current vertex to this seam
  419. thisSeam.push_back(currentVertex);
  420. // remove the previous vertex
  421. auto it = std::find(currentVertexNeighbors->begin(), currentVertexNeighbors->end(), prevVertex);
  422. assert(it != currentVertexNeighbors->end());
  423. currentVertexNeighbors->erase(it);
  424. if (currentVertexNeighbors->size() == 1 && !isStartVertex[currentVertex.v])
  425. {
  426. nextVertex = currentVertexNeighbors->front();
  427. currentVertexNeighbors->pop_front();
  428. }
  429. else
  430. break;
  431. }
  432. verticesPerSeam.push_back(thisSeam);
  433. }
  434. }
  435. return verticesPerSeam;
  436. }
  437. template <typename DerivedV, typename DerivedF>
  438. IGL_INLINE void igl::comiso::VertexIndexing<DerivedV, DerivedF>::InitSeamInfo()
  439. {
  440. auto verticesPerSeam = GetVerticesPerSeam();
  441. Handle_SystemInfo.EdgeSeamInfo.clear();
  442. int integerVar = 0;
  443. // Loop over each seam
  444. for(auto seam : verticesPerSeam){
  445. //choose initial side of the seam such that the start vertex corresponds to Fcut(f, k) and the end vertex corresponds to Fcut(f, (k+1)%3) and not vice versa.
  446. int priorVertexIdx;
  447. if(seam.size() > 2){
  448. auto v1 = seam[1];
  449. auto v2 = seam[2];
  450. if(Fcut(v1.f0, (v1.k0+1) % 3) == Fcut(v2.f0, v2.k0) || Fcut(v1.f0, (v1.k0+1) % 3) == Fcut(v2.f1, v2.k1)){
  451. priorVertexIdx = Fcut(v1.f0, v1.k0);
  452. }
  453. else{
  454. priorVertexIdx = Fcut(v1.f1, v1.k1);
  455. assert(Fcut(v1.f1, (v1.k1+1) % 3) == Fcut(v2.f0, v2.k0) || Fcut(v1.f1, (v1.k1+1) % 3) == Fcut(v2.f1, v2.k1));
  456. }
  457. }
  458. else{
  459. auto v1 = seam[1];
  460. priorVertexIdx = Fcut(v1.f0, v1.k0);
  461. }
  462. // Loop over each vertex of the seam
  463. for(auto it=seam.begin()+1; it != seam.end(); ++it){
  464. auto vertex = *it;
  465. // choose the correct side of the seam
  466. int f,k,ff,kk;
  467. if(priorVertexIdx == Fcut(vertex.f0, vertex.k0)){
  468. f = vertex.f0; ff = vertex.f1;
  469. k = vertex.k0; kk = vertex.k1;
  470. }
  471. else{
  472. f = vertex.f1; ff = vertex.f0;
  473. k = vertex.k1; kk = vertex.k0;
  474. assert(priorVertexIdx == Fcut(vertex.f1, vertex.k1));
  475. }
  476. int vtx0,vtx0p,vtx1,vtx1p;
  477. unsigned char MM;
  478. GetSeamInfo(f,ff,k,vtx0,vtx1,vtx0p,vtx1p,MM);
  479. Handle_SystemInfo.EdgeSeamInfo.push_back(SeamInfo(vtx0,vtx0p,MM,integerVar));
  480. if(it == seam.end() -1){
  481. Handle_SystemInfo.EdgeSeamInfo.push_back(SeamInfo(vtx1,vtx1p,MM,integerVar));
  482. }
  483. priorVertexIdx = vtx1;
  484. }
  485. // use the same integer for each seam
  486. integerVar++;
  487. }
  488. Handle_SystemInfo.num_integer_cuts = integerVar;
  489. #ifndef NDEBUG
  490. int totalNVerticesOnSeams = 0;
  491. for(auto seam : verticesPerSeam){
  492. totalNVerticesOnSeams += seam.size();
  493. }
  494. assert(Handle_SystemInfo.EdgeSeamInfo.size() == totalNVerticesOnSeams);
  495. #endif
  496. }
  497. template <typename DerivedV, typename DerivedF>
  498. IGL_INLINE void igl::comiso::PoissonSolver<DerivedV, DerivedF>::SolvePoisson(Eigen::VectorXd Stiffness,
  499. double vector_field_scale,
  500. double grid_res,
  501. bool direct_round,
  502. int localIter,
  503. bool _integer_rounding,
  504. bool _singularity_rounding,
  505. std::vector<int> roundVertices,
  506. std::vector<std::vector<int> > hardFeatures)
  507. {
  508. Handle_Stiffness = Stiffness;
  509. //initialization of flags and data structures
  510. integer_rounding=_integer_rounding;
  511. ids_to_round.clear();
  512. clearUserConstraint();
  513. // copy the user constraints number
  514. for (size_t i = 0; i < hardFeatures.size(); ++i)
  515. {
  516. addSharpEdgeConstraint(hardFeatures[i][0],hardFeatures[i][1]);
  517. }
  518. ///Initializing Matrix
  519. int t0=clock();
  520. ///initialize the matrix ALLOCATING SPACE
  521. InitMatrix();
  522. if (DEBUGPRINT)
  523. printf("\n ALLOCATED THE MATRIX \n");
  524. ///build the laplacian system
  525. BuildLaplacianMatrix(vector_field_scale);
  526. // add seam constraints
  527. BuildSeamConstraintsExplicitTranslation();
  528. // add user defined constraints
  529. BuildUserDefinedConstraints();
  530. ////add the lagrange multiplier
  531. FixBlockedVertex();
  532. if (DEBUGPRINT)
  533. printf("\n BUILT THE MATRIX \n");
  534. if (integer_rounding)
  535. AddToRoundVertices(roundVertices);
  536. if (_singularity_rounding)
  537. AddSingularityRound();
  538. int t1=clock();
  539. if (DEBUGPRINT) printf("\n time:%d \n",t1-t0);
  540. if (DEBUGPRINT) printf("\n SOLVING \n");
  541. MixedIntegerSolve(grid_res,direct_round,localIter);
  542. int t2=clock();
  543. if (DEBUGPRINT) printf("\n time:%d \n",t2-t1);
  544. if (DEBUGPRINT) printf("\n ASSIGNING COORDS \n");
  545. MapCoords();
  546. int t3=clock();
  547. if (DEBUGPRINT) printf("\n time:%d \n",t3-t2);
  548. if (DEBUGPRINT) printf("\n FINISHED \n");
  549. }
  550. template <typename DerivedV, typename DerivedF>
  551. IGL_INLINE igl::comiso::PoissonSolver<DerivedV, DerivedF>
  552. ::PoissonSolver(const Eigen::PlainObjectBase<DerivedV> &_V,
  553. const Eigen::PlainObjectBase<DerivedF> &_F,
  554. const Eigen::PlainObjectBase<DerivedV> &_Vcut,
  555. const Eigen::PlainObjectBase<DerivedF> &_Fcut,
  556. const Eigen::PlainObjectBase<DerivedF> &_TT,
  557. const Eigen::PlainObjectBase<DerivedF> &_TTi,
  558. const Eigen::PlainObjectBase<DerivedV> &_PD1,
  559. const Eigen::PlainObjectBase<DerivedV> &_PD2,
  560. const Eigen::Matrix<int, Eigen::Dynamic, 1>&_Handle_Singular,
  561. const MeshSystemInfo &_Handle_SystemInfo
  562. ):
  563. V(_V),
  564. F(_F),
  565. Vcut(_Vcut),
  566. Fcut(_Fcut),
  567. TT(_TT),
  568. TTi(_TTi),
  569. PD1(_PD1),
  570. PD2(_PD2),
  571. Handle_Singular(_Handle_Singular),
  572. Handle_SystemInfo(_Handle_SystemInfo)
  573. {
  574. UV = Eigen::MatrixXd(V.rows(),2);
  575. WUV = Eigen::MatrixXd(F.rows(),6);
  576. UV_out = Eigen::MatrixXd(Vcut.rows(),2);
  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();
  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. if (DEBUGPRINT) printf("\n*** SYSTEM VARIABLES *** \n");
  743. if (DEBUGPRINT) printf("* NUM REAL VERTEX VARIABLES %d \n",n_vert_vars);
  744. if (DEBUGPRINT) printf("\n*** INTEGER VARIABLES *** \n");
  745. if (DEBUGPRINT) printf("* NUM INTEGER VARIABLES %d \n",(int)n_integer_vars);
  746. if (DEBUGPRINT) printf("\n*** CONSTRAINTS *** \n ");
  747. if (DEBUGPRINT) printf("* NUM FIXED CONSTRAINTS %d\n",n_fixed_vars);
  748. if (DEBUGPRINT) printf("* NUM CUTS CONSTRAINTS %d\n",num_cut_constraint);
  749. if (DEBUGPRINT) printf("* NUM USER DEFINED CONSTRAINTS %d\n",num_userdefined_constraint);
  750. if (DEBUGPRINT) printf("\n*** TOTAL SIZE *** \n");
  751. if (DEBUGPRINT) printf("* TOTAL VARIABLE SIZE (WITH INTEGER TRASL) %d \n",num_total_vars);
  752. if (DEBUGPRINT) printf("* TOTAL CONSTRAINTS %d \n",num_constraint_equations);
  753. }
  754. template <typename DerivedV, typename DerivedF>
  755. IGL_INLINE void igl::comiso::PoissonSolver<DerivedV, DerivedF>::AllocateSystem()
  756. {
  757. Lhs.resize(n_vert_vars * 2, n_vert_vars * 2);
  758. Constraints.resize(num_constraint_equations, num_total_vars);
  759. rhs.resize(n_vert_vars * 2);
  760. constraints_rhs.resize(num_constraint_equations);
  761. printf("\n INITIALIZED SPARSE MATRIX OF %d x %d \n",n_vert_vars*2, n_vert_vars*2);
  762. printf("\n INITIALIZED SPARSE MATRIX OF %d x %d \n",num_constraint_equations, num_total_vars);
  763. printf("\n INITIALIZED VECTOR OF %d x 1 \n",n_vert_vars*2);
  764. printf("\n INITIALIZED VECTOR OF %d x 1 \n",num_constraint_equations);
  765. }
  766. ///intitialize the whole matrix
  767. template <typename DerivedV, typename DerivedF>
  768. IGL_INLINE void igl::comiso::PoissonSolver<DerivedV, DerivedF>::InitMatrix()
  769. {
  770. FindSizes();
  771. AllocateSystem();
  772. }
  773. ///map back coordinates after that
  774. ///the system has been solved
  775. template <typename DerivedV, typename DerivedF>
  776. IGL_INLINE void igl::comiso::PoissonSolver<DerivedV, DerivedF>::MapCoords()
  777. {
  778. ///map coords to faces
  779. for (unsigned int f=0;f<Fcut.rows();f++)
  780. {
  781. for (int k=0;k<3;k++)
  782. {
  783. //get the index of the variable in the system
  784. int indexUV = Fcut(f,k);
  785. ///then get U and V coords
  786. double U=X[indexUV*2];
  787. double V=X[indexUV*2+1];
  788. WUV(f,k*2 + 0) = U;
  789. WUV(f,k*2 + 1) = V;
  790. }
  791. }
  792. for(int i = 0; i < Vcut.rows(); i++){
  793. UV_out(i,0) = X[i*2];
  794. UV_out(i,1) = X[i*2+1];
  795. }
  796. #if 0
  797. ///initialize the vector of integer variables to return their values
  798. Handle_SystemInfo.IntegerValues.resize(n_integer_vars*2);
  799. int baseIndex = (n_vert_vars)*2;
  800. int endIndex = baseIndex+n_integer_vars*2;
  801. int index = 0;
  802. for (int i=baseIndex; i<endIndex; i++)
  803. {
  804. ///assert that the value is an integer value
  805. double value=X[i];
  806. double diff = value-(int)floor(value+0.5);
  807. assert(diff<0.00000001);
  808. Handle_SystemInfo.IntegerValues[index] = value;
  809. index++;
  810. }
  811. #endif
  812. }
  813. ///END GENERIC SYSTEM FUNCTIONS
  814. ///set the constraints for the inter-range cuts
  815. template <typename DerivedV, typename DerivedF>
  816. IGL_INLINE void igl::comiso::PoissonSolver<DerivedV, DerivedF>::BuildSeamConstraintsExplicitTranslation()
  817. {
  818. ///current constraint row
  819. int constr_row = 0;
  820. for (unsigned int i=0; i<num_cut_constraint; i++)
  821. {
  822. unsigned char interval = Handle_SystemInfo.EdgeSeamInfo[i].MMatch;
  823. if (interval==1)
  824. interval=3;
  825. else
  826. if(interval==3)
  827. interval=1;
  828. int p0 = Handle_SystemInfo.EdgeSeamInfo[i].v0;
  829. int p0p = Handle_SystemInfo.EdgeSeamInfo[i].v0p;
  830. std::complex<double> rot = GetRotationComplex(interval);
  831. ///get the integer variable
  832. int integerVar = n_vert_vars + Handle_SystemInfo.EdgeSeamInfo[i].integerVar;
  833. if (integer_rounding)
  834. {
  835. ids_to_round.push_back(integerVar*2);
  836. ids_to_round.push_back(integerVar*2+1);
  837. }
  838. // cross boundary compatibility conditions
  839. Constraints.coeffRef(constr_row, 2*p0) += rot.real();
  840. Constraints.coeffRef(constr_row, 2*p0+1) += -rot.imag();
  841. Constraints.coeffRef(constr_row+1, 2*p0) += rot.imag();
  842. Constraints.coeffRef(constr_row+1, 2*p0+1) += rot.real();
  843. Constraints.coeffRef(constr_row, 2*p0p) += -1;
  844. Constraints.coeffRef(constr_row+1, 2*p0p+1) += -1;
  845. Constraints.coeffRef(constr_row, 2*integerVar) += 1;
  846. Constraints.coeffRef(constr_row+1, 2*integerVar+1) += 1;
  847. constraints_rhs[constr_row] = 0;
  848. constraints_rhs[constr_row+1] = 0;
  849. constr_row += 2;
  850. }
  851. }
  852. ///set the constraints for the inter-range cuts
  853. template <typename DerivedV, typename DerivedF>
  854. IGL_INLINE void igl::comiso::PoissonSolver<DerivedV, DerivedF>::BuildUserDefinedConstraints()
  855. {
  856. /// the user defined constraints are at the end
  857. int offset_row = num_cut_constraint*2 + n_fixed_vars*2;
  858. ///current constraint row
  859. int constr_row = offset_row;
  860. assert(num_userdefined_constraint == userdefined_constraints.size());
  861. for (unsigned int i=0; i<num_userdefined_constraint; i++)
  862. {
  863. for (unsigned int j=0; j<userdefined_constraints[i].size()-1; ++j)
  864. {
  865. Constraints.coeffRef(constr_row, j) = userdefined_constraints[i][j];
  866. }
  867. constraints_rhs[constr_row] = userdefined_constraints[i][userdefined_constraints[i].size()-1];
  868. constr_row +=1;
  869. }
  870. }
  871. ///call of the mixed integer solver
  872. template <typename DerivedV, typename DerivedF>
  873. IGL_INLINE void igl::comiso::PoissonSolver<DerivedV, DerivedF>::MixedIntegerSolve(double cone_grid_res,
  874. bool direct_round,
  875. int localIter)
  876. {
  877. X = std::vector<double>((n_vert_vars+n_integer_vars)*2);
  878. if (DEBUGPRINT)
  879. printf("\n ALLOCATED X \n");
  880. ///variables part
  881. int ScalarSize = n_vert_vars*2;
  882. int SizeMatrix = (n_vert_vars+n_integer_vars)*2;
  883. ///matrix A
  884. gmm::col_matrix< gmm::wsvector< double > > A(SizeMatrix,SizeMatrix); // lhs matrix variables
  885. ///constraints part
  886. int CsizeX = num_constraint_equations;
  887. int CsizeY = SizeMatrix+1;
  888. gmm::row_matrix< gmm::wsvector< double > > C(CsizeX,CsizeY); // constraints
  889. if (DEBUGPRINT)
  890. printf("\n ALLOCATED QMM STRUCTURES \n");
  891. std::vector<double> B(SizeMatrix,0); // rhs
  892. if (DEBUGPRINT)
  893. printf("\n ALLOCATED RHS STRUCTURES \n");
  894. //// copy LHS
  895. for (int k=0; k < Lhs.outerSize(); ++k){
  896. for (Eigen::SparseMatrix<double>::InnerIterator it(Lhs,k); it; ++it){
  897. int row = it.row();
  898. int col = it.col();
  899. A(row, col) += it.value();
  900. }
  901. }
  902. //// copy Constraints
  903. for (int k=0; k < Constraints.outerSize(); ++k){
  904. for (Eigen::SparseMatrix<double>::InnerIterator it(Constraints,k); it; ++it){
  905. int row = it.row();
  906. int col = it.col();
  907. C(row, col) += it.value();
  908. }
  909. }
  910. if (DEBUGPRINT)
  911. printf("\n SET %d INTEGER VALUES \n",n_integer_vars);
  912. ///add penalization term for integer variables
  913. double penalization = 0.000001;
  914. int offline_index = ScalarSize;
  915. for(unsigned int i = 0; i < (n_integer_vars)*2; ++i)
  916. {
  917. int index=offline_index+i;
  918. A(index,index)=penalization;
  919. }
  920. if (DEBUGPRINT)
  921. printf("\n SET RHS \n");
  922. // copy RHS
  923. for(int i = 0; i < (int)ScalarSize; ++i)
  924. {
  925. B[i] = rhs[i] * cone_grid_res;
  926. }
  927. // copy constraint RHS
  928. if (DEBUGPRINT)
  929. printf("\n SET %d CONSTRAINTS \n",num_constraint_equations);
  930. for(unsigned int i = 0; i < num_constraint_equations; ++i)
  931. {
  932. C(i, SizeMatrix) = -constraints_rhs[i] * cone_grid_res;
  933. }
  934. COMISO::ConstrainedSolver solver;
  935. solver.misolver().set_local_iters(localIter);
  936. solver.misolver().set_direct_rounding(direct_round);
  937. std::sort(ids_to_round.begin(),ids_to_round.end());
  938. std::vector<int>::iterator new_end=std::unique(ids_to_round.begin(),ids_to_round.end());
  939. int dist=distance(ids_to_round.begin(),new_end);
  940. ids_to_round.resize(dist);
  941. solver.solve( C, A, X, B, ids_to_round, 0.0, false, false);
  942. }
  943. template <typename DerivedV, typename DerivedF>
  944. IGL_INLINE void igl::comiso::PoissonSolver<DerivedV, DerivedF>::clearUserConstraint()
  945. {
  946. num_userdefined_constraint = 0;
  947. userdefined_constraints.clear();
  948. }
  949. template <typename DerivedV, typename DerivedF>
  950. IGL_INLINE void igl::comiso::PoissonSolver<DerivedV, DerivedF>::addSharpEdgeConstraint(int fid, int vid)
  951. {
  952. // prepare constraint
  953. std::vector<int> c(Handle_SystemInfo.num_vert_variables*2 + 1);
  954. for (size_t i = 0; i < c.size(); ++i)
  955. {
  956. c[i] = 0;
  957. }
  958. int v1 = Fcut(fid,vid);
  959. int v2 = Fcut(fid,(vid+1)%3);
  960. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> e = Vcut.row(v2) - Vcut.row(v1);
  961. e = e.normalized();
  962. double d1 = fabs(e.dot(PD1.row(fid).normalized()));
  963. double d2 = fabs(e.dot(PD2.row(fid).normalized()));
  964. int offset = 0;
  965. if (d1>d2)
  966. offset = 1;
  967. ids_to_round.push_back((v1 * 2) + offset);
  968. ids_to_round.push_back((v2 * 2) + offset);
  969. // add constraint
  970. c[(v1 * 2) + offset] = 1;
  971. c[(v2 * 2) + offset] = -1;
  972. // add to the user-defined constraints
  973. num_userdefined_constraint++;
  974. userdefined_constraints.push_back(c);
  975. }
  976. template <typename DerivedV, typename DerivedF, typename DerivedU>
  977. IGL_INLINE igl::comiso::MIQ_class<DerivedV, DerivedF, DerivedU>::MIQ_class(const Eigen::PlainObjectBase<DerivedV> &V_,
  978. const Eigen::PlainObjectBase<DerivedF> &F_,
  979. const Eigen::PlainObjectBase<DerivedV> &PD1_combed,
  980. const Eigen::PlainObjectBase<DerivedV> &PD2_combed,
  981. const Eigen::Matrix<int, Eigen::Dynamic, 3> &Handle_MMatch,
  982. const Eigen::Matrix<int, Eigen::Dynamic, 1> &Handle_Singular,
  983. const Eigen::Matrix<int, Eigen::Dynamic, 3> &Handle_Seams,
  984. Eigen::PlainObjectBase<DerivedU> &UV,
  985. Eigen::PlainObjectBase<DerivedF> &FUV,
  986. double GradientSize,
  987. double Stiffness,
  988. bool DirectRound,
  989. int iter,
  990. int localIter,
  991. bool DoRound,
  992. bool SingularityRound,
  993. std::vector<int> roundVertices,
  994. std::vector<std::vector<int> > hardFeatures):
  995. V(V_),
  996. F(F_)
  997. {
  998. igl::cut_mesh(V, F, Handle_Seams, Vcut, Fcut);
  999. igl::local_basis(V,F,B1,B2,B3);
  1000. igl::triangle_triangle_adjacency(V,F,TT,TTi);
  1001. // Prepare indexing for the linear system
  1002. VertexIndexing<DerivedV, DerivedF> VInd(V, F, Vcut, Fcut, TT, TTi, Handle_MMatch, Handle_Singular, Handle_Seams);
  1003. VInd.InitSeamInfo();
  1004. // Eigen::PlainObjectBase<DerivedV> PD1_combed_for_poisson, PD2_combed_for_poisson;
  1005. // // Rotate by 90 degrees CCW
  1006. // PD1_combed_for_poisson.setZero(PD1_combed.rows(),3);
  1007. // PD2_combed_for_poisson.setZero(PD2_combed.rows(),3);
  1008. // for (unsigned i=0; i<PD1_combed.rows();++i)
  1009. // {
  1010. // double n1 = PD1_combed.row(i).norm();
  1011. // double n2 = PD2_combed.row(i).norm();
  1012. //
  1013. // double a1 = atan2(B2.row(i).dot(PD1_combed.row(i)),B1.row(i).dot(PD1_combed.row(i)));
  1014. // double a2 = atan2(B2.row(i).dot(PD2_combed.row(i)),B1.row(i).dot(PD2_combed.row(i)));
  1015. //
  1016. // // a1 += M_PI/2;
  1017. // // a2 += M_PI/2;
  1018. //
  1019. //
  1020. // PD1_combed_for_poisson.row(i) = cos(a1) * B1.row(i) + sin(a1) * B2.row(i);
  1021. // PD2_combed_for_poisson.row(i) = cos(a2) * B1.row(i) + sin(a2) * B2.row(i);
  1022. //
  1023. // PD1_combed_for_poisson.row(i) = PD1_combed_for_poisson.row(i).normalized() * n1;
  1024. // PD2_combed_for_poisson.row(i) = PD2_combed_for_poisson.row(i).normalized() * n2;
  1025. // }
  1026. // Assemble the system and solve
  1027. PoissonSolver<DerivedV, DerivedF> PSolver(V,
  1028. F,
  1029. Vcut,
  1030. Fcut,
  1031. TT,
  1032. TTi,
  1033. PD1_combed,
  1034. PD2_combed,
  1035. Handle_Singular,
  1036. VInd.Handle_SystemInfo);
  1037. Handle_Stiffness = Eigen::VectorXd::Constant(F.rows(),1);
  1038. if (iter > 0) // do stiffening
  1039. {
  1040. for (int i=0;i<iter;i++)
  1041. {
  1042. PSolver.SolvePoisson(Handle_Stiffness, GradientSize,1.f,DirectRound,localIter,DoRound,SingularityRound,roundVertices,hardFeatures);
  1043. int nflips=NumFlips(PSolver.WUV);
  1044. bool folded = updateStiffeningJacobianDistorsion(GradientSize,PSolver.WUV);
  1045. printf("ITERATION %d FLIPS %d \n",i,nflips);
  1046. if (!folded)break;
  1047. }
  1048. }
  1049. else
  1050. {
  1051. PSolver.SolvePoisson(Handle_Stiffness,GradientSize,1.f,DirectRound,localIter,DoRound,SingularityRound,roundVertices,hardFeatures);
  1052. }
  1053. int nflips=NumFlips(PSolver.WUV);
  1054. printf("**** END OPTIMIZING #FLIPS %d ****\n",nflips);
  1055. UV_out = PSolver.UV_out;
  1056. FUV_out = PSolver.Fcut;
  1057. fflush(stdout);
  1058. }
  1059. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1060. IGL_INLINE void igl::comiso::MIQ_class<DerivedV, DerivedF, DerivedU>::extractUV(Eigen::PlainObjectBase<DerivedU> &UV_out,
  1061. Eigen::PlainObjectBase<DerivedF> &FUV_out)
  1062. {
  1063. UV_out = this->UV_out;
  1064. FUV_out = this->FUV_out;
  1065. }
  1066. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1067. IGL_INLINE int igl::comiso::MIQ_class<DerivedV, DerivedF, DerivedU>::NumFlips(const Eigen::MatrixXd& WUV)
  1068. {
  1069. int numFl=0;
  1070. for (unsigned int i=0;i<F.rows();i++)
  1071. {
  1072. if (IsFlipped(i, WUV))
  1073. numFl++;
  1074. }
  1075. return numFl;
  1076. }
  1077. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1078. IGL_INLINE double igl::comiso::MIQ_class<DerivedV, DerivedF, DerivedU>::Distortion(int f, double h, const Eigen::MatrixXd& WUV)
  1079. {
  1080. assert(h > 0);
  1081. Eigen::Vector2d uv0,uv1,uv2;
  1082. uv0 << WUV(f,0), WUV(f,1);
  1083. uv1 << WUV(f,2), WUV(f,3);
  1084. uv2 << WUV(f,4), WUV(f,5);
  1085. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> p0 = V.row(Fcut(f,0));
  1086. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> p1 = V.row(Fcut(f,1));
  1087. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> p2 = V.row(Fcut(f,2));
  1088. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> norm = (p1 - p0).cross(p2 - p0);
  1089. double area2 = norm.norm();
  1090. double area2_inv = 1.0 / area2;
  1091. norm *= area2_inv;
  1092. if (area2 > 0)
  1093. {
  1094. // Singular values of the Jacobian
  1095. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> neg_t0 = norm.cross(p2 - p1);
  1096. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> neg_t1 = norm.cross(p0 - p2);
  1097. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> neg_t2 = norm.cross(p1 - p0);
  1098. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> diffu = (neg_t0 * uv0(0) +neg_t1 *uv1(0) + neg_t2 * uv2(0) )*area2_inv;
  1099. Eigen::Matrix<typename DerivedV::Scalar, 3, 1> diffv = (neg_t0 * uv0(1) + neg_t1*uv1(1) + neg_t2*uv2(1) )*area2_inv;
  1100. // first fundamental form
  1101. double I00 = diffu.dot(diffu); // guaranteed non-neg
  1102. double I01 = diffu.dot(diffv); // I01 = I10
  1103. double I11 = diffv.dot(diffv); // guaranteed non-neg
  1104. // eigenvalues of a 2x2 matrix
  1105. // [a00 a01]
  1106. // [a10 a11]
  1107. // 1/2 * [ (a00 + a11) +/- sqrt((a00 - a11)^2 + 4 a01 a10) ]
  1108. double trI = I00 + I11; // guaranteed non-neg
  1109. double diffDiag = I00 - I11; // guaranteed non-neg
  1110. double sqrtDet = sqrt(std::max(0.0, diffDiag*diffDiag +
  1111. 4 * I01 * I01)); // guaranteed non-neg
  1112. double sig1 = 0.5 * (trI + sqrtDet); // higher singular value
  1113. double sig2 = 0.5 * (trI - sqrtDet); // lower singular value
  1114. // Avoid sig2 < 0 due to numerical error
  1115. if (fabs(sig2) < 1.0e-8)
  1116. sig2 = 0;
  1117. assert(sig1 >= 0);
  1118. assert(sig2 >= 0);
  1119. if (sig2 < 0) {
  1120. printf("Distortion will be NaN! sig1^2 is negative (%lg)\n",
  1121. sig2);
  1122. }
  1123. // The singular values of the Jacobian are the sqrts of the
  1124. // eigenvalues of the first fundamental form.
  1125. sig1 = sqrt(sig1);
  1126. sig2 = sqrt(sig2);
  1127. // distortion
  1128. double tao = IsFlipped(f,WUV) ? -1 : 1;
  1129. double factor = tao / h;
  1130. double lam = fabs(factor * sig1 - 1) + fabs(factor * sig2 - 1);
  1131. return lam;
  1132. }
  1133. else {
  1134. return 10; // something "large"
  1135. }
  1136. }
  1137. ////////////////////////////////////////////////////////////////////////////
  1138. // Approximate the distortion laplacian using a uniform laplacian on
  1139. // the dual mesh:
  1140. // ___________
  1141. // \-1 / \-1 /
  1142. // \ / 3 \ /
  1143. // \-----/
  1144. // \-1 /
  1145. // \ /
  1146. //
  1147. // @param[in] f facet on which to compute distortion laplacian
  1148. // @param[in] h scaling factor applied to cross field
  1149. // @return distortion laplacian for f
  1150. ///////////////////////////////////////////////////////////////////////////
  1151. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1152. IGL_INLINE double igl::comiso::MIQ_class<DerivedV, DerivedF, DerivedU>::LaplaceDistortion(const int f, double h, const Eigen::MatrixXd& WUV)
  1153. {
  1154. double mydist = Distortion(f, h, WUV);
  1155. double lapl=0;
  1156. for (int i=0;i<3;i++)
  1157. {
  1158. if (TT(f,i) != -1)
  1159. lapl += (mydist - Distortion(TT(f,i), h, WUV));
  1160. }
  1161. return lapl;
  1162. }
  1163. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1164. IGL_INLINE bool igl::comiso::MIQ_class<DerivedV, DerivedF, DerivedU>::updateStiffeningJacobianDistorsion(double grad_size, const Eigen::MatrixXd& WUV)
  1165. {
  1166. bool flipped = NumFlips(WUV)>0;
  1167. if (!flipped)
  1168. return false;
  1169. double maxL=0;
  1170. double maxD=0;
  1171. if (flipped)
  1172. {
  1173. const double c = 1.0;
  1174. const double d = 5.0;
  1175. for (unsigned int i = 0; i < Fcut.rows(); ++i)
  1176. {
  1177. double dist=Distortion(i,grad_size,WUV);
  1178. if (dist > maxD)
  1179. maxD=dist;
  1180. double absLap=fabs(LaplaceDistortion(i, grad_size,WUV));
  1181. if (absLap > maxL)
  1182. maxL = absLap;
  1183. double stiffDelta = std::min(c * absLap, d);
  1184. Handle_Stiffness[i]+=stiffDelta;
  1185. }
  1186. }
  1187. printf("Maximum Distorsion %4.4f \n",maxD);
  1188. printf("Maximum Laplacian %4.4f \n",maxL);
  1189. return flipped;
  1190. }
  1191. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1192. IGL_INLINE bool igl::comiso::MIQ_class<DerivedV, DerivedF, DerivedU>::IsFlipped(const Eigen::Vector2d &uv0,
  1193. const Eigen::Vector2d &uv1,
  1194. const Eigen::Vector2d &uv2)
  1195. {
  1196. Eigen::Vector2d e0 = (uv1-uv0);
  1197. Eigen::Vector2d e1 = (uv2-uv0);
  1198. double Area = e0(0)*e1(1) - e0(1)*e1(0);
  1199. return (Area<=0);
  1200. }
  1201. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1202. IGL_INLINE bool igl::comiso::MIQ_class<DerivedV, DerivedF, DerivedU>::IsFlipped(
  1203. const int i, const Eigen::MatrixXd& WUV)
  1204. {
  1205. Eigen::Vector2d uv0,uv1,uv2;
  1206. uv0 << WUV(i,0), WUV(i,1);
  1207. uv1 << WUV(i,2), WUV(i,3);
  1208. uv2 << WUV(i,4), WUV(i,5);
  1209. return (IsFlipped(uv0,uv1,uv2));
  1210. }
  1211. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1212. IGL_INLINE void igl::comiso::miq(
  1213. const Eigen::PlainObjectBase<DerivedV> &V,
  1214. const Eigen::PlainObjectBase<DerivedF> &F,
  1215. const Eigen::PlainObjectBase<DerivedV> &PD1_combed,
  1216. const Eigen::PlainObjectBase<DerivedV> &PD2_combed,
  1217. const Eigen::Matrix<int, Eigen::Dynamic, 3> &Handle_MMatch,
  1218. const Eigen::Matrix<int, Eigen::Dynamic, 1> &Handle_Singular,
  1219. const Eigen::Matrix<int, Eigen::Dynamic, 3> &Handle_Seams,
  1220. Eigen::PlainObjectBase<DerivedU> &UV,
  1221. Eigen::PlainObjectBase<DerivedF> &FUV,
  1222. double GradientSize,
  1223. double Stiffness,
  1224. bool DirectRound,
  1225. int iter,
  1226. int localIter,
  1227. bool DoRound,
  1228. bool SingularityRound,
  1229. std::vector<int> roundVertices,
  1230. std::vector<std::vector<int> > hardFeatures)
  1231. {
  1232. GradientSize = GradientSize/(V.colwise().maxCoeff()-V.colwise().minCoeff()).norm();
  1233. igl::comiso::MIQ_class<DerivedV, DerivedF, DerivedU> miq(V,
  1234. F,
  1235. PD1_combed,
  1236. PD2_combed,
  1237. Handle_MMatch,
  1238. Handle_Singular,
  1239. Handle_Seams,
  1240. UV,
  1241. FUV,
  1242. GradientSize,
  1243. Stiffness,
  1244. DirectRound,
  1245. iter,
  1246. localIter,
  1247. DoRound,
  1248. SingularityRound,
  1249. roundVertices,
  1250. hardFeatures);
  1251. miq.extractUV(UV,FUV);
  1252. }
  1253. template <typename DerivedV, typename DerivedF, typename DerivedU>
  1254. IGL_INLINE void igl::comiso::miq(
  1255. const Eigen::PlainObjectBase<DerivedV> &V,
  1256. const Eigen::PlainObjectBase<DerivedF> &F,
  1257. const Eigen::PlainObjectBase<DerivedV> &PD1,
  1258. const Eigen::PlainObjectBase<DerivedV> &PD2,
  1259. Eigen::PlainObjectBase<DerivedU> &UV,
  1260. Eigen::PlainObjectBase<DerivedF> &FUV,
  1261. double GradientSize,
  1262. double Stiffness,
  1263. bool DirectRound,
  1264. int iter,
  1265. int localIter,
  1266. bool DoRound,
  1267. bool SingularityRound,
  1268. std::vector<int> roundVertices,
  1269. std::vector<std::vector<int> > hardFeatures)
  1270. {
  1271. // Eigen::MatrixXd PD2i = PD2;
  1272. // if (PD2i.size() == 0)
  1273. // {
  1274. // Eigen::MatrixXd B1, B2, B3;
  1275. // igl::local_basis(V,F,B1,B2,B3);
  1276. // PD2i = igl::rotate_vectors(V,Eigen::VectorXd::Constant(1,M_PI/2),B1,B2);
  1277. // }
  1278. Eigen::PlainObjectBase<DerivedV> BIS1, BIS2;
  1279. igl::compute_frame_field_bisectors(V, F, PD1, PD2, BIS1, BIS2);
  1280. Eigen::PlainObjectBase<DerivedV> BIS1_combed, BIS2_combed;
  1281. igl::comb_cross_field(V, F, BIS1, BIS2, BIS1_combed, BIS2_combed);
  1282. Eigen::PlainObjectBase<DerivedF> Handle_MMatch;
  1283. igl::cross_field_missmatch(V, F, BIS1_combed, BIS2_combed, true, Handle_MMatch);
  1284. Eigen::Matrix<int, Eigen::Dynamic, 1> isSingularity, singularityIndex;
  1285. igl::find_cross_field_singularities(V, F, Handle_MMatch, isSingularity, singularityIndex);
  1286. Eigen::Matrix<int, Eigen::Dynamic, 3> Handle_Seams;
  1287. igl::cut_mesh_from_singularities(V, F, Handle_MMatch, Handle_Seams);
  1288. Eigen::PlainObjectBase<DerivedV> PD1_combed, PD2_combed;
  1289. igl::comb_frame_field(V, F, PD1, PD2, BIS1_combed, BIS2_combed, PD1_combed, PD2_combed);
  1290. igl::comiso::miq(V,
  1291. F,
  1292. PD1_combed,
  1293. PD2_combed,
  1294. Handle_MMatch,
  1295. isSingularity,
  1296. Handle_Seams,
  1297. UV,
  1298. FUV,
  1299. GradientSize,
  1300. Stiffness,
  1301. DirectRound,
  1302. iter,
  1303. localIter,
  1304. DoRound,
  1305. SingularityRound,
  1306. roundVertices,
  1307. hardFeatures);
  1308. }
  1309. #ifdef IGL_STATIC_LIBRARY
  1310. // Explicit template specialization
  1311. 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> > > >);
  1312. 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> > > >);
  1313. 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> > > >);
  1314. #endif