linprog.cpp 8.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2015 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 "linprog.h"
  9. #include "find.h"
  10. #include "colon.h"
  11. #include "placeholders.h"
  12. //#define IGL_LINPROG_VERBOSE
  13. IGL_INLINE bool igl::linprog(
  14. const Eigen::VectorXd & c,
  15. const Eigen::MatrixXd & _A,
  16. const Eigen::VectorXd & b,
  17. const int k,
  18. Eigen::VectorXd & x)
  19. {
  20. // This is a very literal translation of
  21. // http://www.mathworks.com/matlabcentral/fileexchange/2166-introduction-to-linear-algebra/content/strang/linprog.m
  22. bool success = true;
  23. // number of constraints
  24. const int m = _A.rows();
  25. // number of original variables
  26. const int n = _A.cols();
  27. // number of iterations
  28. int it = 0;
  29. // maximum number of iterations
  30. //const int MAXIT = 10*m;
  31. const int MAXIT = 100*m;
  32. // residual tolerance
  33. const double tol = 1e-10;
  34. const auto & sign = [](const Eigen::VectorXd & B) -> Eigen::VectorXd
  35. {
  36. Eigen::VectorXd Bsign(B.size());
  37. for(int i = 0;i<B.size();i++)
  38. {
  39. Bsign(i) = B(i)>0?1:(B(i)<0?-1:0);
  40. }
  41. return Bsign;
  42. };
  43. // initial (inverse) basis matrix
  44. Eigen::VectorXd Dv = sign(sign(b).array()+0.5);
  45. Dv.head(k).setConstant(1.);
  46. Eigen::MatrixXd D = Dv.asDiagonal();
  47. // Incorporate slack variables
  48. Eigen::MatrixXd A(_A.rows(),_A.cols()+D.cols());
  49. A<<_A,D;
  50. // Initial basis
  51. Eigen::VectorXi B = igl::colon<int>(n,n+m-1);
  52. // non-basis, may turn out that vector<> would be better here
  53. Eigen::VectorXi N = igl::colon<int>(0,n-1);
  54. int j;
  55. double bmin = b.minCoeff(&j);
  56. int phase;
  57. Eigen::VectorXd xb;
  58. Eigen::VectorXd s;
  59. Eigen::VectorXi J;
  60. if(k>0 && bmin<0)
  61. {
  62. phase = 1;
  63. xb = Eigen::VectorXd::Ones(m);
  64. // super cost
  65. s.resize(n+m+1);
  66. s<<Eigen::VectorXd::Zero(n+k),Eigen::VectorXd::Ones(m-k+1);
  67. N.resize(n+1);
  68. N<<igl::colon<int>(0,n-1),B(j);
  69. J.resize(B.size()-1);
  70. // [0 1 2 3 4]
  71. // ^
  72. // [0 1]
  73. // [3 4]
  74. J.head(j) = B.head(j);
  75. J.tail(B.size()-j-1) = B.tail(B.size()-j-1);
  76. B(j) = n+m;
  77. Eigen::MatrixXd AJ = A(igl::placeholders::all,J);
  78. const Eigen::VectorXd a = b - AJ.rowwise().sum();
  79. {
  80. Eigen::MatrixXd old_A = A;
  81. A.resize(A.rows(),A.cols()+a.cols());
  82. A<<old_A,a;
  83. }
  84. D.col(j) = -a/a(j);
  85. D(j,j) = 1./a(j);
  86. }else if(k==m)
  87. {
  88. phase = 2;
  89. xb = b;
  90. s.resize(c.size()+m);
  91. // cost function
  92. s<<c,Eigen::VectorXd::Zero(m);
  93. }else //k = 0 or bmin >=0
  94. {
  95. phase = 1;
  96. xb = b.array().abs();
  97. s.resize(n+m);
  98. // super cost
  99. s<<Eigen::VectorXd::Zero(n+k),Eigen::VectorXd::Ones(m-k);
  100. }
  101. while(phase<3)
  102. {
  103. double df = -1;
  104. int t = std::numeric_limits<int>::max();
  105. // Lagrange mutipliers fro Ax=b
  106. Eigen::VectorXd yb = D.transpose() * s(B);
  107. while(true)
  108. {
  109. if(MAXIT>0 && it>=MAXIT)
  110. {
  111. #ifdef IGL_LINPROG_VERBOSE
  112. std::cerr<<"linprog: warning! maximum iterations without convergence."<<std::endl;
  113. #endif
  114. success = false;
  115. break;
  116. }
  117. // no freedom for minimization
  118. if(N.size() == 0)
  119. {
  120. break;
  121. }
  122. // reduced costs
  123. Eigen::VectorXd sN = s(N);
  124. Eigen::MatrixXd AN = A(igl::placeholders::all,N);
  125. Eigen::VectorXd r = sN - AN.transpose() * yb;
  126. int q;
  127. // determine new basic variable
  128. double rmin = r.minCoeff(&q);
  129. // optimal! infinity norm
  130. if(rmin>=-tol*(sN.array().abs().maxCoeff()+1))
  131. {
  132. break;
  133. }
  134. // increment iteration count
  135. it++;
  136. // apply Bland's rule to avoid cycling
  137. if(df>=0)
  138. {
  139. if(MAXIT == -1)
  140. {
  141. #ifdef IGL_LINPROG_VERBOSE
  142. std::cerr<<"linprog: warning! degenerate vertex"<<std::endl;
  143. #endif
  144. success = false;
  145. }
  146. igl::find((r.array()<0).eval(),J);
  147. double Nq = N(J).minCoeff();
  148. // again seems like q is assumed to be a scalar though matlab code
  149. // could produce a vector for multiple matches
  150. (N.array()==Nq).cast<int>().maxCoeff(&q);
  151. }
  152. Eigen::VectorXd d = D*A.col(N(q));
  153. Eigen::VectorXi I;
  154. igl::find((d.array()>tol).eval(),I);
  155. if(I.size() == 0)
  156. {
  157. #ifdef IGL_LINPROG_VERBOSE
  158. std::cerr<<"linprog: warning! solution is unbounded"<<std::endl;
  159. #endif
  160. // This seems dubious:
  161. it=-it;
  162. success = false;
  163. break;
  164. }
  165. Eigen::VectorXd xbd = xb(I).array()/d(I).array();
  166. // new use of r
  167. int p;
  168. {
  169. double r;
  170. r = xbd.minCoeff(&p);
  171. p = I(p);
  172. // apply Bland's rule to avoid cycling
  173. if(df>=0)
  174. {
  175. igl::find((xbd.array()==r).eval(),J);
  176. double Bp = B(I(J)).minCoeff();
  177. // idiotic way of finding index in B of Bp
  178. // code down the line seems to assume p is a scalar though the matlab
  179. // code could find a vector of matches)
  180. (B.array()==Bp).cast<int>().maxCoeff(&p);
  181. }
  182. // update x
  183. xb -= r*d;
  184. xb(p) = r;
  185. // change in f
  186. df = r*rmin;
  187. }
  188. // row vector
  189. Eigen::RowVectorXd v = D.row(p)/d(p);
  190. yb += v.transpose() * (s(N(q)) - d.transpose()*s(B));
  191. d(p)-=1;
  192. // update inverse basis matrix
  193. D = D - d*v;
  194. t = B(p);
  195. B(p) = N(q);
  196. if(t>(n+k-1))
  197. {
  198. // remove qth entry from N
  199. Eigen::VectorXi old_N = N;
  200. N.resize(N.size()-1);
  201. N.head(q) = old_N.head(q);
  202. N.head(q) = old_N.head(q);
  203. N.tail(old_N.size()-q-1) = old_N.tail(old_N.size()-q-1);
  204. }else
  205. {
  206. N(q) = t;
  207. }
  208. }
  209. // iterative refinement
  210. xb = (xb+D*(b-A(igl::placeholders::all,B)*xb)).eval();
  211. // must be due to rounding
  212. Eigen::VectorXi I;
  213. igl::find((xb.array()<0).eval(),I);
  214. if(I.size()>0)
  215. {
  216. // so correct
  217. xb(I) = Eigen::VectorXd::Zero(I.size(),1);
  218. }
  219. // B, xb,n,m,res=A(:,B)*xb-b
  220. if(phase == 2 || it<0)
  221. {
  222. break;
  223. }
  224. if(xb.transpose()*s(B) > tol)
  225. {
  226. it = -it;
  227. #ifdef IGL_LINPROG_VERBOSE
  228. std::cerr<<"linprog: warning, no feasible solution"<<std::endl;
  229. #endif
  230. success = false;
  231. break;
  232. }
  233. // re-initialize for Phase 2
  234. phase = phase+1;
  235. s*=1e6*c.array().abs().maxCoeff();
  236. s.head(n) = c;
  237. }
  238. x.setZero(std::max(B.maxCoeff()+1,n));
  239. x(B) = xb;
  240. x = x.head(n).eval();
  241. return success;
  242. }
  243. IGL_INLINE bool igl::linprog(
  244. const Eigen::VectorXd & f,
  245. const Eigen::MatrixXd & A,
  246. const Eigen::VectorXd & b,
  247. const Eigen::MatrixXd & B,
  248. const Eigen::VectorXd & c,
  249. Eigen::VectorXd & x)
  250. {
  251. const int m = A.rows();
  252. const int n = A.cols();
  253. const int p = B.rows();
  254. Eigen::MatrixXd Im = Eigen::MatrixXd::Identity(m,m);
  255. Eigen::MatrixXd AS(m,n+m);
  256. AS<<A,Im;
  257. Eigen::MatrixXd bS = b.array().abs();
  258. for(int i = 0;i<m;i++)
  259. {
  260. const auto & sign = [](double x)->double
  261. {
  262. return (x<0?-1:(x>0?1:0));
  263. };
  264. AS.row(i) *= sign(b(i));
  265. }
  266. Eigen::MatrixXd In = Eigen::MatrixXd::Identity(n,n);
  267. Eigen::MatrixXd P(n+m,2*n+m);
  268. P<< In, -In, Eigen::MatrixXd::Zero(n,m),
  269. Eigen::MatrixXd::Zero(m,2*n), Im;
  270. Eigen::MatrixXd ASP = AS*P;
  271. Eigen::MatrixXd BSP(0,2*n+m);
  272. if(p>0)
  273. {
  274. // B ∈ ℝ^(p × n)
  275. Eigen::MatrixXd BS(p,n+m);
  276. BS<<B,Eigen::MatrixXd::Zero(p,m);
  277. // BS ∈ ℝ^(p × n+m)
  278. BSP = BS*P;
  279. // BSP ∈ ℝ^(p × 2n+m)
  280. }
  281. Eigen::VectorXd fSP = Eigen::VectorXd::Ones(2*n+m);
  282. fSP.head(2*n) = P.block(0,0,n,2*n).transpose()*f;
  283. const Eigen::VectorXd & cc = fSP;
  284. Eigen::MatrixXd AA(m+p,2*n+m);
  285. AA<<ASP,BSP;
  286. Eigen::VectorXd bb(m+p);
  287. bb<<bS,c;
  288. Eigen::VectorXd xxs;
  289. // min ccᵀxxs
  290. // s.t. AA xxs ≤ bb
  291. // xxs ≥ 0
  292. //
  293. // x = x⁺ - x⁻
  294. //
  295. // P
  296. // .--^---.
  297. // [I -I 0 [x⁺ = [x
  298. // 0 0 I] x⁻ s]
  299. // s]
  300. // Pᵀ [xᵀ sᵀ] = xxsᵀ
  301. //
  302. // min [fᵀ -fᵀ 𝟙ᵀ] [x⁺;x⁻;s]
  303. // s.t. AA [x⁺;x⁻;s] ≤ b
  304. // s.t. [x⁺;x⁻;s] ≥ 0
  305. bool ret = linprog(cc,AA,bb,0,xxs);
  306. // x = P(1:n,:) xxs
  307. x = P.block(0,0,n,2*n+m)*xxs;
  308. return ret;
  309. }