miq.cpp 67 KB

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