bone_heat.cpp 3.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2013 Alec Jacobson <[email protected]>
  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 "bone_heat.h"
  9. #include "EmbreeIntersector.h"
  10. #include "bone_visible.h"
  11. #include "../project_to_line_segment.h"
  12. #include "../cotmatrix.h"
  13. #include "../massmatrix.h"
  14. #include "../min.h"
  15. #include <Eigen/Sparse>
  16. bool igl::embree::bone_heat(
  17. const Eigen::MatrixXd & V,
  18. const Eigen::MatrixXi & F,
  19. const Eigen::MatrixXd & C,
  20. const Eigen::VectorXi & P,
  21. const Eigen::MatrixXi & BE,
  22. const Eigen::MatrixXi & CE,
  23. Eigen::MatrixXd & W)
  24. {
  25. assert(CE.rows() == 0 && "Cage edges not supported.");
  26. assert(C.cols() == V.cols() && "V and C should have same #cols");
  27. assert(BE.cols() == 2 && "BE should have #cols=2");
  28. assert(F.cols() == 3 && "F should contain triangles.");
  29. assert(V.cols() == 3 && "V should contain 3D positions.");
  30. const int n = V.rows();
  31. const int np = P.rows();
  32. const int nb = BE.rows();
  33. const int m = np + nb;
  34. // "double sided lighting"
  35. Eigen::MatrixXi FF;
  36. FF.resize(F.rows()*2,F.cols());
  37. FF << F, F.rowwise().reverse();
  38. // Initialize intersector
  39. EmbreeIntersector ei;
  40. ei.init(V.cast<float>(),F.cast<int>());
  41. typedef Eigen::Matrix<bool ,Eigen::Dynamic,1> VectorXb;
  42. typedef Eigen::Matrix<bool ,Eigen::Dynamic ,Eigen::Dynamic> MatrixXb;
  43. MatrixXb vis_mask(n,m);
  44. // Distances
  45. Eigen::MatrixXd D(n,m);
  46. // loop over points
  47. for(int j = 0;j<np;j++)
  48. {
  49. const Eigen::Vector3d p = C.row(P(j));
  50. D.col(j) = (V.rowwise()-p.transpose()).rowwise().norm();
  51. VectorXb vj;
  52. bone_visible(V,F,ei,p,p,vj);
  53. vis_mask.col(j) = vj;
  54. }
  55. // loop over bones
  56. for(int j = 0;j<nb;j++)
  57. {
  58. const Eigen::Vector3d s = C.row(BE(j,0));
  59. const Eigen::Vector3d d = C.row(BE(j,1));
  60. Eigen::VectorXd t,sqrD;
  61. project_to_line_segment(V,s,d,t,sqrD);
  62. D.col(np+j) = sqrD.array().sqrt();
  63. VectorXb vj;
  64. bone_visible(V,F,ei,s,d,vj);
  65. vis_mask.col(np+j) = vj;
  66. }
  67. assert(CE.rows() == 0 && "Cage edges not supported.");
  68. Eigen::MatrixXd PP = Eigen::MatrixXd::Zero(n,m);
  69. Eigen::VectorXd min_D;
  70. Eigen::VectorXd Hdiag = Eigen::VectorXd::Zero(n);
  71. Eigen::VectorXi J;
  72. igl::min(D,2,min_D,J);
  73. for(int i = 0;i<n;i++)
  74. {
  75. PP(i,J(i)) = 1;
  76. if(vis_mask(i,J(i)))
  77. {
  78. double hii = pow(min_D(i),-2.);
  79. Hdiag(i) = (hii>1e10?1e10:hii);
  80. }
  81. }
  82. Eigen::SparseMatrix<double> Q,L,M;
  83. cotmatrix(V,F,L);
  84. massmatrix(V,F,MASSMATRIX_TYPE_DEFAULT,M);
  85. const auto & H = Hdiag.asDiagonal();
  86. Q = (-L+M*H);
  87. Eigen::SimplicialLLT <Eigen::SparseMatrix<double > > llt;
  88. llt.compute(Q);
  89. switch(llt.info())
  90. {
  91. case Eigen::Success:
  92. break;
  93. case Eigen::NumericalIssue:
  94. #ifdef IGL_BONE_HEAT_DEBUG
  95. std::cerr<<"Error: Numerical issue."<<std::endl;
  96. #endif
  97. return false;
  98. default:
  99. #ifdef IGL_BONE_HEAT_DEBUG
  100. std::cerr<<"Error: Other."<<std::endl;
  101. #endif
  102. return false;
  103. }
  104. const auto & rhs = M*H*PP;
  105. W = llt.solve(rhs);
  106. return true;
  107. }