conjugate_frame_fields.cpp 31 KB

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