bone_heat.cpp 3.0 KB

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