polyvector_field_poisson_reconstruction.cpp 6.0 KB

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