miq.cpp 53 KB

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