miq.cpp 55 KB

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