active_set.cpp 13 KB

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