active_set.cpp 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2013 Alec Jacobson <alecjacobson@gmail.com>
  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 "active_set.h"
  9. #include "min_quad_with_fixed.h"
  10. #include "slice.h"
  11. #include "cat.h"
  12. #include "matlab_format.h"
  13. #include <iostream>
  14. #include <limits>
  15. #include <algorithm>
  16. template <
  17. typename AT,
  18. typename DerivedB,
  19. typename Derivedknown,
  20. typename DerivedY,
  21. typename AeqT,
  22. typename DerivedBeq,
  23. typename AieqT,
  24. typename DerivedBieq,
  25. typename Derivedlx,
  26. typename Derivedux,
  27. typename DerivedZ
  28. >
  29. IGL_INLINE igl::SolverStatus igl::active_set(
  30. const Eigen::SparseMatrix<AT>& A,
  31. const Eigen::PlainObjectBase<DerivedB> & B,
  32. const Eigen::PlainObjectBase<Derivedknown> & known,
  33. const Eigen::PlainObjectBase<DerivedY> & Y,
  34. const Eigen::SparseMatrix<AeqT>& Aeq,
  35. const Eigen::PlainObjectBase<DerivedBeq> & Beq,
  36. const Eigen::SparseMatrix<AieqT>& Aieq,
  37. const Eigen::PlainObjectBase<DerivedBieq> & Bieq,
  38. const Eigen::PlainObjectBase<Derivedlx> & p_lx,
  39. const Eigen::PlainObjectBase<Derivedux> & p_ux,
  40. const igl::active_set_params & params,
  41. Eigen::PlainObjectBase<DerivedZ> & Z
  42. )
  43. {
  44. //#define ACTIVE_SET_CPP_DEBUG
  45. using namespace igl;
  46. using namespace Eigen;
  47. using namespace std;
  48. SolverStatus ret = SOLVER_STATUS_ERROR;
  49. const int n = A.rows();
  50. assert(n == A.cols() && "A must be square");
  51. // Discard const qualifiers
  52. //if(B.size() == 0)
  53. //{
  54. // B = Eigen::PlainObjectBase<DerivedB>::Zero(n,1);
  55. //}
  56. assert(n == B.rows() && "B.rows() must match A.rows()");
  57. assert(B.cols() == 1 && "B must be a column vector");
  58. assert(Y.cols() == 1 && "Y must be a column vector");
  59. assert((Aeq.size() == 0 && Beq.size() == 0) || Aeq.cols() == n);
  60. assert((Aeq.size() == 0 && Beq.size() == 0) || Aeq.rows() == Beq.rows());
  61. assert((Aeq.size() == 0 && Beq.size() == 0) || Beq.cols() == 1);
  62. assert((Aieq.size() == 0 && Bieq.size() == 0) || Aieq.cols() == n);
  63. assert((Aieq.size() == 0 && Bieq.size() == 0) || Aieq.rows() == Bieq.rows());
  64. assert((Aieq.size() == 0 && Bieq.size() == 0) || Bieq.cols() == 1);
  65. Eigen::PlainObjectBase<Derivedlx> lx;
  66. Eigen::PlainObjectBase<Derivedux> ux;
  67. if(p_lx.size() == 0)
  68. {
  69. lx = Eigen::PlainObjectBase<Derivedlx>::Constant(
  70. n,1,-numeric_limits<typename Derivedlx::Scalar>::max());
  71. }else
  72. {
  73. lx = p_lx;
  74. }
  75. if(ux.size() == 0)
  76. {
  77. ux = Eigen::PlainObjectBase<Derivedux>::Constant(
  78. n,1,numeric_limits<typename Derivedux::Scalar>::max());
  79. }else
  80. {
  81. ux = p_ux;
  82. }
  83. assert(lx.rows() == n && "lx must have n rows");
  84. assert(ux.rows() == n && "ux must have n rows");
  85. assert(ux.cols() == 1 && "lx must be a column vector");
  86. assert(lx.cols() == 1 && "ux must be a column vector");
  87. assert((ux.array()-lx.array()).minCoeff() > 0 && "ux(i) must be > lx(i)");
  88. if(Z.size() != 0)
  89. {
  90. // Initial guess should have correct size
  91. assert(Z.rows() == n && "Z must have n rows");
  92. assert(Z.cols() == 1 && "Z must be a column vector");
  93. }
  94. assert(known.cols() == 1 && "known must be a column vector");
  95. // Number of knowns
  96. const int nk = known.size();
  97. // Initialize active sets
  98. typedef int BOOL;
  99. #define TRUE 1
  100. #define FALSE 0
  101. Matrix<BOOL,Dynamic,1> as_lx = Matrix<BOOL,Dynamic,1>::Constant(n,1,FALSE);
  102. Matrix<BOOL,Dynamic,1> as_ux = Matrix<BOOL,Dynamic,1>::Constant(n,1,FALSE);
  103. Matrix<BOOL,Dynamic,1> as_ieq = Matrix<BOOL,Dynamic,1>::Constant(Aieq.rows(),1,FALSE);
  104. // Keep track of previous Z for comparison
  105. PlainObjectBase<DerivedZ> old_Z;
  106. old_Z = PlainObjectBase<DerivedZ>::Constant(
  107. n,1,numeric_limits<typename DerivedZ::Scalar>::max());
  108. int iter = 0;
  109. while(true)
  110. {
  111. #ifdef ACTIVE_SET_CPP_DEBUG
  112. cout<<"Iteration: "<<iter<<":"<<endl;
  113. cout<<" pre"<<endl;
  114. #endif
  115. // FIND BREACHES OF CONSTRAINTS
  116. int new_as_lx = 0;
  117. int new_as_ux = 0;
  118. int new_as_ieq = 0;
  119. if(Z.size() > 0)
  120. {
  121. for(int z = 0;z < n;z++)
  122. {
  123. if(Z(z) < lx(z))
  124. {
  125. new_as_lx += (as_lx(z)?0:1);
  126. //new_as_lx++;
  127. as_lx(z) = TRUE;
  128. }
  129. if(Z(z) > ux(z))
  130. {
  131. new_as_ux += (as_ux(z)?0:1);
  132. //new_as_ux++;
  133. as_ux(z) = TRUE;
  134. }
  135. }
  136. if(Aieq.rows() > 0)
  137. {
  138. PlainObjectBase<DerivedZ> AieqZ;
  139. AieqZ = Aieq*Z;
  140. for(int a = 0;a<Aieq.rows();a++)
  141. {
  142. if(AieqZ(a) > Bieq(a))
  143. {
  144. new_as_ieq += (as_ieq(a)?0:1);
  145. as_ieq(a) = TRUE;
  146. }
  147. }
  148. }
  149. #ifdef ACTIVE_SET_CPP_DEBUG
  150. cout<<" new_as_lx: "<<new_as_lx<<endl;
  151. cout<<" new_as_ux: "<<new_as_ux<<endl;
  152. #endif
  153. const double diff = (Z-old_Z).squaredNorm();
  154. #ifdef ACTIVE_SET_CPP_DEBUG
  155. cout<<"diff: "<<diff<<endl;
  156. #endif
  157. if(diff < params.solution_diff_threshold)
  158. {
  159. ret = SOLVER_STATUS_CONVERGED;
  160. break;
  161. }
  162. old_Z = Z;
  163. }
  164. const int as_lx_count = count(as_lx.data(),as_lx.data()+n,TRUE);
  165. const int as_ux_count = count(as_ux.data(),as_ux.data()+n,TRUE);
  166. const int as_ieq_count =
  167. count(as_ieq.data(),as_ieq.data()+as_ieq.size(),TRUE);
  168. #ifndef NDEBUG
  169. {
  170. int count = 0;
  171. for(int a = 0;a<as_ieq.size();a++)
  172. {
  173. if(as_ieq(a))
  174. {
  175. assert(as_ieq(a) == TRUE);
  176. count++;
  177. }
  178. }
  179. assert(as_ieq_count == count);
  180. }
  181. #endif
  182. // PREPARE FIXED VALUES
  183. PlainObjectBase<Derivedknown> known_i;
  184. known_i.resize(nk + as_lx_count + as_ux_count,1);
  185. PlainObjectBase<DerivedY> Y_i;
  186. Y_i.resize(nk + as_lx_count + as_ux_count,1);
  187. {
  188. known_i.block(0,0,known.rows(),known.cols()) = known;
  189. Y_i.block(0,0,Y.rows(),Y.cols()) = Y;
  190. int k = nk;
  191. // Then all lx
  192. for(int z = 0;z < n;z++)
  193. {
  194. if(as_lx(z))
  195. {
  196. known_i(k) = z;
  197. Y_i(k) = lx(z);
  198. k++;
  199. }
  200. }
  201. // Finally all ux
  202. for(int z = 0;z < n;z++)
  203. {
  204. if(as_ux(z))
  205. {
  206. known_i(k) = z;
  207. Y_i(k) = ux(z);
  208. k++;
  209. }
  210. }
  211. assert(k==Y_i.size());
  212. assert(k==known_i.size());
  213. }
  214. //cout<<matlab_format((known_i.array()+1).eval(),"known_i")<<endl;
  215. // PREPARE EQUALITY CONSTRAINTS
  216. VectorXi as_ieq_list(as_ieq_count,1);
  217. // Gather active constraints and resp. rhss
  218. PlainObjectBase<DerivedBeq> Beq_i;
  219. Beq_i.resize(Beq.rows()+as_ieq_count,1);
  220. Beq_i.head(Beq.rows()) = Beq;
  221. {
  222. int k =0;
  223. for(int a=0;a<as_ieq.size();a++)
  224. {
  225. if(as_ieq(a))
  226. {
  227. assert(k<as_ieq_list.size());
  228. as_ieq_list(k)=a;
  229. Beq_i(Beq.rows()+k,0) = Bieq(k,0);
  230. k++;
  231. }
  232. }
  233. assert(k == as_ieq_count);
  234. }
  235. // extract active constraint rows
  236. SparseMatrix<AeqT> Aeq_i,Aieq_i;
  237. slice(Aieq,as_ieq_list,1,Aieq_i);
  238. // Append to equality constraints
  239. cat(1,Aeq,Aieq_i,Aeq_i);
  240. min_quad_with_fixed_data<AT> data;
  241. #ifndef NDEBUG
  242. {
  243. // NO DUPES!
  244. Matrix<BOOL,Dynamic,1> fixed = Matrix<BOOL,Dynamic,1>::Constant(n,1,FALSE);
  245. for(int k = 0;k<known_i.size();k++)
  246. {
  247. assert(!fixed[known_i(k)]);
  248. fixed[known_i(k)] = TRUE;
  249. }
  250. }
  251. #endif
  252. #ifdef ACTIVE_SET_CPP_DEBUG
  253. cout<<" min_quad_with_fixed_precompute"<<endl;
  254. #endif
  255. if(!min_quad_with_fixed_precompute(A,known_i,Aeq_i,params.Auu_pd,data))
  256. {
  257. cerr<<"Error: min_quad_with_fixed precomputation failed."<<endl;
  258. if(iter > 0 && Aeq_i.rows() > Aeq.rows())
  259. {
  260. cerr<<" *Are you sure rows of [Aeq;Aieq] are linearly independent?*"<<
  261. endl;
  262. }
  263. ret = SOLVER_STATUS_ERROR;
  264. break;
  265. }
  266. #ifdef ACTIVE_SET_CPP_DEBUG
  267. cout<<" min_quad_with_fixed_solve"<<endl;
  268. #endif
  269. Eigen::PlainObjectBase<DerivedZ> sol;
  270. if(!min_quad_with_fixed_solve(data,B,Y_i,Beq_i,Z,sol))
  271. {
  272. cerr<<"Error: min_quad_with_fixed solve failed."<<endl;
  273. ret = SOLVER_STATUS_ERROR;
  274. break;
  275. }
  276. //cout<<matlab_format((Aeq*Z-Beq).eval(),"cr")<<endl;
  277. //cout<<matlab_format(Z,"Z")<<endl;
  278. #ifdef ACTIVE_SET_CPP_DEBUG
  279. cout<<" post"<<endl;
  280. #endif
  281. // Compute Lagrange multiplier values for known_i
  282. // This needs to be adjusted slightly if A is not symmetric
  283. assert(data.Auu_sym);
  284. SparseMatrix<AT> Ak;
  285. // Slow
  286. slice(A,known_i,1,Ak);
  287. Eigen::PlainObjectBase<DerivedB> Bk;
  288. slice(B,known_i,Bk);
  289. MatrixXd Lambda_known_i = -(0.5*Ak*Z + 0.5*Bk);
  290. // reverse the lambda values for lx
  291. Lambda_known_i.block(nk,0,as_lx_count,1) =
  292. (-1*Lambda_known_i.block(nk,0,as_lx_count,1)).eval();
  293. // Extract Lagrange multipliers for Aieq_i (always at back of sol)
  294. VectorXd Lambda_Aieq_i(Aieq_i.rows(),1);
  295. for(int l = 0;l<Aieq_i.rows();l++)
  296. {
  297. Lambda_Aieq_i(Aieq_i.rows()-1-l) = sol(sol.rows()-1-l);
  298. }
  299. // Remove from active set
  300. for(int l = 0;l<as_lx_count;l++)
  301. {
  302. if(Lambda_known_i(nk + l) < params.inactive_threshold)
  303. {
  304. as_lx(known_i(nk + l)) = FALSE;
  305. }
  306. }
  307. for(int u = 0;u<as_ux_count;u++)
  308. {
  309. if(Lambda_known_i(nk + as_lx_count + u) <
  310. params.inactive_threshold)
  311. {
  312. as_ux(known_i(nk + as_lx_count + u)) = FALSE;
  313. }
  314. }
  315. for(int a = 0;a<as_ieq_count;a++)
  316. {
  317. if(Lambda_Aieq_i(a) < params.inactive_threshold)
  318. {
  319. as_ieq(as_ieq_list(a)) = FALSE;
  320. }
  321. }
  322. iter++;
  323. //cout<<iter<<endl;
  324. if(params.max_iter>0 && iter>=params.max_iter)
  325. {
  326. ret = SOLVER_STATUS_MAX_ITER;
  327. break;
  328. }
  329. }
  330. return ret;
  331. }
  332. #ifdef IGL_STATIC_LIBRARY
  333. // Explicit template specialization
  334. 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> >&);
  335. #endif