active_set.cpp 13 KB

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