miq.cpp 57 KB

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