cotmatrix.cpp 2.5 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2013 Alec Jacobson <alecjacobson@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 "cotmatrix.h"
  9. // For error printing
  10. #include <cstdio>
  11. #include "cotangent.h"
  12. // Bug in unsupported/Eigen/SparseExtra needs iostream first
  13. #include <iostream>
  14. #include <unsupported/Eigen/SparseExtra>
  15. template <typename DerivedV, typename DerivedF, typename Scalar>
  16. IGL_INLINE void igl::cotmatrix(
  17. const Eigen::PlainObjectBase<DerivedV> & V,
  18. const Eigen::PlainObjectBase<DerivedF> & F,
  19. Eigen::SparseMatrix<Scalar>& L)
  20. {
  21. using namespace igl;
  22. using namespace Eigen;
  23. Eigen::SparseMatrix<double> foo;
  24. SparseMatrix<Scalar, RowMajor> dyn_L (V.rows(), V.rows());
  25. Matrix<int,Dynamic,2> edges;
  26. int simplex_size = F.cols();
  27. // 3 for triangles, 4 for tets
  28. assert(simplex_size == 3 || simplex_size == 4);
  29. if(simplex_size == 3)
  30. {
  31. // This is important! it could decrease the comptuation time by a factor of 2
  32. // Laplacian for a closed 2d manifold mesh will have on average 7 entries per
  33. // row
  34. dyn_L.reserve(7*V.rows());
  35. edges.resize(3,2);
  36. edges <<
  37. 1,2,
  38. 2,0,
  39. 0,1;
  40. }else if(simplex_size == 4)
  41. {
  42. dyn_L.reserve(17*V.rows());
  43. edges.resize(6,2);
  44. edges <<
  45. 1,2,
  46. 2,0,
  47. 0,1,
  48. 3,0,
  49. 3,1,
  50. 3,2;
  51. }else
  52. {
  53. return;
  54. }
  55. // Gather cotangents
  56. Matrix<Scalar,Dynamic,Dynamic> C;
  57. cotangent(V,F,C);
  58. // Loop over triangles
  59. for(int i = 0; i < F.rows(); i++)
  60. {
  61. // loop over edges of element
  62. for(int e = 0;e<edges.rows();e++)
  63. {
  64. int source = F(i,edges(e,0));
  65. int dest = F(i,edges(e,1));
  66. dyn_L.coeffRef(source,dest) += C(i,e);
  67. dyn_L.coeffRef(dest,source) += C(i,e);
  68. dyn_L.coeffRef(source,source) += -C(i,e);
  69. dyn_L.coeffRef(dest,dest) += -C(i,e);
  70. }
  71. }
  72. // Corner indices of this triangle
  73. L = SparseMatrix<Scalar>(dyn_L);
  74. }
  75. #ifndef IGL_HEADER_ONLY
  76. // Explicit template specialization
  77. // generated by autoexplicit.sh
  78. template void igl::cotmatrix<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, double>(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::SparseMatrix<double, 0, int>&);
  79. #endif