cotmatrix.cpp 3.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990
  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. #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. #include <unsupported/Eigen/SparseExtra>
  16. template <typename DerivedV, typename DerivedF, typename Scalar>
  17. IGL_INLINE void igl::cotmatrix(
  18. const Eigen::PlainObjectBase<DerivedV> & V,
  19. const Eigen::PlainObjectBase<DerivedF> & F,
  20. Eigen::SparseMatrix<Scalar>& L)
  21. {
  22. using namespace Eigen;
  23. using namespace std;
  24. L.resize(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. L.reserve(10*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. 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. cotmatrix_entries(V,F,C);
  58. vector<Triplet<Scalar> > IJV;
  59. IJV.reserve(F.rows()*edges.rows()*4);
  60. // Loop over triangles
  61. for(int i = 0; i < F.rows(); i++)
  62. {
  63. // loop over edges of element
  64. for(int e = 0;e<edges.rows();e++)
  65. {
  66. int source = F(i,edges(e,0));
  67. int dest = F(i,edges(e,1));
  68. IJV.push_back(Triplet<Scalar>(source,dest,C(i,e)));
  69. IJV.push_back(Triplet<Scalar>(dest,source,C(i,e)));
  70. IJV.push_back(Triplet<Scalar>(source,source,-C(i,e)));
  71. IJV.push_back(Triplet<Scalar>(dest,dest,-C(i,e)));
  72. }
  73. }
  74. L.setFromTriplets(IJV.begin(),IJV.end());
  75. }
  76. #ifdef IGL_STATIC_LIBRARY
  77. // Explicit template specialization
  78. // generated by autoexplicit.sh
  79. template void igl::cotmatrix<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 4, 0, -1, 4>, double>(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 4, 0, -1, 4> > const&, Eigen::SparseMatrix<double, 0, int>&);
  80. // generated by autoexplicit.sh
  81. template void igl::cotmatrix<Eigen::Matrix<double, -1, 3, 0, -1, 3>, Eigen::Matrix<int, -1, 3, 0, -1, 3>, double>(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 3, 0, -1, 3> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 3, 0, -1, 3> > const&, Eigen::SparseMatrix<double, 0, int>&);
  82. // generated by autoexplicit.sh
  83. 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>&);
  84. #endif