cotmatrix.cpp 8.2 KB

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