cotmatrix.cpp 8.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2013 Alec Jacobson <[email protected]>
  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 "cotmatrix.h"
  9. #include <vector>
  10. // For error printing
  11. #include <cstdio>
  12. #include "cotmatrix_entries.h"
  13. // Bug in unsupported/Eigen/SparseExtra needs iostream first
  14. #include <iostream>
  15. template <typename DerivedV, typename DerivedF, typename Scalar>
  16. IGL_INLINE void igl::cotmatrix(
  17. const Eigen::MatrixBase<DerivedV> & V,
  18. const Eigen::MatrixBase<DerivedF> & F,
  19. Eigen::SparseMatrix<Scalar>& L)
  20. {
  21. L.resize(V.rows(),V.rows());
  22. Eigen::Matrix<int ,Eigen::Dynamic,2> edges;
  23. int simplex_size = F.cols();
  24. // 3 for triangles, 4 for tets
  25. assert(simplex_size == 3 || simplex_size == 4);
  26. if(simplex_size == 3)
  27. {
  28. // This is important! it could decrease the comptuation time by a factor of 2
  29. // Laplacian for a closed 2d manifold mesh will have on average 7 entries per
  30. // row
  31. L.reserve(10*V.rows());
  32. edges.resize(3,2);
  33. edges <<
  34. 1,2,
  35. 2,0,
  36. 0,1;
  37. }else if(simplex_size == 4)
  38. {
  39. L.reserve(17*V.rows());
  40. edges.resize(6,2);
  41. edges <<
  42. 1,2,
  43. 2,0,
  44. 0,1,
  45. 3,0,
  46. 3,1,
  47. 3,2;
  48. }else
  49. {
  50. return;
  51. }
  52. // Gather cotangents
  53. Eigen::Matrix<Scalar ,Eigen::Dynamic ,Eigen::Dynamic> C;
  54. cotmatrix_entries(V,F,C);
  55. std::vector<Eigen::Triplet<Scalar> > IJV;
  56. IJV.reserve(F.rows()*edges.rows()*4);
  57. // Loop over triangles
  58. for(int i = 0; i < F.rows(); i++)
  59. {
  60. // loop over edges of element
  61. for(int e = 0;e<edges.rows();e++)
  62. {
  63. int source = F(i,edges(e,0));
  64. int dest = F(i,edges(e,1));
  65. IJV.push_back(Eigen::Triplet<Scalar>(source,dest,C(i,e)));
  66. IJV.push_back(Eigen::Triplet<Scalar>(dest,source,C(i,e)));
  67. IJV.push_back(Eigen::Triplet<Scalar>(source,source,-C(i,e)));
  68. IJV.push_back(Eigen::Triplet<Scalar>(dest,dest,-C(i,e)));
  69. }
  70. }
  71. L.setFromTriplets(IJV.begin(),IJV.end());
  72. }
  73. #include "massmatrix.h"
  74. #include "cotmatrix_entries.h"
  75. #include "massmatrix.h"
  76. #include <Eigen/Geometry>
  77. #include <Eigen/QR>
  78. template <
  79. typename DerivedV,
  80. typename DerivedI,
  81. typename DerivedC,
  82. typename Scalar>
  83. IGL_INLINE void igl::cotmatrix(
  84. const Eigen::MatrixBase<DerivedV> & V,
  85. const Eigen::MatrixBase<DerivedI> & I,
  86. const Eigen::MatrixBase<DerivedC> & C,
  87. Eigen::SparseMatrix<Scalar>& L,
  88. Eigen::SparseMatrix<Scalar>& M,
  89. Eigen::SparseMatrix<Scalar>& P)
  90. {
  91. typedef Eigen::Matrix<Scalar,1,3> RowVector3S;
  92. typedef Eigen::Matrix<Scalar,Eigen::Dynamic,Eigen::Dynamic> MatrixXS;
  93. typedef Eigen::Matrix<Scalar,Eigen::Dynamic,1> VectorXS;
  94. typedef Eigen::Index Index;
  95. // number of vertices
  96. const Index n = V.rows();
  97. // number of polyfaces
  98. const Index m = C.size()-1;
  99. assert(V.cols() == 2 || V.cols() == 3);
  100. std::vector<Eigen::Triplet<Scalar> > Lfijv;
  101. std::vector<Eigen::Triplet<Scalar> > Mfijv;
  102. std::vector<Eigen::Triplet<Scalar> > Pijv;
  103. // loop over vertices; set identity for original vertices
  104. for(Index i = 0;i<V.rows();i++) { Pijv.emplace_back(i,i,1); }
  105. // loop over faces
  106. for(Index p = 0;p<C.size()-1;p++)
  107. {
  108. // number of faces/vertices in this simple polygon
  109. const Index np = C(p+1)-C(p);
  110. // Working "local" list of vertices; last vertex is new one
  111. // this needs to have 3 columns so Eigen doesn't complain about cross
  112. // products below.
  113. Eigen::Matrix<Scalar,Eigen::Dynamic,3> X = decltype(X)::Zero(np+1,3);
  114. for(Index i = 0;i<np;i++){ X.row(i).head(V.cols()) = V.row(I(C(p)+i)); };
  115. // determine weights definig position of inserted vertex
  116. {
  117. MatrixXS A = decltype(A)::Zero(np+1,np);
  118. // My equation (38) would be A w = b.
  119. VectorXS b = decltype(b)::Zero(np+1);
  120. for(Index k = 0;k<np;k++)
  121. {
  122. const RowVector3S Xkp1mk = X.row((k+1)%np)-X.row(k);
  123. const RowVector3S Xkp1mkck = Xkp1mk.cross(X.row(k));
  124. for(Index i = 0;i<np;i++)
  125. {
  126. b(i) -= 2.*(X.row(i).cross(Xkp1mk)).dot(Xkp1mkck);
  127. for(Index j = 0;j<np;j++)
  128. {
  129. A(i,j) += 2.*(X.row(j).cross(Xkp1mk)).dot(X.row(i).cross(Xkp1mk));
  130. }
  131. }
  132. }
  133. A.row(np).setConstant(1);
  134. b(np) = 1;
  135. const VectorXS w =
  136. Eigen::CompleteOrthogonalDecomposition<Eigen::MatrixXd>(A).solve(b);
  137. X.row(np) = w.transpose()*X.topRows(np);
  138. // scatter w into new row of P
  139. for(Index i = 0;i<np;i++) { Pijv.emplace_back(n+p,I(C(p)+i),w(i)); }
  140. }
  141. // "local" fan of faces. These could be statically cached, but this will
  142. // not be the bottleneck.
  143. Eigen::MatrixXi F(np,3);
  144. for(Index i = 0;i<np;i++)
  145. {
  146. F(i,0) = i;
  147. F(i,1) = (i+1)%np;
  148. F(i,2) = np;
  149. }
  150. // Cotangent contributions
  151. MatrixXS K;
  152. igl::cotmatrix_entries(X,F,K);
  153. // Massmatrix entried
  154. VectorXS Mp;
  155. {
  156. Eigen::SparseMatrix<Scalar> M;
  157. igl::massmatrix(X,F,igl::MASSMATRIX_TYPE_DEFAULT,M);
  158. Mp = M.diagonal();
  159. }
  160. // Scatter into fine Laplacian and mass matrices
  161. const auto J = [&n,&np,&p,&I,&C](Index i)->Index{return i==np?n+p:I(C(p)+i);};
  162. // Should just build Mf as a vector...
  163. for(Index i = 0;i<np+1;i++) { Mfijv.emplace_back(J(i),J(i),Mp(i)); }
  164. // loop over faces
  165. for(Index f = 0;f<np;f++)
  166. {
  167. for(Index c = 0;c<3;c++)
  168. {
  169. const Index i = F(f,(c+1)%3);
  170. const Index j = F(f,(c+2)%3);
  171. // symmetric off-diagonal
  172. Lfijv.emplace_back(J(i),J(j),K(f,c));
  173. Lfijv.emplace_back(J(j),J(i),K(f,c));
  174. // diagonal
  175. Lfijv.emplace_back(J(i),J(i),-K(f,c));
  176. Lfijv.emplace_back(J(j),J(j),-K(f,c));
  177. }
  178. }
  179. }
  180. P.resize(n+m,n);
  181. P.setFromTriplets(Pijv.begin(),Pijv.end());
  182. Eigen::SparseMatrix<Scalar> Lf(n+m,n+m);
  183. Lf.setFromTriplets(Lfijv.begin(),Lfijv.end());
  184. Eigen::SparseMatrix<Scalar> Mf(n+m,n+m);
  185. Mf.setFromTriplets(Mfijv.begin(),Mfijv.end());
  186. L = P.transpose() * Lf * P;
  187. // "unlumped" M
  188. const Eigen::SparseMatrix<Scalar> PTMP = P.transpose() * Mf * P;
  189. // Lump M
  190. const VectorXS Mdiag = PTMP * VectorXS::Ones(n,1);
  191. M = Eigen::SparseMatrix<Scalar>(Mdiag.asDiagonal());
  192. MatrixXS Vf = P*V;
  193. Eigen::MatrixXi Ff(I.size(),3);
  194. {
  195. Index f = 0;
  196. for(Index p = 0;p<C.size()-1;p++)
  197. {
  198. const Index np = C(p+1)-C(p);
  199. for(Index c = 0;c<np;c++)
  200. {
  201. Ff(f,0) = I(C(p)+c);
  202. Ff(f,1) = I(C(p)+(c+1)%np);
  203. Ff(f,2) = V.rows()+p;
  204. f++;
  205. }
  206. }
  207. assert(f == Ff.rows());
  208. }
  209. }
  210. #ifdef IGL_STATIC_LIBRARY
  211. // Explicit template instantiation
  212. // generated by autoexplicit.sh
  213. template void igl::cotmatrix<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, double>(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::SparseMatrix<double, 0, int>&, Eigen::SparseMatrix<double, 0, int>&, Eigen::SparseMatrix<double, 0, int>&);
  214. // generated by autoexplicit.sh
  215. template void igl::cotmatrix<Eigen::Matrix<double, -1, -1, 1, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, double>(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 1, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::SparseMatrix<double, 0, int>&);
  216. template void igl::cotmatrix<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 4, 0, -1, 4>, double>(Eigen::MatrixBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 4, 0, -1, 4> > const&, Eigen::SparseMatrix<double, 0, int>&);
  217. template void igl::cotmatrix<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, double>(Eigen::MatrixBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, Eigen::SparseMatrix<double, 0, int>&);
  218. template void igl::cotmatrix<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, double>(Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::SparseMatrix<double, 0, int>&);
  219. #endif