linprog.cpp 7.6 KB

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