miq.cpp 53 KB

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