cotmatrix.cpp 3.1 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889
  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. 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 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. #ifdef IGL_STATIC_LIBRARY
  76. // Explicit template specialization
  77. // generated by autoexplicit.sh
  78. 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>&);
  79. // generated by autoexplicit.sh
  80. 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>&);
  81. // generated by autoexplicit.sh
  82. 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>&);
  83. #endif