heat_geodesics.cpp 5.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2018 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 "heat_geodesics.h"
  9. #include "grad.h"
  10. #include "doublearea.h"
  11. #include "cotmatrix.h"
  12. #include "intrinsic_delaunay_cotmatrix.h"
  13. #include "massmatrix.h"
  14. #include "massmatrix_intrinsic.h"
  15. #include "grad_intrinsic.h"
  16. #include "boundary_facets.h"
  17. #include "unique.h"
  18. #include "slice.h"
  19. #include "avg_edge_length.h"
  20. template < typename DerivedV, typename DerivedF, typename Scalar >
  21. IGL_INLINE bool igl::heat_geodesics_precompute(
  22. const Eigen::MatrixBase<DerivedV> & V,
  23. const Eigen::MatrixBase<DerivedF> & F,
  24. HeatGeodesicsData<Scalar> & data)
  25. {
  26. // default t value
  27. const Scalar h = avg_edge_length(V,F);
  28. const Scalar t = h*h;
  29. return heat_geodesics_precompute(V,F,t,data);
  30. }
  31. template < typename DerivedV, typename DerivedF, typename Scalar >
  32. IGL_INLINE bool igl::heat_geodesics_precompute(
  33. const Eigen::MatrixBase<DerivedV> & V,
  34. const Eigen::MatrixBase<DerivedF> & F,
  35. const Scalar t,
  36. HeatGeodesicsData<Scalar> & data)
  37. {
  38. typedef Eigen::Matrix<Scalar,Eigen::Dynamic,1> VectorXS;
  39. typedef Eigen::Matrix<Scalar,Eigen::Dynamic,Eigen::Dynamic> MatrixXS;
  40. Eigen::SparseMatrix<Scalar> L,M;
  41. Eigen::Matrix<Scalar,Eigen::Dynamic,3> l_intrinsic;
  42. DerivedF F_intrinsic;
  43. VectorXS dblA;
  44. if(data.use_intrinsic_delaunay)
  45. {
  46. igl::intrinsic_delaunay_cotmatrix(V,F,L,l_intrinsic,F_intrinsic);
  47. igl::massmatrix_intrinsic(l_intrinsic,F_intrinsic,MASSMATRIX_TYPE_DEFAULT,M);
  48. igl::doublearea(l_intrinsic,0,dblA);
  49. igl::grad_intrinsic(l_intrinsic,F_intrinsic,data.Grad);
  50. }else
  51. {
  52. igl::cotmatrix(V,F,L);
  53. igl::massmatrix(V,F,MASSMATRIX_TYPE_DEFAULT,M);
  54. igl::doublearea(V,F,dblA);
  55. igl::grad(V,F,data.Grad);
  56. }
  57. // div
  58. assert(F.cols() == 3 && "Only triangles are supported");
  59. // number of gradient components
  60. data.ng = data.Grad.rows() / F.rows();
  61. assert(data.ng == 3 || data.ng == 2);
  62. data.Div = -0.25*data.Grad.transpose()*dblA.colwise().replicate(data.ng).asDiagonal();
  63. Eigen::SparseMatrix<Scalar> Q = M - t*L;
  64. Eigen::MatrixXi O;
  65. igl::boundary_facets(F,O);
  66. igl::unique(O,data.b);
  67. {
  68. Eigen::SparseMatrix<Scalar> _;
  69. if(!igl::min_quad_with_fixed_precompute(
  70. Q,Eigen::VectorXi(),_,true,data.Neumann))
  71. {
  72. return false;
  73. }
  74. // Only need if there's a boundary
  75. if(data.b.size()>0)
  76. {
  77. if(!igl::min_quad_with_fixed_precompute(Q,data.b,_,true,data.Dirichlet))
  78. {
  79. return false;
  80. }
  81. }
  82. const Eigen::SparseMatrix<double> Aeq = M.diagonal().transpose().sparseView();
  83. L *= -0.5;
  84. if(!igl::min_quad_with_fixed_precompute(
  85. L,Eigen::VectorXi(),Aeq,true,data.Poisson))
  86. {
  87. return false;
  88. }
  89. }
  90. return true;
  91. }
  92. template < typename Scalar, typename Derivedgamma, typename DerivedD>
  93. IGL_INLINE void igl::heat_geodesics_solve(
  94. const HeatGeodesicsData<Scalar> & data,
  95. const Eigen::MatrixBase<Derivedgamma> & gamma,
  96. Eigen::PlainObjectBase<DerivedD> & D)
  97. {
  98. // number of mesh vertices
  99. const int n = data.Grad.cols();
  100. // Set up delta at gamma
  101. DerivedD u0 = DerivedD::Zero(n,1);
  102. for(int g = 0;g<gamma.size();g++)
  103. {
  104. u0(gamma(g)) = 1;
  105. }
  106. // Neumann solution
  107. DerivedD u;
  108. igl::min_quad_with_fixed_solve(
  109. data.Neumann,u0,DerivedD(),DerivedD(),u);
  110. if(data.b.size()>0)
  111. {
  112. // Average Dirichelt and Neumann solutions
  113. DerivedD uD;
  114. igl::min_quad_with_fixed_solve(
  115. data.Dirichlet,u0,DerivedD::Zero(data.b.size()).eval(),DerivedD(),uD);
  116. u += uD;
  117. u *= 0.5;
  118. }
  119. DerivedD grad_u = data.Grad*u;
  120. const int m = data.Grad.rows()/data.ng;
  121. for(int i = 0;i<m;i++)
  122. {
  123. Scalar norm = 0;
  124. for(int d = 0;d<data.ng;d++)
  125. {
  126. norm += grad_u(d*m+i)*grad_u(d*m+i);
  127. }
  128. norm = sqrt(norm);
  129. if(norm == 0)
  130. {
  131. for(int d = 0;d<data.ng;d++) { grad_u(d*m+i) = 0; }
  132. }else
  133. {
  134. for(int d = 0;d<data.ng;d++) { grad_u(d*m+i) /= norm; }
  135. }
  136. }
  137. const DerivedD div_X = -data.Div*grad_u;
  138. const DerivedD Beq = (DerivedD(1,1)<<0).finished();
  139. igl::min_quad_with_fixed_solve(data.Poisson,(-2.0*div_X).eval(),DerivedD(),Beq,D);
  140. DerivedD Dgamma;
  141. igl::slice(D,gamma,Dgamma);
  142. D.array() -= Dgamma.mean();
  143. if(D.mean() < 0)
  144. {
  145. D = -D;
  146. }
  147. }
  148. #ifdef IGL_STATIC_LIBRARY
  149. // Explicit template instantiation
  150. template void igl::heat_geodesics_solve<double, Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >(igl::HeatGeodesicsData<double> const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  151. template bool igl::heat_geodesics_precompute<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&, double, igl::HeatGeodesicsData<double>&);
  152. template bool igl::heat_geodesics_precompute<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&, igl::HeatGeodesicsData<double>&);
  153. #endif