conjugate_frame_fields.cpp 31 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2014 Olga Diamanti <olga.diam@gmail.com>
  4. //
  5. // This Source Code Form is subject to the terms of the Mozilla Public License
  6. // v. 2.0. If a copy of the MPL was not distributed with this file, You can
  7. // obtain one at http://mozilla.org/MPL/2.0/.
  8. #include <igl/conjugate_frame_fields.h>
  9. #include <igl/edge_topology.h>
  10. #include <igl/local_basis.h>
  11. #include <igl/nchoosek.h>
  12. #include <igl/sparse.h>
  13. #include <igl/speye.h>
  14. #include <igl/slice.h>
  15. #include <igl/polyroots.h>
  16. #include <igl/colon.h>
  17. #include <igl/false_barycentric_subdivision.h>
  18. #include <igl/principal_curvature.h>
  19. #include <Eigen/Sparse>
  20. #include <iostream>
  21. namespace igl {
  22. template <typename DerivedV, typename DerivedF>
  23. class ConjugateFFSolverData
  24. {
  25. public:
  26. const Eigen::PlainObjectBase<DerivedV> &V; int numV;
  27. const Eigen::PlainObjectBase<DerivedF> &F; int numF;
  28. Eigen::MatrixXi EV; int numE;
  29. Eigen::MatrixXi F2E;
  30. Eigen::MatrixXi E2F;
  31. Eigen::VectorXd K;
  32. Eigen::VectorXi isBorderEdge;
  33. int numInteriorEdges;
  34. Eigen::Matrix<int,Eigen::Dynamic,2> E2F_int;
  35. Eigen::VectorXi indInteriorToFull;
  36. Eigen::VectorXi indFullToInterior;
  37. Eigen::PlainObjectBase<DerivedV> B1, B2, FN;
  38. Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic,1> kmin, kmax;
  39. Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic,2> dmin, dmax;
  40. Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic,3> dmin3, dmax3;
  41. Eigen::VectorXd nonPlanarityMeasure;
  42. Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar> > planarityWeight;
  43. //conjugacy matrix
  44. std::vector<Eigen::Matrix<typename DerivedV::Scalar, 4,4> > H;
  45. //conjugacy matrix eigenvectors and (scaled) eigenvalues
  46. std::vector<Eigen::Matrix<typename DerivedV::Scalar, 4,4> > UH;
  47. std::vector<Eigen::Matrix<typename DerivedV::Scalar, 4,1> > s;
  48. //laplacians
  49. Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar>> DDA, DDB;
  50. private:
  51. IGL_INLINE void computeCurvatureAndPrincipals();
  52. IGL_INLINE void precomputeConjugacyStuff();
  53. IGL_INLINE void computeLaplacians();
  54. IGL_INLINE void computek();
  55. IGL_INLINE void computeCoefficientLaplacian(int n, Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar> > &D);
  56. IGL_INLINE void precomputeInteriorEdges();
  57. public:
  58. IGL_INLINE ConjugateFFSolverData(const Eigen::PlainObjectBase<DerivedV> &_V,
  59. const Eigen::PlainObjectBase<DerivedF> &_F);
  60. IGL_INLINE void evaluateConjugacy(const Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 2> &pvU,
  61. const Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 2> &pvV,
  62. Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 1> &conjValues) const ;
  63. };
  64. template <typename DerivedV, typename DerivedF, typename DerivedO>
  65. class ConjugateFFSolver
  66. {
  67. public:
  68. IGL_INLINE ConjugateFFSolver(const ConjugateFFSolverData<DerivedV, DerivedF> &_data,
  69. int _maxIter = 50,
  70. const typename DerivedV::Scalar &_lambdaOrtho = .1,
  71. const typename DerivedV::Scalar &_lambdaInit = 100,
  72. const typename DerivedV::Scalar &_lambdaMultFactor = 1.01,
  73. bool _doHardConstraints = true);
  74. IGL_INLINE bool solve(const Eigen::VectorXi &isConstrained,
  75. const Eigen::PlainObjectBase<DerivedO> &initialSolution,
  76. Eigen::PlainObjectBase<DerivedO> &output,
  77. typename DerivedV::Scalar *lambdaOut = NULL);
  78. private:
  79. const ConjugateFFSolverData<DerivedV, DerivedF> &data;
  80. //polyVF data
  81. Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> Acoeff, Bcoeff;
  82. Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 2> pvU, pvV;
  83. typename DerivedV::Scalar lambda;
  84. //parameters
  85. typename DerivedV::Scalar lambdaOrtho;
  86. typename DerivedV::Scalar lambdaInit,lambdaMultFactor;
  87. int maxIter;
  88. bool doHardConstraints;
  89. IGL_INLINE void localStep();
  90. IGL_INLINE void getPolyCoeffsForLocalSolve(const Eigen::Matrix<typename DerivedV::Scalar, 4, 1> &s,
  91. const Eigen::Matrix<typename DerivedV::Scalar, 4, 1> &z,
  92. Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 1> &polyCoeff);
  93. IGL_INLINE void globalStep(const Eigen::Matrix<int, Eigen::Dynamic, 1> &isConstrained,
  94. const Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> &Ak,
  95. const Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> &Bk);
  96. IGL_INLINE void minQuadWithKnownMini(const Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar> > &Q,
  97. const Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar> > &f,
  98. const Eigen::VectorXi isConstrained,
  99. const Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> &xknown,
  100. Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> &x);
  101. IGL_INLINE void setFieldFromCoefficients();
  102. IGL_INLINE void setCoefficientsFromField();
  103. };
  104. }
  105. //Implementation
  106. /***************************** Data ***********************************/
  107. template <typename DerivedV, typename DerivedF>
  108. IGL_INLINE igl::ConjugateFFSolverData<DerivedV, DerivedF>::
  109. ConjugateFFSolverData(const Eigen::PlainObjectBase<DerivedV> &_V,
  110. const Eigen::PlainObjectBase<DerivedF> &_F):
  111. V(_V),
  112. numV(_V.rows()),
  113. F(_F),
  114. numF(_F.rows())
  115. {
  116. igl::edge_topology(V,F,EV,F2E,E2F);
  117. numE = EV.rows();
  118. precomputeInteriorEdges();
  119. igl::local_basis(V,F,B1,B2,FN);
  120. computek();
  121. computeLaplacians();
  122. computeCurvatureAndPrincipals();
  123. precomputeConjugacyStuff();
  124. };
  125. template <typename DerivedV, typename DerivedF>
  126. IGL_INLINE void igl::ConjugateFFSolverData<DerivedV, DerivedF>::computeCurvatureAndPrincipals()
  127. {
  128. Eigen::MatrixXd VCBary;
  129. Eigen::MatrixXi FCBary;
  130. VCBary.setZero(numV+numF,3);
  131. FCBary.setZero(3*numF,3);
  132. igl::false_barycentric_subdivision(V, F, VCBary, FCBary);
  133. Eigen::MatrixXd dmax3_,dmin3_;
  134. igl::principal_curvature(VCBary, FCBary, dmax3_, dmin3_, kmax, kmin, 5,true);
  135. dmax3 = dmax3_.bottomRows(numF);
  136. dmin3 = dmin3_.bottomRows(numF);
  137. kmax = kmax.bottomRows(numF);
  138. kmin = kmin.bottomRows(numF);
  139. // kmax = dmax3.rowwise().norm();
  140. // kmin = dmin3.rowwise().norm();
  141. dmin3.rowwise().normalize();
  142. dmax3.rowwise().normalize();
  143. dmax.setZero(numF,2);
  144. dmin.setZero(numF,2);
  145. for (int i= 0; i <numF; ++i)
  146. {
  147. if(kmin[i] != kmin[i] || kmax[i] != kmax[i] || (dmin3.row(i).array() != dmin3.row(i).array()).any() || (dmax3.row(i).array() != dmax3.row(i).array()).any())
  148. {
  149. kmin[i] = 0;
  150. kmax[i] = 0;
  151. dmin3.row(i) = B1.row(i);
  152. dmax3.row(i) = B2.row(i);
  153. }
  154. else
  155. {
  156. dmax3.row(i) = (dmax3.row(i) - (dmax3.row(i).dot(FN.row(i)))*FN.row(i)).normalized();
  157. dmin3.row(i) = dmin3.row(i) - (dmin3.row(i).dot(FN.row(i)))*FN.row(i);
  158. dmin3.row(i) = (dmin3.row(i) - (dmin3.row(i).dot(dmax3.row(i)))*dmax3.row(i)).normalized();
  159. if ((dmin3.row(i).cross(dmax3.row(i))).dot(FN.row(i))<0)
  160. dmin3.row(i) = -dmin3.row(i);
  161. }
  162. dmax.row(i) << dmax3.row(i).dot(B1.row(i)), dmax3.row(i).dot(B2.row(i));
  163. dmax.row(i).normalize();
  164. dmin.row(i) << dmin3.row(i).dot(B1.row(i)), dmin3.row(i).dot(B2.row(i));
  165. dmin.row(i).normalize();
  166. }
  167. nonPlanarityMeasure = kmax.cwiseAbs().array()*kmin.cwiseAbs().array();
  168. typename DerivedV::Scalar minP = nonPlanarityMeasure.minCoeff();
  169. typename DerivedV::Scalar maxP = nonPlanarityMeasure.maxCoeff();
  170. nonPlanarityMeasure = (nonPlanarityMeasure.array()-minP)/(maxP-minP);
  171. Eigen::VectorXi I = igl::colon<typename DerivedF::Scalar>(0, numF-1);
  172. igl::sparse(I, I, nonPlanarityMeasure, numF, numF, planarityWeight);
  173. }
  174. template <typename DerivedV, typename DerivedF>
  175. IGL_INLINE void igl::ConjugateFFSolverData<DerivedV, DerivedF>::precomputeConjugacyStuff()
  176. {
  177. H.resize(numF);
  178. UH.resize(numF);
  179. s.resize(numF);
  180. for (int i = 0; i<numF; ++i)
  181. {
  182. //compute conjugacy matrix
  183. typename DerivedV::Scalar e1x = dmin(i,0), e1y = dmin(i,1), e2x = dmax(i,0), e2y = dmax(i,1), k1 = kmin[i], k2 = kmax[i];
  184. H[i]<<
  185. 0, 0, k1*e1x*e1x, k1*e1x*e1y,
  186. 0, 0, k1*e1x*e1y, k1*e1y*e1y,
  187. k2*e2x*e2x, k2*e2x*e2y, 0, 0,
  188. k2*e2x*e2y, k2*e2y*e2y, 0, 0;
  189. Eigen::Matrix<typename DerivedV::Scalar, 4, 4> Ht = H[i].transpose();
  190. H[i] = .5*(H[i]+Ht);
  191. Eigen::EigenSolver<Eigen::Matrix<typename DerivedV::Scalar, 4, 4> > es(H[i]);
  192. s[i] = es.eigenvalues().real();//ok to do this because H symmetric
  193. //scale
  194. s[i] = s[i]/(s[i].cwiseAbs().minCoeff());
  195. UH[i] = es.eigenvectors().real();
  196. }
  197. }
  198. template <typename DerivedV, typename DerivedF>
  199. IGL_INLINE void igl::ConjugateFFSolverData<DerivedV, DerivedF>::computeLaplacians()
  200. {
  201. computeCoefficientLaplacian(2, DDA);
  202. computeCoefficientLaplacian(4, DDB);
  203. }
  204. template<typename DerivedV, typename DerivedF>
  205. IGL_INLINE void igl::ConjugateFFSolverData<DerivedV, DerivedF>::
  206. precomputeInteriorEdges()
  207. {
  208. // Flag border edges
  209. numInteriorEdges = 0;
  210. isBorderEdge.setZero(numE,1);
  211. indFullToInterior = -1.*Eigen::VectorXi::Ones(numE,1);
  212. for(unsigned i=0; i<numE; ++i)
  213. {
  214. if ((E2F(i,0) == -1) || ((E2F(i,1) == -1)))
  215. isBorderEdge[i] = 1;
  216. else
  217. {
  218. indFullToInterior[i] = numInteriorEdges;
  219. numInteriorEdges++;
  220. }
  221. }
  222. E2F_int.resize(numInteriorEdges, 2);
  223. indInteriorToFull.setZero(numInteriorEdges,1);
  224. int ii = 0;
  225. for (int k=0; k<numE; ++k)
  226. {
  227. if (isBorderEdge[k])
  228. continue;
  229. E2F_int.row(ii) = E2F.row(k);
  230. indInteriorToFull[ii] = k;
  231. ii++;
  232. }
  233. }
  234. template<typename DerivedV, typename DerivedF>
  235. IGL_INLINE void igl::ConjugateFFSolverData<DerivedV, DerivedF>::
  236. computeCoefficientLaplacian(int n, Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar> > &D)
  237. {
  238. std::vector<Eigen::Triplet<std::complex<typename DerivedV::Scalar> >> tripletList;
  239. // For every non-border edge
  240. for (unsigned eid=0; eid<numE; ++eid)
  241. {
  242. if (!isBorderEdge[eid])
  243. {
  244. int fid0 = E2F(eid,0);
  245. int fid1 = E2F(eid,1);
  246. tripletList.push_back(Eigen::Triplet<std::complex<typename DerivedV::Scalar> >(fid0,
  247. fid0,
  248. std::complex<typename DerivedV::Scalar>(1.)));
  249. tripletList.push_back(Eigen::Triplet<std::complex<typename DerivedV::Scalar> >(fid1,
  250. fid1,
  251. std::complex<typename DerivedV::Scalar>(1.)));
  252. tripletList.push_back(Eigen::Triplet<std::complex<typename DerivedV::Scalar> >(fid0,
  253. fid1,
  254. -1.*std::polar(1.,-1.*n*K[eid])));
  255. tripletList.push_back(Eigen::Triplet<std::complex<typename DerivedV::Scalar> >(fid1,
  256. fid0,
  257. -1.*std::polar(1.,1.*n*K[eid])));
  258. }
  259. }
  260. D.resize(numF,numF);
  261. D.setFromTriplets(tripletList.begin(), tripletList.end());
  262. }
  263. template<typename DerivedV, typename DerivedF>
  264. IGL_INLINE void igl::ConjugateFFSolverData<DerivedV, DerivedF>::
  265. computek()
  266. {
  267. K.setZero(numE);
  268. // For every non-border edge
  269. for (unsigned eid=0; eid<numE; ++eid)
  270. {
  271. if (!isBorderEdge[eid])
  272. {
  273. int fid0 = E2F(eid,0);
  274. int fid1 = E2F(eid,1);
  275. Eigen::Matrix<typename DerivedV::Scalar, 1, 3> N0 = FN.row(fid0);
  276. Eigen::Matrix<typename DerivedV::Scalar, 1, 3> N1 = FN.row(fid1);
  277. // find common edge on triangle 0 and 1
  278. int fid0_vc = -1;
  279. int fid1_vc = -1;
  280. for (unsigned i=0;i<3;++i)
  281. {
  282. if (F2E(fid0,i) == eid)
  283. fid0_vc = i;
  284. if (F2E(fid1,i) == eid)
  285. fid1_vc = i;
  286. }
  287. assert(fid0_vc != -1);
  288. assert(fid1_vc != -1);
  289. Eigen::Matrix<typename DerivedV::Scalar, 1, 3> common_edge = V.row(F(fid0,(fid0_vc+1)%3)) - V.row(F(fid0,fid0_vc));
  290. common_edge.normalize();
  291. // Map the two triangles in a new space where the common edge is the x axis and the N0 the z axis
  292. Eigen::Matrix<typename DerivedV::Scalar, 3, 3> P;
  293. Eigen::Matrix<typename DerivedV::Scalar, 1, 3> o = V.row(F(fid0,fid0_vc));
  294. Eigen::Matrix<typename DerivedV::Scalar, 1, 3> tmp = -N0.cross(common_edge);
  295. P << common_edge, tmp, N0;
  296. // P.transposeInPlace();
  297. Eigen::Matrix<typename DerivedV::Scalar, 3, 3> V0;
  298. V0.row(0) = V.row(F(fid0,0)) -o;
  299. V0.row(1) = V.row(F(fid0,1)) -o;
  300. V0.row(2) = V.row(F(fid0,2)) -o;
  301. V0 = (P*V0.transpose()).transpose();
  302. Eigen::Matrix<typename DerivedV::Scalar, 3, 3> V1;
  303. V1.row(0) = V.row(F(fid1,0)) -o;
  304. V1.row(1) = V.row(F(fid1,1)) -o;
  305. V1.row(2) = V.row(F(fid1,2)) -o;
  306. V1 = (P*V1.transpose()).transpose();
  307. // compute rotation R such that R * N1 = N0
  308. // i.e. map both triangles to the same plane
  309. double alpha = -atan2(V1((fid1_vc+2)%3,2),V1((fid1_vc+2)%3,1));
  310. Eigen::Matrix<typename DerivedV::Scalar, 3, 3> R;
  311. R << 1, 0, 0,
  312. 0, cos(alpha), -sin(alpha) ,
  313. 0, sin(alpha), cos(alpha);
  314. V1 = (R*V1.transpose()).transpose();
  315. // measure the angle between the reference frames
  316. // k_ij is the angle between the triangle on the left and the one on the right
  317. Eigen::Matrix<typename DerivedV::Scalar, 1, 3> ref0 = V0.row(1) - V0.row(0);
  318. Eigen::Matrix<typename DerivedV::Scalar, 1, 3> ref1 = V1.row(1) - V1.row(0);
  319. ref0.normalize();
  320. ref1.normalize();
  321. double ktemp = atan2(ref1(1),ref1(0)) - atan2(ref0(1),ref0(0));
  322. // just to be sure, rotate ref0 using angle ktemp...
  323. Eigen::Matrix<typename DerivedV::Scalar, 2, 2> R2;
  324. R2 << cos(ktemp), -sin(ktemp), sin(ktemp), cos(ktemp);
  325. Eigen::Matrix<typename DerivedV::Scalar, 1, 2> tmp1 = R2*(ref0.head(2)).transpose();
  326. K[eid] = ktemp;
  327. }
  328. }
  329. }
  330. template<typename DerivedV, typename DerivedF>
  331. IGL_INLINE void igl::ConjugateFFSolverData<DerivedV, DerivedF>::
  332. evaluateConjugacy(const Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 2> &pvU,
  333. const Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 2> &pvV,
  334. Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 1> &conjValues) const
  335. {
  336. conjValues.resize(numF,1);
  337. for (int j =0; j<numF; ++j)
  338. {
  339. Eigen::Matrix<typename DerivedV::Scalar, 4, 1> x; x<<pvU.row(j).transpose(), pvV.row(j).transpose();
  340. conjValues[j] = x.transpose()*H[j]*x;
  341. }
  342. }
  343. /***************************** Solver ***********************************/
  344. template <typename DerivedV, typename DerivedF, typename DerivedO>
  345. IGL_INLINE igl::ConjugateFFSolver<DerivedV, DerivedF, DerivedO>::
  346. ConjugateFFSolver(const ConjugateFFSolverData<DerivedV, DerivedF> &_data,
  347. int _maxIter,
  348. const typename DerivedV::Scalar &_lambdaOrtho,
  349. const typename DerivedV::Scalar &_lambdaInit,
  350. const typename DerivedV::Scalar &_lambdaMultFactor,
  351. bool _doHardConstraints):
  352. data(_data),
  353. lambdaOrtho(_lambdaOrtho),
  354. lambdaInit(_lambdaInit),
  355. maxIter(_maxIter),
  356. lambdaMultFactor(_lambdaMultFactor),
  357. doHardConstraints(_doHardConstraints)
  358. {
  359. Acoeff.resize(data.numF,1);
  360. Bcoeff.resize(data.numF,1);
  361. pvU.setZero(data.numF, 2);
  362. pvV.setZero(data.numF, 2);
  363. };
  364. template<typename DerivedV, typename DerivedF, typename DerivedO>
  365. IGL_INLINE void igl::ConjugateFFSolver<DerivedV, DerivedF, DerivedO>::
  366. getPolyCoeffsForLocalSolve(const Eigen::Matrix<typename DerivedV::Scalar, 4, 1> &s,
  367. const Eigen::Matrix<typename DerivedV::Scalar, 4, 1> &z,
  368. Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 1> &polyCoeff)
  369. {
  370. typename DerivedV::Scalar s0 = s(0);
  371. typename DerivedV::Scalar s1 = s(1);
  372. typename DerivedV::Scalar s2 = s(2);
  373. typename DerivedV::Scalar s3 = s(3);
  374. typename DerivedV::Scalar z0 = z(0);
  375. typename DerivedV::Scalar z1 = z(1);
  376. typename DerivedV::Scalar z2 = z(2);
  377. typename DerivedV::Scalar z3 = z(3);
  378. polyCoeff.resize(7,1);
  379. polyCoeff(0) = s0*s0* s1*s1* s2*s2* s3* z3*z3 + s0*s0* s1*s1* s2* s3*s3* z2*z2 + s0*s0* s1* s2*s2* s3*s3* z1*z1 + s0* s1*s1* s2*s2* s3*s3* z0*z0 ;
  380. polyCoeff(1) = 2* s0*s0* s1*s1* s2* s3* z2*z2 + 2* s0*s0* s1*s1* s2* s3* z3*z3 + 2* s0*s0* s1* s2*s2* s3* z1*z1 + 2* s0*s0* s1* s2*s2* s3* z3*z3 + 2* s0*s0* s1* s2* s3*s3* z1*z1 + 2* s0*s0* s1* s2* s3*s3* z2*z2 + 2* s0* s1*s1* s2*s2* s3* z0*z0 + 2* s0* s1*s1* s2*s2* s3* z3*z3 + 2* s0* s1*s1* s2* s3*s3* z0*z0 + 2* s0* s1*s1* s2* s3*s3* z2*z2 + 2* s0* s1* s2*s2* s3*s3* z0*z0 + 2* s0* s1* s2*s2* s3*s3* z1*z1 ;
  381. polyCoeff(2) = s0*s0* s1*s1* s2* z2*z2 + s0*s0* s1*s1* s3* z3*z3 + s0*s0* s1* s2*s2* z1*z1 + 4* s0*s0* s1* s2* s3* z1*z1 + 4* s0*s0* s1* s2* s3* z2*z2 + 4* s0*s0* s1* s2* s3* z3*z3 + s0*s0* s1* s3*s3* z1*z1 + s0*s0* s2*s2* s3* z3*z3 + s0*s0* s2* s3*s3* z2*z2 + s0* s1*s1* s2*s2* z0*z0 + 4* s0* s1*s1* s2* s3* z0*z0 + 4* s0* s1*s1* s2* s3* z2*z2 + 4* s0* s1*s1* s2* s3* z3*z3 + s0* s1*s1* s3*s3* z0*z0 + 4* s0* s1* s2*s2* s3* z0*z0 + 4* s0* s1* s2*s2* s3* z1*z1 + 4* s0* s1* s2*s2* s3* z3*z3 + 4* s0* s1* s2* s3*s3* z0*z0 + 4* s0* s1* s2* s3*s3* z1*z1 + 4* s0* s1* s2* s3*s3* z2*z2 + s0* s2*s2* s3*s3* z0*z0 + s1*s1* s2*s2* s3* z3*z3 + s1*s1* s2* s3*s3* z2*z2 + s1* s2*s2* s3*s3* z1*z1;
  382. polyCoeff(3) = 2* s0*s0* s1* s2* z1*z1 + 2* s0*s0* s1* s2* z2*z2 + 2* s0*s0* s1* s3* z1*z1 + 2* s0*s0* s1* s3* z3*z3 + 2* s0*s0* s2* s3* z2*z2 + 2* s0*s0* s2* s3* z3*z3 + 2* s0* s1*s1* s2* z0*z0 + 2* s0* s1*s1* s2* z2*z2 + 2* s0* s1*s1* s3* z0*z0 + 2* s0* s1*s1* s3* z3*z3 + 2* s0* s1* s2*s2* z0*z0 + 2* s0* s1* s2*s2* z1*z1 + 8* s0* s1* s2* s3* z0*z0 + 8* s0* s1* s2* s3* z1*z1 + 8* s0* s1* s2* s3* z2*z2 + 8* s0* s1* s2* s3* z3*z3 + 2* s0* s1* s3*s3* z0*z0 + 2* s0* s1* s3*s3* z1*z1 + 2* s0* s2*s2* s3* z0*z0 + 2* s0* s2*s2* s3* z3*z3 + 2* s0* s2* s3*s3* z0*z0 + 2* s0* s2* s3*s3* z2*z2 + 2* s1*s1* s2* s3* z2*z2 + 2* s1*s1* s2* s3* z3*z3 + 2* s1* s2*s2* s3* z1*z1 + 2* s1* s2*s2* s3* z3*z3 + 2* s1* s2* s3*s3* z1*z1 + 2* s1* s2* s3*s3* z2*z2 ;
  383. polyCoeff(4) = s0*s0* s1* z1*z1 + s0*s0* s2* z2*z2 + s0*s0* s3* z3*z3 + s0* s1*s1* z0*z0 + 4* s0* s1* s2* z0*z0 + 4* s0* s1* s2* z1*z1 + 4* s0* s1* s2* z2*z2 + 4* s0* s1* s3* z0*z0 + 4* s0* s1* s3* z1*z1 + 4* s0* s1* s3* z3*z3 + s0* s2*s2* z0*z0 + 4* s0* s2* s3* z0*z0 + 4* s0* s2* s3* z2*z2 + 4* s0* s2* s3* z3*z3 + s0* s3*s3* z0*z0 + s1*s1* s2* z2*z2 + s1*s1* s3* z3*z3 + s1* s2*s2* z1*z1 + 4* s1* s2* s3* z1*z1 + 4* s1* s2* s3* z2*z2 + 4* s1* s2* s3* z3*z3 + s1* s3*s3* z1*z1 + s2*s2* s3* z3*z3 + s2* s3*s3* z2*z2;
  384. polyCoeff(5) = 2* s0* s1* z0*z0 + 2* s0* s1* z1*z1 + 2* s0* s2* z0*z0 + 2* s0* s2* z2*z2 + 2* s0* s3* z0*z0 + 2* s0* s3* z3*z3 + 2* s1* s2* z1*z1 + 2* s1* s2* z2*z2 + 2* s1* s3* z1*z1 + 2* s1* s3* z3*z3 + 2* s2* s3* z2*z2 + 2* s2* s3* z3*z3 ;
  385. polyCoeff(6) = s0* z0*z0 + s1* z1*z1 + s2* z2*z2 + s3* z3*z3;
  386. }
  387. template<typename DerivedV, typename DerivedF, typename DerivedO>
  388. IGL_INLINE void igl::ConjugateFFSolver<DerivedV, DerivedF, DerivedO>::
  389. localStep()
  390. {
  391. for (int j =0; j<data.numF; ++j)
  392. {
  393. Eigen::Matrix<typename DerivedV::Scalar, 4, 1> xproj; xproj << pvU.row(j).transpose(),pvV.row(j).transpose();
  394. Eigen::Matrix<typename DerivedV::Scalar, 4, 1> z = data.UH[j].transpose()*xproj;
  395. Eigen::Matrix<typename DerivedV::Scalar, 4, 1> x;
  396. Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 1> polyCoeff;
  397. getPolyCoeffsForLocalSolve(data.s[j], z, polyCoeff);
  398. Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> roots;
  399. igl::polyRoots<typename DerivedV::Scalar, typename DerivedV::Scalar> (polyCoeff, roots );
  400. // find closest real root to xproj
  401. typename DerivedV::Scalar minDist = 1e10;
  402. for (int i =0; i< 6; ++i)
  403. {
  404. if (fabs(imag(roots[i]))>1e-10)
  405. continue;
  406. Eigen::Matrix<typename DerivedV::Scalar, 4, 4> D = ((Eigen::Matrix<typename DerivedV::Scalar, 4, 1>::Ones()+real(roots(i))*data.s[j]).array().inverse()).matrix().asDiagonal();
  407. Eigen::Matrix<typename DerivedV::Scalar, 4, 1> candidate = data.UH[j]*D*z;
  408. typename DerivedV::Scalar dist = (candidate-xproj).norm();
  409. if (dist<minDist)
  410. {
  411. minDist = dist;
  412. x = candidate;
  413. }
  414. }
  415. pvU.row(j) << x(0),x(1);
  416. pvV.row(j) << x(2),x(3);
  417. }
  418. }
  419. template<typename DerivedV, typename DerivedF, typename DerivedO>
  420. IGL_INLINE void igl::ConjugateFFSolver<DerivedV, DerivedF, DerivedO>::
  421. setCoefficientsFromField()
  422. {
  423. for (int i = 0; i <data.numF; ++i)
  424. {
  425. std::complex<typename DerivedV::Scalar> u(pvU(i,0),pvU(i,1));
  426. std::complex<typename DerivedV::Scalar> v(pvV(i,0),pvV(i,1));
  427. Acoeff(i) = u*u+v*v;
  428. Bcoeff(i) = u*u*v*v;
  429. }
  430. }
  431. template<typename DerivedV, typename DerivedF, typename DerivedO>
  432. IGL_INLINE void igl::ConjugateFFSolver<DerivedV, DerivedF, DerivedO>::
  433. globalStep(const Eigen::Matrix<int, Eigen::Dynamic, 1> &isConstrained,
  434. const Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> &Ak,
  435. const Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> &Bk)
  436. {
  437. setCoefficientsFromField();
  438. Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar> > I;
  439. igl::speye(data.numF, data.numF, I);
  440. Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar> > QA = data.DDA+lambda*data.planarityWeight+lambdaOrtho*I;
  441. Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar> > fA = (-2*lambda*data.planarityWeight*Acoeff).sparseView();
  442. Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar> > QB = data.DDB+lambda*data.planarityWeight;
  443. Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar> > fB = (-2*lambda*data.planarityWeight*Bcoeff).sparseView();
  444. if(doHardConstraints)
  445. {
  446. minQuadWithKnownMini(QA, fA, isConstrained, Ak, Acoeff);
  447. minQuadWithKnownMini(QB, fB, isConstrained, Bk, Bcoeff);
  448. }
  449. else
  450. {
  451. Eigen::Matrix<int, Eigen::Dynamic, 1>isknown_; isknown_.setZero(data.numF,1);
  452. Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> xknown_; xknown_.setZero(0,1);
  453. minQuadWithKnownMini(QA, fA, isknown_, xknown_, Acoeff);
  454. minQuadWithKnownMini(QB, fB, isknown_, xknown_, Bcoeff);
  455. }
  456. setFieldFromCoefficients();
  457. }
  458. template<typename DerivedV, typename DerivedF, typename DerivedO>
  459. IGL_INLINE void igl::ConjugateFFSolver<DerivedV, DerivedF, DerivedO>::
  460. setFieldFromCoefficients()
  461. {
  462. for (int i = 0; i <data.numF; ++i)
  463. {
  464. // poly coefficients: 1, 0, -Acoeff, 0, Bcoeff
  465. // matlab code from roots (given there are no trailing zeros in the polynomial coefficients)
  466. Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> polyCoeff(5,1);
  467. polyCoeff<<1., 0., -Acoeff(i), 0., Bcoeff(i);
  468. Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> roots;
  469. polyRoots<std::complex<typename DerivedV::Scalar>>(polyCoeff,roots);
  470. std::complex<typename DerivedV::Scalar> u = roots[0];
  471. int maxi = -1;
  472. float maxd = -1;
  473. for (int k =1; k<4; ++k)
  474. {
  475. float dist = abs(roots[k]+u);
  476. if (dist>maxd)
  477. {
  478. maxd = dist;
  479. maxi = k;
  480. }
  481. }
  482. std::complex<typename DerivedV::Scalar> v = roots[maxi];
  483. pvU(i,0) = real(u); pvU(i,1) = imag(u);
  484. pvV(i,0) = real(v); pvV(i,1) = imag(v);
  485. }
  486. }
  487. template<typename DerivedV, typename DerivedF, typename DerivedO>
  488. IGL_INLINE void igl::ConjugateFFSolver<DerivedV, DerivedF, DerivedO>::
  489. minQuadWithKnownMini(const Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar> > &Q,
  490. const Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar> > &f,
  491. const Eigen::VectorXi isConstrained,
  492. const Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> &xknown,
  493. Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> &x)
  494. {
  495. int N = Q.rows();
  496. int nc = xknown.rows();
  497. Eigen::VectorXi known; known.setZero(nc,1);
  498. Eigen::VectorXi unknown; unknown.setZero(N-nc,1);
  499. int indk = 0, indu = 0;
  500. for (int i = 0; i<N; ++i)
  501. if (isConstrained[i])
  502. {
  503. known[indk] = i;
  504. indk++;
  505. }
  506. else
  507. {
  508. unknown[indu] = i;
  509. indu++;
  510. }
  511. Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar>> Quu, Quk;
  512. igl::slice(Q,unknown, unknown, Quu);
  513. igl::slice(Q,unknown, known, Quk);
  514. std::vector<typename Eigen::Triplet<std::complex<typename DerivedV::Scalar> > > tripletList;
  515. Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar> > fu(N-nc,1);
  516. igl::slice(f,unknown, Eigen::VectorXi::Zero(1,1), fu);
  517. Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar> > rhs = (Quk*xknown).sparseView()+.5*fu;
  518. Eigen::SparseLU< Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar>>> solver;
  519. solver.compute(-Quu);
  520. if(solver.info()!=Eigen::Success)
  521. {
  522. std::cerr<<"Decomposition failed!"<<std::endl;
  523. return;
  524. }
  525. Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar>> b = solver.solve(rhs);
  526. if(solver.info()!=Eigen::Success)
  527. {
  528. std::cerr<<"Solving failed!"<<std::endl;
  529. return;
  530. }
  531. indk = 0, indu = 0;
  532. x.setZero(N,1);
  533. for (int i = 0; i<N; ++i)
  534. if (isConstrained[i])
  535. x[i] = xknown[indk++];
  536. else
  537. x[i] = b.coeff(indu++,0);
  538. }
  539. template<typename DerivedV, typename DerivedF, typename DerivedO>
  540. IGL_INLINE bool igl::ConjugateFFSolver<DerivedV, DerivedF, DerivedO>::
  541. solve(const Eigen::VectorXi &isConstrained,
  542. const Eigen::PlainObjectBase<DerivedO> &initialSolution,
  543. Eigen::PlainObjectBase<DerivedO> &output,
  544. typename DerivedV::Scalar *lambdaOut)
  545. {
  546. int numConstrained = isConstrained.sum();
  547. // coefficient values
  548. Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> Ak, Bk;
  549. pvU.resize(data.numF,2);
  550. pvV.resize(data.numF,2);
  551. for (int fi = 0; fi <data.numF; ++fi)
  552. {
  553. const Eigen::Matrix<typename DerivedV::Scalar, 1, 3> &b1 = data.B1.row(fi);
  554. const Eigen::Matrix<typename DerivedV::Scalar, 1, 3> &b2 = data.B2.row(fi);
  555. const Eigen::Matrix<typename DerivedV::Scalar, 1, 3> &u3 = initialSolution.block(fi,0,1,3);
  556. const Eigen::Matrix<typename DerivedV::Scalar, 1, 3> &v3 = initialSolution.block(fi,3,1,3);
  557. pvU.row(fi)<< u3.dot(b1), u3.dot(b2);
  558. pvV.row(fi)<< v3.dot(b1), v3.dot(b2);
  559. }
  560. setCoefficientsFromField();
  561. Ak.resize(numConstrained,1);
  562. Bk.resize(numConstrained,1);
  563. int ind = 0;
  564. for (int i = 0; i <data.numF; ++i)
  565. {
  566. if(isConstrained[i])
  567. {
  568. Ak(ind) = Acoeff[i];
  569. Bk(ind) = Bcoeff[i];
  570. ind ++;
  571. }
  572. }
  573. typename DerivedV::Scalar smoothnessValue;
  574. Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 1> conjValues;
  575. typename DerivedV::Scalar meanConj;
  576. typename DerivedV::Scalar maxConj;
  577. data.evaluateConjugacy(pvU, pvV, conjValues);
  578. meanConj = conjValues.cwiseAbs().mean();
  579. maxConj = conjValues.cwiseAbs().maxCoeff();
  580. printf("Initial max non-conjugacy: %.5g\n",maxConj);
  581. smoothnessValue = (Acoeff.adjoint()*data.DDA*Acoeff + Bcoeff.adjoint()*data.DDB*Bcoeff).real()[0];
  582. printf("\n\nInitial smoothness: %.5g\n",smoothnessValue);
  583. lambda = lambdaInit;
  584. bool doit = false;
  585. for (int iter = 0; iter<maxIter; ++iter)
  586. {
  587. printf("\n\n--- Iteration %d ---\n",iter);
  588. typename DerivedV::Scalar oldMeanConj = meanConj;
  589. localStep();
  590. globalStep(isConstrained, Ak, Bk);
  591. smoothnessValue = (Acoeff.adjoint()*data.DDA*Acoeff + Bcoeff.adjoint()*data.DDB*Bcoeff).real()[0];
  592. printf("Smoothness: %.5g\n",smoothnessValue);
  593. data.evaluateConjugacy(pvU, pvV, conjValues);
  594. meanConj = conjValues.cwiseAbs().mean();
  595. maxConj = conjValues.cwiseAbs().maxCoeff();
  596. printf("Mean/Max non-conjugacy: %.5g, %.5g\n",meanConj,maxConj);
  597. typename DerivedV::Scalar diffMeanConj = fabs(oldMeanConj-meanConj);
  598. if (diffMeanConj<1e-4)
  599. doit = true;
  600. if (doit)
  601. lambda = lambda*lambdaMultFactor;
  602. printf(" %d %.5g %.5g\n",iter, smoothnessValue,maxConj);
  603. }
  604. output.setZero(data.numF,6);
  605. for (int fi=0; fi<data.numF; ++fi)
  606. {
  607. const Eigen::Matrix<typename DerivedV::Scalar, 1, 3> &b1 = data.B1.row(fi);
  608. const Eigen::Matrix<typename DerivedV::Scalar, 1, 3> &b2 = data.B2.row(fi);
  609. output.block(fi,0, 1, 3) = pvU(fi,0)*b1 + pvU(fi,1)*b2;
  610. output.block(fi,3, 1, 3) = pvV(fi,0)*b1 + pvV(fi,1)*b2;
  611. }
  612. if (lambdaOut)
  613. *lambdaOut = lambda;
  614. return true;
  615. }
  616. template <typename DerivedV, typename DerivedF, typename DerivedO>
  617. IGL_INLINE void igl::conjugate_frame_fields(const Eigen::PlainObjectBase<DerivedV> &V,
  618. const Eigen::PlainObjectBase<DerivedF> &F,
  619. const Eigen::VectorXi &isConstrained,
  620. const Eigen::PlainObjectBase<DerivedO> &initialSolution,
  621. Eigen::PlainObjectBase<DerivedO> &output,
  622. int maxIter,
  623. const typename DerivedV::Scalar &lambdaOrtho,
  624. const typename DerivedV::Scalar &lambdaInit,
  625. const typename DerivedV::Scalar &lambdaMultFactor,
  626. bool doHardConstraints)
  627. {
  628. igl::ConjugateFFSolverData<DerivedV, DerivedF> csdata(V, F);
  629. igl::ConjugateFFSolver<DerivedV, DerivedF, DerivedO> cs(csdata, maxIter, lambdaOrtho, lambdaInit, lambdaMultFactor, doHardConstraints);
  630. cs.solve(isConstrained, initialSolution, output);
  631. }
  632. template <typename DerivedV, typename DerivedF, typename DerivedO>
  633. IGL_INLINE void igl::conjugate_frame_fields(const igl::ConjugateFFSolverData<DerivedV, DerivedF> &csdata,
  634. const Eigen::VectorXi &isConstrained,
  635. const Eigen::PlainObjectBase<DerivedO> &initialSolution,
  636. Eigen::PlainObjectBase<DerivedO> &output,
  637. int maxIter,
  638. const typename DerivedV::Scalar &lambdaOrtho,
  639. const typename DerivedV::Scalar &lambdaInit,
  640. const typename DerivedV::Scalar &lambdaMultFactor,
  641. bool doHardConstraints,
  642. typename DerivedV::Scalar *lambdaOut)
  643. {
  644. igl::ConjugateFFSolver<DerivedV, DerivedF, DerivedO> cs(csdata, maxIter, lambdaOrtho, lambdaInit, lambdaMultFactor, doHardConstraints);
  645. cs.solve(isConstrained, initialSolution, output, lambdaOut);
  646. }
  647. #ifdef IGL_STATIC_LIBRARY
  648. // Explicit template specialization
  649. #endif