active_set.cpp 7.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233
  1. #include "active_set.h"
  2. #include "min_quad_with_fixed.h"
  3. #include "slice.h"
  4. #include "matlab_format.h"
  5. #include <iostream>
  6. #include <limits>
  7. #include <algorithm>
  8. template <
  9. typename AT,
  10. typename DerivedB,
  11. typename Derivedknown,
  12. typename DerivedY,
  13. typename AeqT,
  14. typename DerivedBeq,
  15. typename AieqT,
  16. typename DerivedBieq,
  17. typename Derivedlx,
  18. typename Derivedux,
  19. typename DerivedZ
  20. >
  21. IGL_INLINE igl::SolverStatus igl::active_set(
  22. const Eigen::SparseMatrix<AT>& A,
  23. const Eigen::PlainObjectBase<DerivedB> & B,
  24. const Eigen::PlainObjectBase<Derivedknown> & known,
  25. const Eigen::PlainObjectBase<DerivedY> & Y,
  26. const Eigen::SparseMatrix<AeqT>& Aeq,
  27. const Eigen::PlainObjectBase<DerivedBeq> & Beq,
  28. const Eigen::SparseMatrix<AieqT>& Aieq,
  29. const Eigen::PlainObjectBase<DerivedBieq> & Bieq,
  30. const Eigen::PlainObjectBase<Derivedlx> & lx,
  31. const Eigen::PlainObjectBase<Derivedux> & ux,
  32. const igl::active_set_params & params,
  33. Eigen::PlainObjectBase<DerivedZ> & Z
  34. )
  35. {
  36. using namespace igl;
  37. using namespace Eigen;
  38. using namespace std;
  39. // Linear inequality constraints are not supported yet
  40. assert(Aieq.size() == 0);
  41. assert(Bieq.size() == 0);
  42. SolverStatus ret = SOLVER_STATUS_ERROR;
  43. const int n = A.rows();
  44. assert(n == A.cols());
  45. // Discard const qualifiers
  46. //if(B.size() == 0)
  47. //{
  48. // B = Eigen::PlainObjectBase<DerivedB>::Zero(n,1);
  49. //}
  50. assert(n == B.rows());
  51. assert(B.cols() == 1);
  52. assert(Y.cols() == 1);
  53. assert((Aeq.size() == 0 && Beq.size() == 0) || Aeq.cols() == n);
  54. assert((Aeq.size() == 0 && Beq.size() == 0) || Aeq.rows() == Beq.rows());
  55. assert((Aeq.size() == 0 && Beq.size() == 0) || Beq.cols() == 1);
  56. assert((Aieq.size() == 0 && Bieq.size() == 0) || Aieq.cols() == n);
  57. assert((Aieq.size() == 0 && Bieq.size() == 0) || Aieq.rows() == Bieq.rows());
  58. assert((Aieq.size() == 0 && Bieq.size() == 0) || Bieq.cols() == 1);
  59. // Discard const qualifiers
  60. //if(lx.size() == 0)
  61. //{
  62. // lx = Eigen::PlainObjectBase<Derivedlx>::Constant(
  63. // n,1,numeric_limits<typename Derivedlx::Scalar>::min());
  64. //}
  65. //if(ux.size() == 0)
  66. //{
  67. // ux = Eigen::PlainObjectBase<Derivedux>::Constant(
  68. // n,1,numeric_limits<typename Derivedux::Scalar>::max());
  69. //}
  70. assert(lx.rows() == n);
  71. assert(ux.rows() == n);
  72. assert(ux.cols() == 1);
  73. assert(lx.cols() == 1);
  74. assert((ux.array()-lx.array()).minCoeff() > 0);
  75. if(Z.size() != 0)
  76. {
  77. // Initial guess should have correct size
  78. assert(Z.rows() == n);
  79. assert(Z.cols() == 1);
  80. }
  81. assert(known.cols() == 1);
  82. // Number of knowns
  83. const int nk = known.size();
  84. // Initialize active sets
  85. typedef bool BOOL;
  86. #define TRUE true
  87. #define FALSE false
  88. Matrix<BOOL,Dynamic,1> as_lx = Matrix<BOOL,Dynamic,1>::Constant(n,1,FALSE);
  89. Matrix<BOOL,Dynamic,1> as_ux = Matrix<BOOL,Dynamic,1>::Constant(n,1,FALSE);
  90. Matrix<BOOL,Dynamic,1> as_ieq(Aieq.rows(),1);
  91. // Keep track of previous Z for comparison
  92. Eigen::PlainObjectBase<DerivedZ> old_Z;
  93. old_Z = PlainObjectBase<DerivedZ>::Constant(
  94. n,1,numeric_limits<typename DerivedZ::Scalar>::max());
  95. int iter = 0;
  96. while(true)
  97. {
  98. // FIND BREACHES OF CONSTRAINTS
  99. int new_as_lx = 0;
  100. int new_as_ux = 0;
  101. if(Z.size() > 0)
  102. {
  103. for(int z = 0;z < n;z++)
  104. {
  105. if(Z(z) < lx(z))
  106. {
  107. new_as_lx += (as_lx(z)?0:1);
  108. //new_as_lx++;
  109. as_lx(z) = TRUE;
  110. }
  111. if(Z(z) > ux(z))
  112. {
  113. new_as_ux += (as_ux(z)?0:1);
  114. //new_as_ux++;
  115. as_ux(z) = TRUE;
  116. }
  117. }
  118. //cout<<"new_as_lx: "<<new_as_lx<<endl;
  119. //cout<<"new_as_ux: "<<new_as_ux<<endl;
  120. const double diff = (Z-old_Z).squaredNorm();
  121. //cout<<"diff: "<<diff<<endl;
  122. if(diff < params.solution_diff_threshold)
  123. {
  124. ret = SOLVER_STATUS_CONVERGED;
  125. break;
  126. }
  127. old_Z = Z;
  128. }
  129. const int as_lx_count = count(as_lx.data(),as_lx.data()+n,TRUE);
  130. const int as_ux_count = count(as_ux.data(),as_ux.data()+n,TRUE);
  131. // PREPARE FIXED VALUES
  132. Eigen::PlainObjectBase<Derivedknown> known_i;
  133. known_i.resize(nk + as_lx_count + as_ux_count,1);
  134. Eigen::PlainObjectBase<DerivedY> Y_i;
  135. Y_i.resize(nk + as_lx_count + as_ux_count,1);
  136. {
  137. known_i.block(0,0,known.rows(),known.cols()) = known;
  138. Y_i.block(0,0,Y.rows(),Y.cols()) = Y;
  139. int k = nk;
  140. // Then all lx
  141. for(int z = 0;z < n;z++)
  142. {
  143. if(as_lx(z))
  144. {
  145. known_i(k) = z;
  146. Y_i(k) = lx(z);
  147. k++;
  148. }
  149. }
  150. // Finally all ux
  151. for(int z = 0;z < n;z++)
  152. {
  153. if(as_ux(z))
  154. {
  155. known_i(k) = z;
  156. Y_i(k) = ux(z);
  157. k++;
  158. }
  159. }
  160. assert(k==Y_i.size());
  161. assert(k==known_i.size());
  162. }
  163. //cout<<matlab_format((known_i.array()+1).eval(),"known_i")<<endl;
  164. min_quad_with_fixed_data<AT> data;
  165. if(!min_quad_with_fixed_precompute(A,known_i,Aeq,params.Auu_pd,data))
  166. {
  167. cerr<<"Error: min_quad_with_fixed precomputation failed."<<endl;
  168. ret = SOLVER_STATUS_ERROR;
  169. break;
  170. }
  171. if(!min_quad_with_fixed_solve(data,B,Y_i,Beq,Z))
  172. {
  173. cerr<<"Error: min_quad_with_fixed solve failed."<<endl;
  174. ret = SOLVER_STATUS_ERROR;
  175. break;
  176. }
  177. // Compute Lagrange multiplier values for known_i
  178. // This needs to be adjusted slightly if A is not symmetric
  179. assert(data.Auu_sym);
  180. SparseMatrix<AT> Ak;
  181. // Slow
  182. slice(A,known_i,1,Ak);
  183. Eigen::PlainObjectBase<DerivedB> Bk;
  184. slice(B,known_i,Bk);
  185. MatrixXd Lambda_known_i = -(Ak*Z + 0.5*Bk);
  186. // reverse the lambda values for lx
  187. Lambda_known_i.block(nk,0,as_lx_count,1) =
  188. (-1*Lambda_known_i.block(nk,0,as_lx_count,1)).eval();
  189. // Remove from active set
  190. for(int l = 0;l<as_lx_count;l++)
  191. {
  192. if(Lambda_known_i(nk + l) < params.inactive_threshold)
  193. {
  194. as_lx(known_i(nk + l)) = FALSE;
  195. }
  196. }
  197. for(int u = 0;u<as_ux_count;u++)
  198. {
  199. if(Lambda_known_i(nk + as_lx_count + u) <
  200. params.inactive_threshold)
  201. {
  202. as_ux(known_i(nk + as_lx_count + u)) = FALSE;
  203. }
  204. }
  205. iter++;
  206. //cout<<iter<<endl;
  207. if(params.max_iter>0 && iter>=params.max_iter)
  208. {
  209. ret = SOLVER_STATUS_MAX_ITER;
  210. break;
  211. }
  212. }
  213. return ret;
  214. }
  215. #ifndef IGL_HEADER_ONLY
  216. // Explicit template specialization
  217. template igl::SolverStatus igl::active_set<double, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, double, Eigen::Matrix<double, -1, 1, 0, -1, 1>, double, 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::Matrix<double, -1, 1, 0, -1, 1> >(Eigen::SparseMatrix<double, 0, int> const&, 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::SparseMatrix<double, 0, int> const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::SparseMatrix<double, 0, int> const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -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> > const&, igl::active_set_params const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  218. #endif