heat_geodesics.cpp 5.1 KB

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