bbw.cpp 5.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204
  1. #include "bbw.h"
  2. #include <igl/cotmatrix.h>
  3. #include <igl/massmatrix.h>
  4. #include <igl/invert_diag.h>
  5. #include <igl/speye.h>
  6. #include <igl/slice_into.h>
  7. #include <igl/normalize_row_sums.h>
  8. #include <igl/min_quad_with_fixed.h>
  9. #define EIGEN_YES_I_KNOW_SPARSE_MODULE_IS_NOT_STABLE_YET
  10. #include <Eigen/Sparse>
  11. #include <iostream>
  12. #include <cstdio>
  13. igl::BBWData::BBWData():
  14. partition_unity(false),
  15. W0(),
  16. #ifndef IGL_NO_MOSEK
  17. mosek_data(),
  18. #endif
  19. active_set_params(),
  20. qp_solver(QP_SOLVER_IGL_ACTIVE_SET),
  21. verbosity(0)
  22. {
  23. // We know that the Bilaplacian is positive semi-definite
  24. active_set_params.Auu_pd = true;
  25. }
  26. void igl::BBWData::print()
  27. {
  28. using namespace std;
  29. cout<<"partition_unity: "<<partition_unity<<endl;
  30. cout<<"W0=["<<endl<<W0<<endl<<"];"<<endl;
  31. cout<<"qp_solver: "<<QPSolverNames[qp_solver]<<endl;
  32. }
  33. template <
  34. typename DerivedV,
  35. typename DerivedEle,
  36. typename Derivedb,
  37. typename Derivedbc,
  38. typename DerivedW>
  39. IGL_INLINE bool igl::bbw(
  40. const Eigen::PlainObjectBase<DerivedV> & V,
  41. const Eigen::PlainObjectBase<DerivedEle> & Ele,
  42. const Eigen::PlainObjectBase<Derivedb> & b,
  43. const Eigen::PlainObjectBase<Derivedbc> & bc,
  44. igl::BBWData & data,
  45. Eigen::PlainObjectBase<DerivedW> & W
  46. )
  47. {
  48. using namespace igl;
  49. using namespace std;
  50. using namespace Eigen;
  51. // number of domain vertices
  52. int n = V.rows();
  53. // number of handles
  54. int m = bc.cols();
  55. SparseMatrix<typename DerivedW::Scalar> L;
  56. cotmatrix(V,Ele,L);
  57. MassMatrixType mmtype = MASSMATRIX_VORONOI;
  58. if(Ele.cols() == 4)
  59. {
  60. mmtype = MASSMATRIX_BARYCENTRIC;
  61. }
  62. SparseMatrix<typename DerivedW::Scalar> M;
  63. SparseMatrix<typename DerivedW::Scalar> Mi;
  64. massmatrix(V,Ele,mmtype,M);
  65. invert_diag(M,Mi);
  66. // Biharmonic operator
  67. SparseMatrix<typename DerivedW::Scalar> Q = L.transpose() * Mi * L;
  68. W.derived().resize(n,m);
  69. if(data.partition_unity)
  70. {
  71. // Not yet implemented
  72. assert(false);
  73. }else
  74. {
  75. // No linear terms
  76. VectorXd c = VectorXd::Zero(n);
  77. // No linear constraints
  78. SparseMatrix<typename DerivedW::Scalar> A(0,n),Aeq(0,n),Aieq(0,n);
  79. VectorXd uc(0,1),Beq(0,1),Bieq(0,1),lc(0,1);
  80. // Upper and lower box constraints (Constant bounds)
  81. VectorXd ux = VectorXd::Ones(n);
  82. VectorXd lx = VectorXd::Zero(n);
  83. active_set_params eff_params = data.active_set_params;
  84. switch(data.qp_solver)
  85. {
  86. case QP_SOLVER_IGL_ACTIVE_SET:
  87. {
  88. if(data.verbosity >= 1)
  89. {
  90. cout<<"BBW: Computing initial weights for "<<m<<" handle"<<
  91. (m!=1?"s":"")<<"."<<endl;
  92. }
  93. min_quad_with_fixed_data<typename DerivedW::Scalar > mqwf;
  94. min_quad_with_fixed_precompute(Q,b,Aeq,true,mqwf);
  95. min_quad_with_fixed_solve(mqwf,c,bc,Beq,W);
  96. // decrement
  97. eff_params.max_iter--;
  98. bool error = false;
  99. // Loop over handles
  100. #pragma omp parallel for
  101. for(int i = 0;i<m;i++)
  102. {
  103. // Quicker exit for openmp
  104. if(error)
  105. {
  106. continue;
  107. }
  108. if(data.verbosity >= 1)
  109. {
  110. cout<<"BBW: Computing weight for handle "<<i+1<<" out of "<<m<<
  111. "."<<endl;
  112. }
  113. VectorXd bci = bc.col(i);
  114. VectorXd Wi;
  115. // use initial guess
  116. Wi = W.col(i);
  117. SolverStatus ret = active_set(
  118. Q,c,b,bci,Aeq,Beq,Aieq,Bieq,lx,ux,eff_params,Wi);
  119. switch(ret)
  120. {
  121. case SOLVER_STATUS_CONVERGED:
  122. break;
  123. case SOLVER_STATUS_MAX_ITER:
  124. cerr<<"active_set: max iter without convergence."<<endl;
  125. break;
  126. case SOLVER_STATUS_ERROR:
  127. default:
  128. cerr<<"active_set error."<<endl;
  129. error = true;
  130. }
  131. W.col(i) = Wi;
  132. }
  133. if(error)
  134. {
  135. return false;
  136. }
  137. break;
  138. }
  139. case QP_SOLVER_MOSEK:
  140. {
  141. #ifdef IGL_NO_MOSEK
  142. assert(false && "Use another QPSolver. Recompile without IGL_NO_MOSEK defined.");
  143. cerr<<"Use another QPSolver. Recompile without IGL_NO_MOSEK defined."<<endl;
  144. return false;
  145. #else
  146. // Loop over handles
  147. for(int i = 0;i<m;i++)
  148. {
  149. if(data.verbosity >= 1)
  150. {
  151. cout<<"BBW: Computing weight for handle "<<i+1<<" out of "<<m<<
  152. "."<<endl;
  153. }
  154. VectorXd bci = bc.col(i);
  155. VectorXd Wi;
  156. // impose boundary conditions via bounds
  157. slice_into(bci,b,ux);
  158. slice_into(bci,b,lx);
  159. bool r = mosek_quadprog(Q,c,0,A,lc,uc,lx,ux,data.mosek_data,Wi);
  160. if(!r)
  161. {
  162. return false;
  163. }
  164. W.col(i) = Wi;
  165. }
  166. #endif
  167. break;
  168. }
  169. default:
  170. {
  171. assert(false && "Unknown qp_solver");
  172. return false;
  173. }
  174. }
  175. #ifndef NDEBUG
  176. const double min_rowsum = W.rowwise().sum().array().abs().minCoeff();
  177. if(min_rowsum < 0.1)
  178. {
  179. cerr<<"bbw.cpp: Warning, minimum row sum is very low. Consider more "
  180. "active set iterations or enforcing partition of unity."<<endl;
  181. }
  182. #endif
  183. }
  184. return true;
  185. }
  186. #ifndef IGL_HEADER_ONLY
  187. // Explicit template specialization
  188. 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::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<int, -1, 1, 0, -1, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, igl::BBWData&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  189. #endif