miq.cpp 58 KB

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