miq.cpp 55 KB

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