bbw.cpp 2.9 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2016 Alec Jacobson <[email protected]>
  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 "bbw.h"
  9. #include "mosek_quadprog.h"
  10. #include "../harmonic.h"
  11. #include <Eigen/Sparse>
  12. #include <iostream>
  13. #include <cstdio>
  14. template <
  15. typename DerivedV,
  16. typename DerivedEle,
  17. typename Derivedb,
  18. typename Derivedbc,
  19. typename DerivedW>
  20. IGL_INLINE bool igl::mosek::bbw(
  21. const Eigen::MatrixBase<DerivedV> & V,
  22. const Eigen::MatrixBase<DerivedEle> & Ele,
  23. const Eigen::MatrixBase<Derivedb> & b,
  24. const Eigen::MatrixBase<Derivedbc> & bc,
  25. igl::BBWData & data,
  26. igl::mosek::MosekData & mosek_data,
  27. Eigen::PlainObjectBase<DerivedW> & W
  28. )
  29. {
  30. assert(!data.partition_unity && "partition_unity not implemented yet");
  31. // number of domain vertices
  32. int n = V.rows();
  33. // number of handles
  34. int m = bc.cols();
  35. // Build biharmonic operator
  36. Eigen::SparseMatrix<typename DerivedV::Scalar> Q;
  37. harmonic(V,Ele,2,Q);
  38. W.derived().resize(n,m);
  39. // No linear terms
  40. Eigen::VectorXd c = Eigen::VectorXd::Zero(n);
  41. // No linear constraints
  42. Eigen::SparseMatrix<typename DerivedW::Scalar> A(0,n);
  43. Eigen::VectorXd uc(0,1),lc(0,1);
  44. // Upper and lower box constraints (Constant bounds)
  45. Eigen::VectorXd ux = Eigen::VectorXd::Ones(n);
  46. Eigen::VectorXd lx = Eigen::VectorXd::Zero(n);
  47. // Loop over handles
  48. for(int i = 0;i<m;i++)
  49. {
  50. if(data.verbosity >= 1)
  51. {
  52. std::cout<<"BBW: Computing weight for handle "<<i+1<<" out of "<<m<<
  53. "."<<std::endl;
  54. }
  55. Eigen::VectorXd bci = bc.col(i);
  56. Eigen::VectorXd Wi;
  57. // impose boundary conditions via bounds
  58. ux(b) = bci;
  59. lx(b) = bci;
  60. bool r = mosek_quadprog(Q,c,0,A,lc,uc,lx,ux,mosek_data,Wi);
  61. if(!r)
  62. {
  63. return false;
  64. }
  65. W.col(i) = Wi;
  66. }
  67. #ifndef NDEBUG
  68. const double min_rowsum = W.rowwise().sum().array().abs().minCoeff();
  69. if(min_rowsum < 0.1)
  70. {
  71. std::cerr<<"bbw.cpp: Warning, minimum row sum is very low. Consider more "
  72. "active set iterations or enforcing partition of unity."<<std::endl;
  73. }
  74. #endif
  75. return true;
  76. }
  77. #ifdef IGL_STATIC_LIBRARY
  78. // Explicit template instantiation
  79. template bool igl::mosek::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&, igl::mosek::MosekData&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  80. #endif