polyvector_field_poisson_reconstruction.cpp 4.0 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697
  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, typename DerivedE>
  11. IGL_INLINE double 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. Eigen::PlainObjectBase<DerivedS> &sol3D_recon,
  17. Eigen::PlainObjectBase<DerivedE> &max_error )
  18. {
  19. Eigen::SparseMatrix<typename DerivedV::Scalar> gradMatrix;
  20. igl::grad(Vcut, Fcut, gradMatrix);
  21. Eigen::VectorXd FAreas;
  22. igl::doublearea(Vcut, Fcut, FAreas);
  23. FAreas = FAreas.array() * .5;
  24. int nf = FAreas.rows();
  25. Eigen::SparseMatrix<typename DerivedV::Scalar> M,M1;
  26. Eigen::VectorXi II = igl::colon<int>(0, nf-1);
  27. igl::sparse(II, II, FAreas, M1);
  28. igl::repdiag(M1, 3, M) ;
  29. int half_degree = sol3D_combed.cols()/3;
  30. sol3D_recon.setZero(sol3D_combed.rows(),sol3D_combed.cols());
  31. int numF = Fcut.rows();
  32. scalars.setZero(Vcut.rows(),half_degree);
  33. Eigen::SparseMatrix<typename DerivedV::Scalar> Q = gradMatrix.transpose()* M *gradMatrix;
  34. //fix one point at Ik=fix, value at fixed xk=0
  35. int fix = 0;
  36. Eigen::VectorXi Ik(1);Ik<<fix;
  37. Eigen::VectorXd xk(1);xk<<0;
  38. //unknown indices
  39. Eigen::VectorXi Iu(Vcut.rows()-1,1);
  40. Iu<<igl::colon<int>(0, fix-1), igl::colon<int>(fix+1,Vcut.rows()-1);
  41. Eigen::SparseMatrix<typename DerivedV::Scalar> Quu, Quk;
  42. igl::slice(Q, Iu, Iu, Quu);
  43. igl::slice(Q, Iu, Ik, Quk);
  44. Eigen::SimplicialLDLT<Eigen::SparseMatrix<typename DerivedV::Scalar> > solver;
  45. solver.compute(Quu);
  46. Eigen::VectorXd vec; vec.setZero(3*numF,1);
  47. for (int i =0; i<half_degree; ++i)
  48. {
  49. vec<<sol3D_combed.col(i*3+0),sol3D_combed.col(i*3+1),sol3D_combed.col(i*3+2);
  50. Eigen::VectorXd b = gradMatrix.transpose()* M * vec;
  51. Eigen::VectorXd bu = igl::slice(b, Iu);
  52. Eigen::VectorXd rhs = bu-Quk*xk;
  53. Eigen::MatrixXd yu = solver.solve(rhs);
  54. Eigen::VectorXi index = i*Eigen::VectorXi::Ones(Iu.rows(),1);
  55. igl::slice_into(yu, Iu, index, scalars);scalars(Ik[0],i)=xk[0];
  56. }
  57. // evaluate gradient of found scalar function
  58. for (int i =0; i<half_degree; ++i)
  59. {
  60. Eigen::VectorXd vec_poisson = gradMatrix*scalars.col(i);
  61. sol3D_recon.col(i*3+0) = vec_poisson.segment(0*numF, numF);
  62. sol3D_recon.col(i*3+1) = vec_poisson.segment(1*numF, numF);
  63. sol3D_recon.col(i*3+2) = vec_poisson.segment(2*numF, numF);
  64. }
  65. max_error.setZero(numF,1);
  66. for (int i =0; i<half_degree; ++i)
  67. {
  68. Eigen::VectorXd diff = (sol3D_recon.block(0, i*3, numF, 3)-sol3D_combed.block(0, i*3, numF, 3)).rowwise().norm();
  69. diff = diff.array() / sol3D_combed.block(0, i*3, numF, 3).rowwise().norm().array();
  70. max_error = max_error.cwiseMax(diff.cast<typename DerivedE::Scalar>());
  71. }
  72. return max_error.mean();
  73. }
  74. #ifdef IGL_STATIC_LIBRARY
  75. // Explicit template specialization
  76. 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> >&);
  77. #endif