bbw.cpp 3.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114
  1. #define VERBOSE
  2. #include "bbw.h"
  3. #include <igl/cotmatrix.h>
  4. #include <igl/massmatrix.h>
  5. #include <igl/invert_diag.h>
  6. #include <igl/speye.h>
  7. #include <igl/slice_into.h>
  8. #include <igl/normalize_row_sums.h>
  9. #include <igl/verbose.h>
  10. #define EIGEN_YES_I_KNOW_SPARSE_MODULE_IS_NOT_STABLE_YET
  11. #include <Eigen/Sparse>
  12. #include <iostream>
  13. #include <cstdio>
  14. igl::BBWData::BBWData():
  15. partition_unity(false)
  16. {}
  17. void igl::BBWData::print()
  18. {
  19. using namespace std;
  20. cout<<"partition_unity: "<<partition_unity<<endl;
  21. cout<<"W0=["<<endl<<W0<<endl<<"];"<<endl;
  22. }
  23. template <
  24. typename DerivedV,
  25. typename DerivedEle,
  26. typename Derivedb,
  27. typename Derivedbc,
  28. typename DerivedW>
  29. IGL_INLINE bool igl::bbw(
  30. const Eigen::MatrixBase<DerivedV> & V,
  31. const Eigen::MatrixBase<DerivedEle> & Ele,
  32. const Eigen::MatrixBase<Derivedb> & b,
  33. const Eigen::MatrixBase<Derivedbc> & bc,
  34. igl::BBWData & data,
  35. Eigen::MatrixBase<DerivedW> & W
  36. )
  37. {
  38. using namespace igl;
  39. using namespace std;
  40. using namespace Eigen;
  41. // number of domain vertices
  42. int n = V.rows();
  43. // number of handles
  44. int m = bc.cols();
  45. SparseMatrix<typename DerivedW::Scalar> L;
  46. cotmatrix(V,Ele,L);
  47. MassMatrixType mmtype = MASSMATRIX_VORONOI;
  48. if(Ele.cols() == 4)
  49. {
  50. mmtype = MASSMATRIX_BARYCENTRIC;
  51. }
  52. SparseMatrix<typename DerivedW::Scalar> M;
  53. SparseMatrix<typename DerivedW::Scalar> Mi;
  54. massmatrix(V,Ele,mmtype,M);
  55. invert_diag(M,Mi);
  56. // Biharmonic operator
  57. SparseMatrix<typename DerivedW::Scalar> Q = L.transpose() * Mi * L;
  58. W.derived().resize(n,m);
  59. if(data.partition_unity)
  60. {
  61. // Not yet implemented
  62. assert(false);
  63. }else
  64. {
  65. // No linear terms
  66. VectorXd c = VectorXd::Zero(n);
  67. // No linear constraints
  68. SparseMatrix<typename DerivedW::Scalar> A(0,n);
  69. VectorXd uc(0,1);
  70. VectorXd lc(0,1);
  71. // Upper and lower box constraints (Constant bounds)
  72. VectorXd ux = VectorXd::Ones(n);
  73. VectorXd lx = VectorXd::Zero(n);
  74. // Loop over handles
  75. for(int i = 0;i<m;i++)
  76. {
  77. verbose("\n^%s: Computing weight for handle %d out of %d.\n\n",
  78. __FUNCTION__,i+1,m);
  79. // impose boundary conditions
  80. VectorXd bci = bc.col(i);
  81. slice_into(bci,b,ux);
  82. slice_into(bci,b,lx);
  83. VectorXd Wi;
  84. bool r = igl::mosek_quadprog(Q,c,0,A,lc,uc,lx,ux,data.mosek_data,Wi);
  85. if(!r)
  86. {
  87. return false;
  88. }
  89. W.col(i) = Wi;
  90. }
  91. // Need to normalize
  92. igl::normalize_row_sums(W,W);
  93. }
  94. return true;
  95. }
  96. #ifndef IGL_HEADER_ONLY
  97. // Explicit template specialization
  98. template bool igl::bbw<Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<int, -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::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, igl::BBWData&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  99. #endif