polyvector_field_poisson_reconstruction.cpp 5.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116
  1. #include <igl/polyvector_field_poisson_reconstruction.h>
  2. #include <igl/grad.h>
  3. #include <igl/doublearea.h>
  4. #include <igl/sparse.h>
  5. #include <igl/repdiag.h>
  6. #include <igl/slice.h>
  7. #include <igl/slice_into.h>
  8. #include <igl/colon.h>
  9. #include <Eigen/Sparse>
  10. template <typename DerivedV, typename DerivedF, typename DerivedSF, typename DerivedS>
  11. IGL_INLINE void igl::polyvector_field_poisson_reconstruction(
  12. const Eigen::PlainObjectBase<DerivedV> &Vcut,
  13. const Eigen::PlainObjectBase<DerivedF> &Fcut,
  14. const Eigen::PlainObjectBase<DerivedS> &sol3D_combed,
  15. Eigen::PlainObjectBase<DerivedSF> &scalars)
  16. {
  17. Eigen::SparseMatrix<typename DerivedV::Scalar> gradMatrix;
  18. igl::grad(Vcut, Fcut, gradMatrix);
  19. Eigen::VectorXd FAreas;
  20. igl::doublearea(Vcut, Fcut, FAreas);
  21. FAreas = FAreas.array() * .5;
  22. int nf = FAreas.rows();
  23. Eigen::SparseMatrix<typename DerivedV::Scalar> M,M1;
  24. Eigen::VectorXi II = igl::colon<int>(0, nf-1);
  25. igl::sparse(II, II, FAreas, M1);
  26. igl::repdiag(M1, 3, M) ;
  27. int half_degree = sol3D_combed.cols()/3;
  28. int numF = Fcut.rows();
  29. scalars.setZero(Vcut.rows(),half_degree);
  30. Eigen::SparseMatrix<typename DerivedV::Scalar> Q = gradMatrix.transpose()* M *gradMatrix;
  31. //fix one point at Ik=fix, value at fixed xk=0
  32. int fix = 0;
  33. Eigen::VectorXi Ik(1);Ik<<fix;
  34. Eigen::VectorXd xk(1);xk<<0;
  35. //unknown indices
  36. Eigen::VectorXi Iu(Vcut.rows()-1,1);
  37. Iu<<igl::colon<int>(0, fix-1), igl::colon<int>(fix+1,Vcut.rows()-1);
  38. Eigen::SparseMatrix<typename DerivedV::Scalar> Quu, Quk;
  39. igl::slice(Q, Iu, Iu, Quu);
  40. igl::slice(Q, Iu, Ik, Quk);
  41. Eigen::SimplicialLDLT<Eigen::SparseMatrix<typename DerivedV::Scalar> > solver;
  42. solver.compute(Quu);
  43. Eigen::VectorXd vec; vec.setZero(3*numF,1);
  44. for (int i =0; i<half_degree; ++i)
  45. {
  46. vec<<sol3D_combed.col(i*3+0),sol3D_combed.col(i*3+1),sol3D_combed.col(i*3+2);
  47. Eigen::VectorXd b = gradMatrix.transpose()* M * vec;
  48. Eigen::VectorXd bu = igl::slice(b, Iu);
  49. Eigen::VectorXd rhs = bu-Quk*xk;
  50. Eigen::VectorXd yu = solver.solve(rhs);
  51. Eigen::VectorXd y(Vcut.rows(),1);
  52. igl::slice_into(yu, Iu, 1, y);y(Ik[0])=xk[0];
  53. scalars.col(i) = y;
  54. }
  55. }
  56. template <typename DerivedV, typename DerivedF, typename DerivedSF, typename DerivedS, typename DerivedE>
  57. IGL_INLINE double igl::polyvector_field_poisson_reconstruction(
  58. const Eigen::PlainObjectBase<DerivedV> &Vcut,
  59. const Eigen::PlainObjectBase<DerivedF> &Fcut,
  60. const Eigen::PlainObjectBase<DerivedS> &sol3D_combed,
  61. Eigen::PlainObjectBase<DerivedSF> &scalars,
  62. Eigen::PlainObjectBase<DerivedS> &sol3D_recon,
  63. Eigen::PlainObjectBase<DerivedE> &max_error )
  64. {
  65. igl::polyvector_field_poisson_reconstruction(Vcut, Fcut, sol3D_combed, scalars);
  66. Eigen::SparseMatrix<typename DerivedV::Scalar> gradMatrix;
  67. igl::grad(Vcut, Fcut, gradMatrix);
  68. int numF = Fcut.rows();
  69. int half_degree = sol3D_combed.cols()/3;
  70. // evaluate gradient of found scalar function
  71. sol3D_recon.setZero(sol3D_combed.rows(),sol3D_combed.cols());
  72. for (int i =0; i<half_degree; ++i)
  73. {
  74. Eigen::VectorXd vec_poisson = gradMatrix*scalars.col(i);
  75. sol3D_recon.col(i*3+0) = vec_poisson.segment(0*numF, numF);
  76. sol3D_recon.col(i*3+1) = vec_poisson.segment(1*numF, numF);
  77. sol3D_recon.col(i*3+2) = vec_poisson.segment(2*numF, numF);
  78. }
  79. max_error.setZero(numF,1);
  80. for (int i =0; i<half_degree; ++i)
  81. {
  82. Eigen::VectorXd diff = (sol3D_recon.block(0, i*3, numF, 3)-sol3D_combed.block(0, i*3, numF, 3)).rowwise().norm();
  83. diff = diff.array() / sol3D_combed.block(0, i*3, numF, 3).rowwise().norm().array();
  84. max_error = max_error.cwiseMax(diff.cast<typename DerivedE::Scalar>());
  85. }
  86. return max_error.mean();
  87. }
  88. #ifdef IGL_STATIC_LIBRARY
  89. // Explicit template specialization
  90. template double igl::polyvector_field_poisson_reconstruction<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<float, -1, 1, 0, -1, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<float, -1, 1, 0, -1, 1> >&);
  91. template double igl::polyvector_field_poisson_reconstruction<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >(Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  92. #endif