conjugate_frame_fields.cpp 20 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438
  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/speye.h>
  10. #include <igl/slice.h>
  11. #include <igl/polyroots.h>
  12. #include <Eigen/Sparse>
  13. #include <iostream>
  14. namespace igl {
  15. template <typename DerivedV, typename DerivedF, typename DerivedO>
  16. class ConjugateFFSolver
  17. {
  18. public:
  19. IGL_INLINE ConjugateFFSolver(const ConjugateFFSolverData<DerivedV, DerivedF> &_data,
  20. int _maxIter = 50,
  21. const typename DerivedV::Scalar &_lambdaOrtho = .1,
  22. const typename DerivedV::Scalar &_lambdaInit = 100,
  23. const typename DerivedV::Scalar &_lambdaMultFactor = 1.01,
  24. bool _doHardConstraints = true);
  25. IGL_INLINE bool solve(const Eigen::VectorXi &isConstrained,
  26. const Eigen::PlainObjectBase<DerivedO> &initialSolution,
  27. Eigen::PlainObjectBase<DerivedO> &output,
  28. typename DerivedV::Scalar *lambdaOut = NULL);
  29. private:
  30. const ConjugateFFSolverData<DerivedV, DerivedF> &data;
  31. //polyVF data
  32. Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> Acoeff, Bcoeff;
  33. Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 2> pvU, pvV;
  34. typename DerivedV::Scalar lambda;
  35. //parameters
  36. typename DerivedV::Scalar lambdaOrtho;
  37. typename DerivedV::Scalar lambdaInit,lambdaMultFactor;
  38. int maxIter;
  39. bool doHardConstraints;
  40. IGL_INLINE void localStep();
  41. IGL_INLINE void getPolyCoeffsForLocalSolve(const Eigen::Matrix<typename DerivedV::Scalar, 4, 1> &s,
  42. const Eigen::Matrix<typename DerivedV::Scalar, 4, 1> &z,
  43. Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 1> &polyCoeff);
  44. IGL_INLINE void globalStep(const Eigen::Matrix<int, Eigen::Dynamic, 1> &isConstrained,
  45. const Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> &Ak,
  46. const Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> &Bk);
  47. IGL_INLINE void minQuadWithKnownMini(const Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar> > &Q,
  48. const Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar> > &f,
  49. const Eigen::VectorXi isConstrained,
  50. const Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> &xknown,
  51. Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> &x);
  52. IGL_INLINE void setFieldFromCoefficients();
  53. IGL_INLINE void setCoefficientsFromField();
  54. };
  55. }
  56. //Implementation
  57. /***************************** Solver ***********************************/
  58. template <typename DerivedV, typename DerivedF, typename DerivedO>
  59. IGL_INLINE igl::ConjugateFFSolver<DerivedV, DerivedF, DerivedO>::
  60. ConjugateFFSolver(const ConjugateFFSolverData<DerivedV, DerivedF> &_data,
  61. int _maxIter,
  62. const typename DerivedV::Scalar &_lambdaOrtho,
  63. const typename DerivedV::Scalar &_lambdaInit,
  64. const typename DerivedV::Scalar &_lambdaMultFactor,
  65. bool _doHardConstraints):
  66. data(_data),
  67. lambdaOrtho(_lambdaOrtho),
  68. lambdaInit(_lambdaInit),
  69. maxIter(_maxIter),
  70. lambdaMultFactor(_lambdaMultFactor),
  71. doHardConstraints(_doHardConstraints)
  72. {
  73. Acoeff.resize(data.numF,1);
  74. Bcoeff.resize(data.numF,1);
  75. pvU.setZero(data.numF, 2);
  76. pvV.setZero(data.numF, 2);
  77. };
  78. template<typename DerivedV, typename DerivedF, typename DerivedO>
  79. IGL_INLINE void igl::ConjugateFFSolver<DerivedV, DerivedF, DerivedO>::
  80. getPolyCoeffsForLocalSolve(const Eigen::Matrix<typename DerivedV::Scalar, 4, 1> &s,
  81. const Eigen::Matrix<typename DerivedV::Scalar, 4, 1> &z,
  82. Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 1> &polyCoeff)
  83. {
  84. typename DerivedV::Scalar s0 = s(0);
  85. typename DerivedV::Scalar s1 = s(1);
  86. typename DerivedV::Scalar s2 = s(2);
  87. typename DerivedV::Scalar s3 = s(3);
  88. typename DerivedV::Scalar z0 = z(0);
  89. typename DerivedV::Scalar z1 = z(1);
  90. typename DerivedV::Scalar z2 = z(2);
  91. typename DerivedV::Scalar z3 = z(3);
  92. polyCoeff.resize(7,1);
  93. 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 ;
  94. 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 ;
  95. 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;
  96. 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 ;
  97. 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;
  98. 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 ;
  99. polyCoeff(6) = s0* z0*z0 + s1* z1*z1 + s2* z2*z2 + s3* z3*z3;
  100. }
  101. template<typename DerivedV, typename DerivedF, typename DerivedO>
  102. IGL_INLINE void igl::ConjugateFFSolver<DerivedV, DerivedF, DerivedO>::
  103. localStep()
  104. {
  105. for (int j =0; j<data.numF; ++j)
  106. {
  107. Eigen::Matrix<typename DerivedV::Scalar, 4, 1> xproj; xproj << pvU.row(j).transpose(),pvV.row(j).transpose();
  108. Eigen::Matrix<typename DerivedV::Scalar, 4, 1> z = data.UH[j].transpose()*xproj;
  109. Eigen::Matrix<typename DerivedV::Scalar, 4, 1> x;
  110. Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 1> polyCoeff;
  111. getPolyCoeffsForLocalSolve(data.s[j], z, polyCoeff);
  112. Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> roots;
  113. igl::polyRoots<typename DerivedV::Scalar, typename DerivedV::Scalar> (polyCoeff, roots );
  114. // find closest real root to xproj
  115. typename DerivedV::Scalar minDist = 1e10;
  116. for (int i =0; i< 6; ++i)
  117. {
  118. if (fabs(imag(roots[i]))>1e-10)
  119. continue;
  120. 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();
  121. Eigen::Matrix<typename DerivedV::Scalar, 4, 1> candidate = data.UH[j]*D*z;
  122. typename DerivedV::Scalar dist = (candidate-xproj).norm();
  123. if (dist<minDist)
  124. {
  125. minDist = dist;
  126. x = candidate;
  127. }
  128. }
  129. pvU.row(j) << x(0),x(1);
  130. pvV.row(j) << x(2),x(3);
  131. }
  132. }
  133. template<typename DerivedV, typename DerivedF, typename DerivedO>
  134. IGL_INLINE void igl::ConjugateFFSolver<DerivedV, DerivedF, DerivedO>::
  135. setCoefficientsFromField()
  136. {
  137. for (int i = 0; i <data.numF; ++i)
  138. {
  139. std::complex<typename DerivedV::Scalar> u(pvU(i,0),pvU(i,1));
  140. std::complex<typename DerivedV::Scalar> v(pvV(i,0),pvV(i,1));
  141. Acoeff(i) = u*u+v*v;
  142. Bcoeff(i) = u*u*v*v;
  143. }
  144. }
  145. template<typename DerivedV, typename DerivedF, typename DerivedO>
  146. IGL_INLINE void igl::ConjugateFFSolver<DerivedV, DerivedF, DerivedO>::
  147. globalStep(const Eigen::Matrix<int, Eigen::Dynamic, 1> &isConstrained,
  148. const Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> &Ak,
  149. const Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> &Bk)
  150. {
  151. setCoefficientsFromField();
  152. Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar> > I;
  153. igl::speye(data.numF, data.numF, I);
  154. Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar> > QA = data.DDA+lambda*data.planarityWeight+lambdaOrtho*I;
  155. Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar> > fA = (-2*lambda*data.planarityWeight*Acoeff).sparseView();
  156. Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar> > QB = data.DDB+lambda*data.planarityWeight;
  157. Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar> > fB = (-2*lambda*data.planarityWeight*Bcoeff).sparseView();
  158. if(doHardConstraints)
  159. {
  160. minQuadWithKnownMini(QA, fA, isConstrained, Ak, Acoeff);
  161. minQuadWithKnownMini(QB, fB, isConstrained, Bk, Bcoeff);
  162. }
  163. else
  164. {
  165. Eigen::Matrix<int, Eigen::Dynamic, 1>isknown_; isknown_.setZero(data.numF,1);
  166. Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> xknown_; xknown_.setZero(0,1);
  167. minQuadWithKnownMini(QA, fA, isknown_, xknown_, Acoeff);
  168. minQuadWithKnownMini(QB, fB, isknown_, xknown_, Bcoeff);
  169. }
  170. setFieldFromCoefficients();
  171. }
  172. template<typename DerivedV, typename DerivedF, typename DerivedO>
  173. IGL_INLINE void igl::ConjugateFFSolver<DerivedV, DerivedF, DerivedO>::
  174. setFieldFromCoefficients()
  175. {
  176. for (int i = 0; i <data.numF; ++i)
  177. {
  178. // poly coefficients: 1, 0, -Acoeff, 0, Bcoeff
  179. // matlab code from roots (given there are no trailing zeros in the polynomial coefficients)
  180. Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> polyCoeff(5,1);
  181. polyCoeff<<1., 0., -Acoeff(i), 0., Bcoeff(i);
  182. Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> roots;
  183. polyRoots<std::complex<typename DerivedV::Scalar>>(polyCoeff,roots);
  184. std::complex<typename DerivedV::Scalar> u = roots[0];
  185. int maxi = -1;
  186. float maxd = -1;
  187. for (int k =1; k<4; ++k)
  188. {
  189. float dist = abs(roots[k]+u);
  190. if (dist>maxd)
  191. {
  192. maxd = dist;
  193. maxi = k;
  194. }
  195. }
  196. std::complex<typename DerivedV::Scalar> v = roots[maxi];
  197. pvU(i,0) = real(u); pvU(i,1) = imag(u);
  198. pvV(i,0) = real(v); pvV(i,1) = imag(v);
  199. }
  200. }
  201. template<typename DerivedV, typename DerivedF, typename DerivedO>
  202. IGL_INLINE void igl::ConjugateFFSolver<DerivedV, DerivedF, DerivedO>::
  203. minQuadWithKnownMini(const Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar> > &Q,
  204. const Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar> > &f,
  205. const Eigen::VectorXi isConstrained,
  206. const Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> &xknown,
  207. Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> &x)
  208. {
  209. int N = Q.rows();
  210. int nc = xknown.rows();
  211. Eigen::VectorXi known; known.setZero(nc,1);
  212. Eigen::VectorXi unknown; unknown.setZero(N-nc,1);
  213. int indk = 0, indu = 0;
  214. for (int i = 0; i<N; ++i)
  215. if (isConstrained[i])
  216. {
  217. known[indk] = i;
  218. indk++;
  219. }
  220. else
  221. {
  222. unknown[indu] = i;
  223. indu++;
  224. }
  225. Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar>> Quu, Quk;
  226. igl::slice(Q,unknown, unknown, Quu);
  227. igl::slice(Q,unknown, known, Quk);
  228. std::vector<typename Eigen::Triplet<std::complex<typename DerivedV::Scalar> > > tripletList;
  229. Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar> > fu(N-nc,1);
  230. igl::slice(f,unknown, Eigen::VectorXi::Zero(1,1), fu);
  231. Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar> > rhs = (Quk*xknown).sparseView()+.5*fu;
  232. Eigen::SparseLU< Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar>>> solver;
  233. solver.compute(-Quu);
  234. if(solver.info()!=Eigen::Success)
  235. {
  236. std::cerr<<"Decomposition failed!"<<std::endl;
  237. return;
  238. }
  239. Eigen::SparseMatrix<std::complex<typename DerivedV::Scalar>> b = solver.solve(rhs);
  240. if(solver.info()!=Eigen::Success)
  241. {
  242. std::cerr<<"Solving failed!"<<std::endl;
  243. return;
  244. }
  245. indk = 0, indu = 0;
  246. x.setZero(N,1);
  247. for (int i = 0; i<N; ++i)
  248. if (isConstrained[i])
  249. x[i] = xknown[indk++];
  250. else
  251. x[i] = b.coeff(indu++,0);
  252. }
  253. template<typename DerivedV, typename DerivedF, typename DerivedO>
  254. IGL_INLINE bool igl::ConjugateFFSolver<DerivedV, DerivedF, DerivedO>::
  255. solve(const Eigen::VectorXi &isConstrained,
  256. const Eigen::PlainObjectBase<DerivedO> &initialSolution,
  257. Eigen::PlainObjectBase<DerivedO> &output,
  258. typename DerivedV::Scalar *lambdaOut)
  259. {
  260. int numConstrained = isConstrained.sum();
  261. // coefficient values
  262. Eigen::Matrix<std::complex<typename DerivedV::Scalar>, Eigen::Dynamic, 1> Ak, Bk;
  263. pvU.resize(data.numF,2);
  264. pvV.resize(data.numF,2);
  265. for (int fi = 0; fi <data.numF; ++fi)
  266. {
  267. const Eigen::Matrix<typename DerivedV::Scalar, 1, 3> &b1 = data.B1.row(fi);
  268. const Eigen::Matrix<typename DerivedV::Scalar, 1, 3> &b2 = data.B2.row(fi);
  269. const Eigen::Matrix<typename DerivedV::Scalar, 1, 3> &u3 = initialSolution.block(fi,0,1,3);
  270. const Eigen::Matrix<typename DerivedV::Scalar, 1, 3> &v3 = initialSolution.block(fi,3,1,3);
  271. pvU.row(fi)<< u3.dot(b1), u3.dot(b2);
  272. pvV.row(fi)<< v3.dot(b1), v3.dot(b2);
  273. }
  274. setCoefficientsFromField();
  275. Ak.resize(numConstrained,1);
  276. Bk.resize(numConstrained,1);
  277. int ind = 0;
  278. for (int i = 0; i <data.numF; ++i)
  279. {
  280. if(isConstrained[i])
  281. {
  282. Ak(ind) = Acoeff[i];
  283. Bk(ind) = Bcoeff[i];
  284. ind ++;
  285. }
  286. }
  287. typename DerivedV::Scalar smoothnessValue;
  288. Eigen::Matrix<typename DerivedV::Scalar, Eigen::Dynamic, 1> conjValues;
  289. typename DerivedV::Scalar meanConj;
  290. typename DerivedV::Scalar maxConj;
  291. data.evaluateConjugacy(pvU, pvV, conjValues);
  292. meanConj = conjValues.cwiseAbs().mean();
  293. maxConj = conjValues.cwiseAbs().maxCoeff();
  294. printf("Initial max non-conjugacy: %.5g\n",maxConj);
  295. smoothnessValue = (Acoeff.adjoint()*data.DDA*Acoeff + Bcoeff.adjoint()*data.DDB*Bcoeff).real()[0];
  296. printf("\n\nInitial smoothness: %.5g\n",smoothnessValue);
  297. lambda = lambdaInit;
  298. bool doit = false;
  299. for (int iter = 0; iter<maxIter; ++iter)
  300. {
  301. printf("\n\n--- Iteration %d ---\n",iter);
  302. typename DerivedV::Scalar oldMeanConj = meanConj;
  303. localStep();
  304. globalStep(isConstrained, Ak, Bk);
  305. smoothnessValue = (Acoeff.adjoint()*data.DDA*Acoeff + Bcoeff.adjoint()*data.DDB*Bcoeff).real()[0];
  306. printf("Smoothness: %.5g\n",smoothnessValue);
  307. data.evaluateConjugacy(pvU, pvV, conjValues);
  308. meanConj = conjValues.cwiseAbs().mean();
  309. maxConj = conjValues.cwiseAbs().maxCoeff();
  310. printf("Mean/Max non-conjugacy: %.5g, %.5g\n",meanConj,maxConj);
  311. typename DerivedV::Scalar diffMeanConj = fabs(oldMeanConj-meanConj);
  312. if (diffMeanConj<1e-4)
  313. doit = true;
  314. if (doit)
  315. lambda = lambda*lambdaMultFactor;
  316. printf(" %d %.5g %.5g\n",iter, smoothnessValue,maxConj);
  317. }
  318. output.setZero(data.numF,6);
  319. for (int fi=0; fi<data.numF; ++fi)
  320. {
  321. const Eigen::Matrix<typename DerivedV::Scalar, 1, 3> &b1 = data.B1.row(fi);
  322. const Eigen::Matrix<typename DerivedV::Scalar, 1, 3> &b2 = data.B2.row(fi);
  323. output.block(fi,0, 1, 3) = pvU(fi,0)*b1 + pvU(fi,1)*b2;
  324. output.block(fi,3, 1, 3) = pvV(fi,0)*b1 + pvV(fi,1)*b2;
  325. }
  326. if (lambdaOut)
  327. *lambdaOut = lambda;
  328. return true;
  329. }
  330. template <typename DerivedV, typename DerivedF, typename DerivedO>
  331. IGL_INLINE void igl::conjugate_frame_fields(const Eigen::PlainObjectBase<DerivedV> &V,
  332. const Eigen::PlainObjectBase<DerivedF> &F,
  333. const Eigen::VectorXi &isConstrained,
  334. const Eigen::PlainObjectBase<DerivedO> &initialSolution,
  335. Eigen::PlainObjectBase<DerivedO> &output,
  336. int maxIter,
  337. const typename DerivedV::Scalar &lambdaOrtho,
  338. const typename DerivedV::Scalar &lambdaInit,
  339. const typename DerivedV::Scalar &lambdaMultFactor,
  340. bool doHardConstraints)
  341. {
  342. igl::ConjugateFFSolverData<DerivedV, DerivedF> csdata(V, F);
  343. igl::ConjugateFFSolver<DerivedV, DerivedF, DerivedO> cs(csdata, maxIter, lambdaOrtho, lambdaInit, lambdaMultFactor, doHardConstraints);
  344. cs.solve(isConstrained, initialSolution, output);
  345. }
  346. template <typename DerivedV, typename DerivedF, typename DerivedO>
  347. IGL_INLINE void igl::conjugate_frame_fields(const igl::ConjugateFFSolverData<DerivedV, DerivedF> &csdata,
  348. const Eigen::VectorXi &isConstrained,
  349. const Eigen::PlainObjectBase<DerivedO> &initialSolution,
  350. Eigen::PlainObjectBase<DerivedO> &output,
  351. int maxIter,
  352. const typename DerivedV::Scalar &lambdaOrtho,
  353. const typename DerivedV::Scalar &lambdaInit,
  354. const typename DerivedV::Scalar &lambdaMultFactor,
  355. bool doHardConstraints,
  356. typename DerivedV::Scalar *lambdaOut)
  357. {
  358. igl::ConjugateFFSolver<DerivedV, DerivedF, DerivedO> cs(csdata, maxIter, lambdaOrtho, lambdaInit, lambdaMultFactor, doHardConstraints);
  359. cs.solve(isConstrained, initialSolution, output, lambdaOut);
  360. }
  361. #ifdef IGL_STATIC_LIBRARY
  362. // Explicit template specialization
  363. template void igl::conjugate_frame_fields<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(igl::ConjugateFFSolverData<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, 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> >&, int, Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar const&, Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar const&, Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar const&, bool, Eigen::Matrix<double, -1, -1, 0, -1, -1>::Scalar*);
  364. #endif