bbw.cpp 4.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145
  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 "min_quad_with_fixed.h"
  10. #include "harmonic.h"
  11. #include "parallel_for.h"
  12. #include <Eigen/Sparse>
  13. #include <iostream>
  14. #include <mutex>
  15. #include <cstdio>
  16. igl::BBWData::BBWData():
  17. partition_unity(false),
  18. W0(),
  19. active_set_params(),
  20. verbosity(0)
  21. {
  22. // We know that the Bilaplacian is positive semi-definite
  23. active_set_params.Auu_pd = true;
  24. }
  25. void igl::BBWData::print()
  26. {
  27. std::cout<<"partition_unity: "<<partition_unity<<std::endl;
  28. std::cout<<"W0=["<<std::endl<<W0<<std::endl<<"];"<<std::endl;
  29. }
  30. template <
  31. typename DerivedV,
  32. typename DerivedEle,
  33. typename Derivedb,
  34. typename Derivedbc,
  35. typename DerivedW>
  36. IGL_INLINE bool igl::bbw(
  37. const Eigen::MatrixBase<DerivedV> & V,
  38. const Eigen::MatrixBase<DerivedEle> & Ele,
  39. const Eigen::MatrixBase<Derivedb> & b,
  40. const Eigen::MatrixBase<Derivedbc> & bc,
  41. igl::BBWData & data,
  42. Eigen::PlainObjectBase<DerivedW> & W
  43. )
  44. {
  45. assert(!data.partition_unity && "partition_unity not implemented yet");
  46. // number of domain vertices
  47. int n = V.rows();
  48. // number of handles
  49. int m = bc.cols();
  50. // Build biharmonic operator
  51. Eigen::SparseMatrix<typename DerivedV::Scalar> Q;
  52. harmonic(V,Ele,2,Q);
  53. W.derived().resize(n,m);
  54. // No linear terms
  55. Eigen::VectorXd c = Eigen::VectorXd::Zero(n);
  56. // No linear constraints
  57. Eigen::SparseMatrix<typename DerivedW::Scalar> A(0,n),Aeq(0,n),Aieq(0,n);
  58. Eigen::VectorXd Beq(0,1),Bieq(0,1);
  59. // Upper and lower box constraints (Constant bounds)
  60. Eigen::VectorXd ux = Eigen::VectorXd::Ones(n);
  61. Eigen::VectorXd lx = Eigen::VectorXd::Zero(n);
  62. active_set_params eff_params = data.active_set_params;
  63. if(data.verbosity >= 1)
  64. {
  65. std::cout<<"BBW: max_iter: "<<data.active_set_params.max_iter<<std::endl;
  66. std::cout<<"BBW: eff_max_iter: "<<eff_params.max_iter<<std::endl;
  67. }
  68. if(data.verbosity >= 1)
  69. {
  70. std::cout<<"BBW: Computing initial weights for "<<m<<" handle"<<
  71. (m!=1?"s":"")<<"."<<std::endl;
  72. }
  73. min_quad_with_fixed_data<typename DerivedW::Scalar > mqwf;
  74. min_quad_with_fixed_precompute(Q,b,Aeq,true,mqwf);
  75. min_quad_with_fixed_solve(mqwf,c,bc,Beq,W);
  76. // decrement
  77. eff_params.max_iter--;
  78. bool error = false;
  79. // Loop over handles
  80. std::mutex critical;
  81. const auto & optimize_weight = [&](const int i)
  82. {
  83. // Quicker exit for paralle_for
  84. if(error)
  85. {
  86. return;
  87. }
  88. if(data.verbosity >= 1)
  89. {
  90. std::lock_guard<std::mutex> lock(critical);
  91. std::cout<<"BBW: Computing weight for handle "<<i+1<<" out of "<<m<<
  92. "."<<std::endl;
  93. }
  94. Eigen::VectorXd bci = bc.col(i);
  95. Eigen::VectorXd Wi;
  96. // use initial guess
  97. Wi = W.col(i);
  98. SolverStatus ret = active_set(
  99. Q,c,b,bci,Aeq,Beq,Aieq,Bieq,lx,ux,eff_params,Wi);
  100. switch(ret)
  101. {
  102. case SOLVER_STATUS_CONVERGED:
  103. break;
  104. case SOLVER_STATUS_MAX_ITER:
  105. #ifdef IGL_BBW_DEBUG
  106. cerr<<"active_set: max iter without convergence."<<endl;
  107. #endif
  108. break;
  109. case SOLVER_STATUS_ERROR:
  110. default:
  111. #ifdef IGL_BBW_DEBUG
  112. cerr<<"active_set error."<<endl;
  113. #endif
  114. error = true;
  115. }
  116. W.col(i) = Wi;
  117. };
  118. parallel_for(m,optimize_weight,2);
  119. if(error)
  120. {
  121. return false;
  122. }
  123. #ifndef NDEBUG
  124. const double min_rowsum = W.rowwise().sum().array().abs().minCoeff();
  125. if(min_rowsum < 0.1)
  126. {
  127. std::cerr<<"bbw.cpp: Warning, minimum row sum is very low. Consider more "
  128. "active set iterations or enforcing partition of unity."<<std::endl;
  129. }
  130. #endif
  131. return true;
  132. }
  133. #ifdef IGL_STATIC_LIBRARY
  134. // Explicit template instantiation
  135. 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::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1>>&);
  136. #endif