active_set.cpp 13 KB

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