active_set.cpp 9.4 KB

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