├── CMakeLists.txt ├── example └── sdlp_example.cpp ├── README.md └── include └── sdlp └── sdlp.hpp /CMakeLists.txt: -------------------------------------------------------------------------------- 1 | cmake_minimum_required(VERSION 2.8.3) 2 | 3 | project(sdlp) 4 | 5 | set(CMAKE_CXX_FLAGS "-std=c++11") 6 | set(CMAKE_BUILD_TYPE "Release") 7 | set(CMAKE_CXX_FLAGS_RELEASE "-O3 -Wall -fPIC") 8 | 9 | find_package(Eigen3 REQUIRED) 10 | 11 | include_directories( 12 | include 13 | ${EIGEN3_INCLUDE_DIRS} 14 | ) 15 | 16 | add_executable(${PROJECT_NAME}_example example/sdlp_example.cpp) 17 | -------------------------------------------------------------------------------- /example/sdlp_example.cpp: -------------------------------------------------------------------------------- 1 | #include 2 | 3 | #include "sdlp/sdlp.hpp" 4 | 5 | using namespace std; 6 | using namespace Eigen; 7 | 8 | int main(int argc, char **argv) 9 | { 10 | int m = 2 * 7; 11 | Eigen::Matrix x; // decision variables 12 | Eigen::Matrix c; // objective coefficients 13 | Eigen::Matrix A(m, 7); // constraint matrix 14 | Eigen::VectorXd b(m); // constraint bound 15 | 16 | c << 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0; 17 | A << 1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 18 | 0.0, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 19 | 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0, 20 | 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 21 | 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0, 22 | 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 23 | 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 24 | -1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 25 | 0.0, -1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 26 | 0.0, 0.0, -1.0, 0.0, 0.0, 0.0, 0.0, 27 | 0.0, 0.0, 0.0, -1.0, 0.0, 0.0, 0.0, 28 | 0.0, 0.0, 0.0, 0.0, -1.0, 0.0, 0.0, 29 | 0.0, 0.0, 0.0, 0.0, 0.0, -1.0, 0.0, 30 | 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, -1.0; 31 | b << 6.0, 5.0, 4.0, 3.0, 2.0, 1.0, 0.0, 1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0; 32 | 33 | double minobj = sdlp::linprog<7>(c, A, b, x); 34 | 35 | std::cout << "prob:\n" 36 | << std::endl; 37 | std::cout << " min x1 + ... + x7," << std::endl; 38 | std::cout << " s.t. x1 <= 6, x2 <= 5, ..., x7 <= 0," << std::endl; 39 | std::cout << " x1 >= -1, x2 >= -2, ..., x7 >= -7.\n" 40 | << std::endl; 41 | std::cout << "optimal sol: " << x.transpose() << std::endl; 42 | std::cout << "optimal obj: " << minobj << std::endl; 43 | 44 | return 0; 45 | } -------------------------------------------------------------------------------- /README.md: -------------------------------------------------------------------------------- 1 | # SDLP 2 | 3 | Seidel's LP Algorithm: Linear-Complexity Linear Programming (LP) for Small-Dimensions 4 | 5 | # About 6 | 7 | 1. This solver is super efficient for small-dimensional LP with any constraint number, mostly encountered in computational geometry. It enjoys [__linear complexity about the constraint number__](https://dl.acm.org/doi/10.1145/2422.322418). 8 | 9 | 2. The speed is __at least an order of magnitude__ faster than [GLPK](https://en.wikipedia.org/wiki/GNU_Linear_Programming_Kit) in small-dimensional LP (<10) with a large constraints number (>100). 10 | 11 | 3. This solver is adapted from the [linear-fractional programming (LFP)](https://en.wikipedia.org/wiki/Linear-fractional_programming) from Mike Hohmeyer at UC Berkeley based on [Raimund Seidel's algorithm](https://link.springer.com/article/10.1007/BF02574699). Kernel functions are reorganized. Previously-existed bugs are fixed here. An easy-to-use interface for LP via [Eigen](http://eigen.tuxfamily.org) is also added. 12 | 13 | 4. Only [a header file](https://github.com/ZJU-FAST-Lab/SDLP/blob/main/include/sdlp/sdlp.hpp) is all you need. 14 | 15 | # Interface 16 | 17 | To solve a linear programming: 18 | 19 | min cTx, 20 | s.t. Ax<=b, 21 | 22 | where x and c are d-dimensional vectors, b an m-dimensional vector and A an mxd matrix. It is assumed that d is small (<10) while m can be arbitrary value (1<= m <= 1e+8). 23 | 24 | Only one function is all you need: 25 | 26 | template 27 | double linprog(const Eigen::Matrix &c, 28 | const Eigen::Matrix &A, 29 | const Eigen::Matrix &b, 30 | Eigen::Matrix &x); 31 | 32 | Input: 33 | 34 | c: objective coefficient 35 | A: constraint matrix 36 | b: constraint bound 37 | 38 | Output: 39 | 40 | x: optimal solution if solved 41 | return: finite value if solved 42 | -infinity if unbounded 43 | infinity if infeasible 44 | 45 | # Reference 46 | 47 | 1. Megiddo, N., 1984. Linear programming in linear time when the dimension is fixed. Journal of the ACM (JACM), 31(1), pp.114-127. 48 | 2. Seidel, R., 1991. Small-dimensional linear programming and convex hulls made easy. Discrete & Computational Geometry, 6(3), pp.423-434. 49 | -------------------------------------------------------------------------------- /include/sdlp/sdlp.hpp: -------------------------------------------------------------------------------- 1 | /* 2 | * Copyright (c) 1990 Michael E. Hohmeyer, 3 | * hohmeyer@icemcfd.com 4 | * Permission is granted to modify and re-distribute this code in any manner 5 | * as long as this notice is preserved. All standard disclaimers apply. 6 | * 7 | * R. Seidel's algorithm for solving LPs (linear programs.) 8 | */ 9 | 10 | /* 11 | * Copyright (c) 2021 Zhepei Wang, 12 | * wangzhepei@live.com 13 | * 1. Bug fix in "move_to_front" function that "prev[m]" is illegally accessed 14 | * while "prev" originally has only m ints. It is fixed by allocating a 15 | * "prev" with m + 1 ints. 16 | * 2. Add Eigen interface. 17 | * 3. Resursive template. 18 | * Permission is granted to modify and re-distribute this code in any manner 19 | * as long as this notice is preserved. All standard disclaimers apply. 20 | * 21 | * Ref: Seidel, R. (1991), "Small-dimensional linear programming and convex 22 | * hulls made easy", Discrete & Computational Geometry 6 (1): 423–434, 23 | * doi:10.1007/BF02574699 24 | */ 25 | 26 | #ifndef SDLP_HPP 27 | #define SDLP_HPP 28 | 29 | #include 30 | #include 31 | #include 32 | 33 | namespace sdlp 34 | { 35 | constexpr double eps = 1.0e-12; 36 | 37 | enum 38 | { 39 | /* minimum attained */ 40 | MINIMUM = 0, 41 | /* no feasible region */ 42 | INFEASIBLE, 43 | /* unbounded solution */ 44 | UNBOUNDED, 45 | /* only a vertex in the solution set */ 46 | AMBIGUOUS, 47 | }; 48 | 49 | inline double dot2(const double a[2], 50 | const double b[2]) 51 | { 52 | return a[0] * b[0] + a[1] * b[1]; 53 | } 54 | 55 | inline double cross2(const double a[2], 56 | const double b[2]) 57 | { 58 | return a[0] * b[1] - a[1] * b[0]; 59 | } 60 | 61 | inline bool unit2(const double a[2], 62 | double b[2]) 63 | { 64 | const double mag = std::sqrt(a[0] * a[0] + 65 | a[1] * a[1]); 66 | if (mag < 2.0 * eps) 67 | { 68 | return true; 69 | } 70 | b[0] = a[0] / mag; 71 | b[1] = a[1] / mag; 72 | return false; 73 | } 74 | 75 | /* unitize a d + 1 dimensional point */ 76 | template 77 | inline bool unit(double *a) 78 | { 79 | double mag = 0.0; 80 | for (int i = 0; i <= d; i++) 81 | { 82 | mag += a[i] * a[i]; 83 | } 84 | if (mag < (d + 1) * eps * eps) 85 | { 86 | return true; 87 | } 88 | mag = 1.0 / std::sqrt(mag); 89 | for (int i = 0; i <= d; i++) 90 | { 91 | a[i] *= mag; 92 | } 93 | return false; 94 | } 95 | 96 | /* optimize the unconstrained objective */ 97 | template 98 | inline int lp_no_con(const double *n_vec, 99 | const double *d_vec, 100 | double *opt) 101 | { 102 | double n_dot_d = 0.0; 103 | double d_dot_d = 0.0; 104 | for (int i = 0; i <= d; i++) 105 | { 106 | n_dot_d += n_vec[i] * d_vec[i]; 107 | d_dot_d += d_vec[i] * d_vec[i]; 108 | } 109 | if (d_dot_d < eps * eps) 110 | { 111 | n_dot_d = 0.0; 112 | d_dot_d = 1.0; 113 | } 114 | for (int i = 0; i <= d; i++) 115 | { 116 | opt[i] = -n_vec[i] + 117 | d_vec[i] * n_dot_d / d_dot_d; 118 | } 119 | /* normalize the optimal point */ 120 | if (unit(opt)) 121 | { 122 | opt[d] = 1.0; 123 | return AMBIGUOUS; 124 | } 125 | else 126 | { 127 | return MINIMUM; 128 | } 129 | } 130 | 131 | /* returns the plane index that is in i's place */ 132 | inline int move_to_front(const int i, 133 | int *next, 134 | int *prev) 135 | { 136 | if (i == 0 || i == next[0]) 137 | { 138 | return i; 139 | } 140 | const int previ = prev[i]; 141 | /* remove i from it's current position */ 142 | next[prev[i]] = next[i]; 143 | prev[next[i]] = prev[i]; 144 | /* put i at the front */ 145 | next[i] = next[0]; 146 | prev[i] = 0; 147 | prev[next[i]] = i; 148 | next[0] = i; 149 | return previ; 150 | } 151 | 152 | inline void lp_min_lin_rat(const bool degen, 153 | const double cw_vec[2], 154 | const double ccw_vec[2], 155 | const double n_vec[2], 156 | const double d_vec[2], 157 | double opt[2]) 158 | { 159 | /* linear rational function case */ 160 | const double d_cw = dot2(cw_vec, d_vec); 161 | const double d_ccw = dot2(ccw_vec, d_vec); 162 | const double n_cw = dot2(cw_vec, n_vec); 163 | const double n_ccw = dot2(ccw_vec, n_vec); 164 | if (degen) 165 | { 166 | /* if degenerate simply compare values */ 167 | if (n_cw / d_cw < n_ccw / d_ccw) 168 | { 169 | opt[0] = cw_vec[0]; 170 | opt[1] = cw_vec[1]; 171 | } 172 | else 173 | { 174 | opt[0] = ccw_vec[0]; 175 | opt[1] = ccw_vec[1]; 176 | } 177 | /* check CW/CCW bounds are not near a poles */ 178 | } 179 | else if (std::fabs(d_cw) > 2.0 * eps && 180 | std::fabs(d_ccw) > 2.0 * eps) 181 | { 182 | /* the valid region does not contain a poles */ 183 | if (d_cw * d_ccw > 0.0) 184 | { 185 | /* find which end has the minimum value */ 186 | if (n_cw / d_cw < n_ccw / d_ccw) 187 | { 188 | opt[0] = cw_vec[0]; 189 | opt[1] = cw_vec[1]; 190 | } 191 | else 192 | { 193 | opt[0] = ccw_vec[0]; 194 | opt[1] = ccw_vec[1]; 195 | } 196 | } 197 | else 198 | { 199 | /* the valid region does contain a poles */ 200 | if (d_cw > 0.0) 201 | { 202 | opt[0] = -d_vec[1]; 203 | opt[1] = d_vec[0]; 204 | } 205 | else 206 | { 207 | opt[0] = d_vec[1]; 208 | opt[1] = -d_vec[0]; 209 | } 210 | } 211 | } 212 | else if (std::fabs(d_cw) > 2.0 * eps) 213 | { 214 | /* CCW bound is near a pole */ 215 | if (n_ccw * d_cw > 0.0) 216 | { 217 | /* CCW bound is a positive pole */ 218 | opt[0] = cw_vec[0]; 219 | opt[1] = cw_vec[1]; 220 | } 221 | else 222 | { 223 | /* CCW bound is a negative pole */ 224 | opt[0] = ccw_vec[0]; 225 | opt[1] = ccw_vec[1]; 226 | } 227 | } 228 | else if (std::fabs(d_ccw) > 2.0 * eps) 229 | { 230 | /* CW bound is near a pole */ 231 | if (n_cw * d_ccw > 2.0 * eps) 232 | { 233 | /* CW bound is at a positive pole */ 234 | opt[0] = ccw_vec[0]; 235 | opt[1] = ccw_vec[1]; 236 | } 237 | else 238 | { 239 | /* CW bound is at a negative pole */ 240 | opt[0] = cw_vec[0]; 241 | opt[1] = cw_vec[1]; 242 | } 243 | } 244 | else 245 | { 246 | /* both bounds are near poles */ 247 | if (cross2(d_vec, n_vec) > 0.0) 248 | { 249 | opt[0] = cw_vec[0]; 250 | opt[1] = cw_vec[1]; 251 | } 252 | else 253 | { 254 | opt[0] = ccw_vec[0]; 255 | opt[1] = ccw_vec[1]; 256 | } 257 | } 258 | } 259 | 260 | inline int wedge(const double (*halves)[2], 261 | const int m, 262 | int *next, 263 | int *prev, 264 | double cw_vec[2], 265 | double ccw_vec[2], 266 | bool *degen) 267 | { 268 | int i; 269 | double d_cw, d_ccw; 270 | bool offensive; 271 | 272 | *degen = false; 273 | for (i = 0; i != m; i = next[i]) 274 | { 275 | if (!unit2(halves[i], ccw_vec)) 276 | { 277 | /* CW */ 278 | cw_vec[0] = ccw_vec[1]; 279 | cw_vec[1] = -ccw_vec[0]; 280 | /* CCW */ 281 | ccw_vec[0] = -cw_vec[0]; 282 | ccw_vec[1] = -cw_vec[1]; 283 | break; 284 | } 285 | } 286 | if (i == m) 287 | { 288 | return UNBOUNDED; 289 | } 290 | i = 0; 291 | while (i != m) 292 | { 293 | offensive = false; 294 | d_cw = dot2(cw_vec, halves[i]); 295 | d_ccw = dot2(ccw_vec, halves[i]); 296 | if (d_ccw >= 2.0 * eps) 297 | { 298 | if (d_cw <= -2.0 * eps) 299 | { 300 | cw_vec[0] = halves[i][1]; 301 | cw_vec[1] = -halves[i][0]; 302 | unit2(cw_vec, cw_vec); 303 | offensive = true; 304 | } 305 | } 306 | else if (d_cw >= 2.0 * eps) 307 | { 308 | if (d_ccw <= -2.0 * eps) 309 | { 310 | ccw_vec[0] = -halves[i][1]; 311 | ccw_vec[1] = halves[i][0]; 312 | unit2(ccw_vec, ccw_vec); 313 | offensive = true; 314 | } 315 | } 316 | else if (d_ccw <= -2.0 * eps && 317 | d_cw <= -2.0 * eps) 318 | { 319 | return INFEASIBLE; 320 | } 321 | else if (d_cw <= -2.0 * eps || 322 | d_ccw <= -2.0 * eps || 323 | cross2(cw_vec, halves[i]) < 0.0) 324 | { 325 | /* degenerate */ 326 | if (d_cw <= -2.0 * eps) 327 | { 328 | unit2(ccw_vec, cw_vec); 329 | } 330 | else if (d_ccw <= -2.0 * eps) 331 | { 332 | unit2(cw_vec, ccw_vec); 333 | } 334 | *degen = true; 335 | offensive = true; 336 | } 337 | /* place this offensive plane in second place */ 338 | if (offensive) 339 | { 340 | i = move_to_front(i, next, prev); 341 | } 342 | i = next[i]; 343 | if (*degen) 344 | { 345 | break; 346 | } 347 | } 348 | if (*degen) 349 | { 350 | while (i != m) 351 | { 352 | d_cw = dot2(cw_vec, halves[i]); 353 | d_ccw = dot2(ccw_vec, halves[i]); 354 | if (d_cw < -2.0 * eps) 355 | { 356 | if (d_ccw < -2.0 * eps) 357 | { 358 | return INFEASIBLE; 359 | } 360 | else 361 | { 362 | cw_vec[0] = ccw_vec[0]; 363 | cw_vec[1] = ccw_vec[1]; 364 | } 365 | } 366 | else if (d_ccw < -2.0 * eps) 367 | { 368 | ccw_vec[0] = cw_vec[0]; 369 | ccw_vec[1] = cw_vec[1]; 370 | } 371 | i = next[i]; 372 | } 373 | } 374 | return MINIMUM; 375 | } 376 | 377 | /* return the minimum on the projective line */ 378 | inline int lp_base_case(const double (*halves)[2], /* halves --- half lines */ 379 | const int m, /* m --- terminal marker */ 380 | const double n_vec[2], /* n_vec --- numerator funciton */ 381 | const double d_vec[2], /* d_vec --- denominator function */ 382 | double opt[2], /* opt --- optimum */ 383 | int *next, /* next, prev --- double linked list of indices */ 384 | int *prev) 385 | { 386 | double cw_vec[2], ccw_vec[2]; 387 | bool degen; 388 | int status; 389 | 390 | /* find the feasible region of the line */ 391 | status = wedge(halves, m, next, prev, cw_vec, ccw_vec, °en); 392 | 393 | if (status == INFEASIBLE) 394 | { 395 | return status; 396 | } 397 | /* no non-trivial constraints one the plane: return the unconstrained optimum */ 398 | if (status == UNBOUNDED) 399 | { 400 | return lp_no_con<1>(n_vec, d_vec, opt); 401 | } 402 | 403 | if (std::fabs(cross2(n_vec, d_vec)) < 2.0 * eps * eps) 404 | { 405 | if (dot2(n_vec, n_vec) < 2.0 * eps * eps || 406 | dot2(d_vec, d_vec) > 2.0 * eps * eps) 407 | { 408 | /* numerator is zero or numerator and denominator are linearly dependent */ 409 | opt[0] = cw_vec[0]; 410 | opt[1] = cw_vec[1]; 411 | status = AMBIGUOUS; 412 | } 413 | else 414 | { 415 | /* numerator is non-zero and denominator is zero minimize linear functional on circle */ 416 | if (!degen && 417 | cross2(cw_vec, n_vec) <= 0.0 && 418 | cross2(n_vec, ccw_vec) <= 0.0) 419 | { 420 | /* optimum is in interior of feasible region */ 421 | opt[0] = -n_vec[0]; 422 | opt[1] = -n_vec[1]; 423 | } 424 | else if (dot2(n_vec, cw_vec) > dot2(n_vec, ccw_vec)) 425 | { 426 | /* optimum is at CCW boundary */ 427 | opt[0] = ccw_vec[0]; 428 | opt[1] = ccw_vec[1]; 429 | } 430 | else 431 | { 432 | /* optimum is at CW boundary */ 433 | opt[0] = cw_vec[0]; 434 | opt[1] = cw_vec[1]; 435 | } 436 | status = MINIMUM; 437 | } 438 | } 439 | else 440 | { 441 | /* niether numerator nor denominator is zero */ 442 | lp_min_lin_rat(degen, cw_vec, ccw_vec, n_vec, d_vec, opt); 443 | status = MINIMUM; 444 | } 445 | return status; 446 | } 447 | 448 | /* find the largest coefficient in a plane */ 449 | template 450 | inline void findimax(const double *pln, 451 | int *imax) 452 | { 453 | *imax = 0; 454 | double rmax = std::fabs(pln[0]); 455 | for (int i = 1; i <= d; i++) 456 | { 457 | const double ab = std::fabs(pln[i]); 458 | if (ab > rmax) 459 | { 460 | *imax = i; 461 | rmax = ab; 462 | } 463 | } 464 | } 465 | 466 | template 467 | inline void vector_up(const double *equation, 468 | const int ivar, 469 | const double *low_vector, 470 | double *vector) 471 | { 472 | vector[ivar] = 0.0; 473 | for (int i = 0; i <= d; i++) 474 | { 475 | if (i != ivar) 476 | { 477 | const int j = i < ivar ? i : i - 1; 478 | vector[i] = low_vector[j]; 479 | vector[ivar] -= equation[i] * low_vector[j]; 480 | } 481 | } 482 | vector[ivar] /= equation[ivar]; 483 | } 484 | 485 | template 486 | inline void vector_down(const double *elim_eqn, 487 | const int ivar, 488 | const double *old_vec, 489 | double *new_vec) 490 | { 491 | double ve = 0.0; 492 | double ee = 0.0; 493 | for (int i = 0; i <= d; i++) 494 | { 495 | ve += old_vec[i] * elim_eqn[i]; 496 | ee += elim_eqn[i] * elim_eqn[i]; 497 | } 498 | const double fac = ve / ee; 499 | for (int i = 0; i <= d; i++) 500 | { 501 | if (i != ivar) 502 | { 503 | new_vec[i < ivar ? i : i - 1] = 504 | old_vec[i] - elim_eqn[i] * fac; 505 | } 506 | } 507 | } 508 | 509 | template 510 | inline void plane_down(const double *elim_eqn, 511 | const int ivar, 512 | const double *old_plane, 513 | double *new_plane) 514 | { 515 | const double crit = old_plane[ivar] / elim_eqn[ivar]; 516 | for (int i = 0; i <= d; i++) 517 | { 518 | if (i != ivar) 519 | { 520 | new_plane[i < ivar ? i : i - 1] = 521 | old_plane[i] - elim_eqn[i] * crit; 522 | } 523 | } 524 | } 525 | 526 | template 527 | inline int linfracprog(const double *halves, /* halves --- half spaces */ 528 | const int max_size, /* max_size --- size of halves array */ 529 | const int m, /* m --- terminal marker */ 530 | const double *n_vec, /* n_vec --- numerator vector */ 531 | const double *d_vec, /* d_vec --- denominator vector */ 532 | double *opt, /* opt --- optimum */ 533 | double *work, /* work --- work space (see below) */ 534 | int *next, /* next --- array of indices into halves */ 535 | int *prev) /* prev --- array of indices into halves */ 536 | /* 537 | ** 538 | ** half-spaces are in the form 539 | ** halves[i][0]*x[0] + halves[i][1]*x[1] + 540 | ** ... + halves[i][d-1]*x[d-1] + halves[i][d]*x[d] >= 0 541 | ** 542 | ** coefficients should be normalized 543 | ** half-spaces should be in random order 544 | ** the order of the half spaces is 0, next[0] next[next[0]] ... 545 | ** and prev[next[i]] = i 546 | ** 547 | ** halves: (max_size)x(d+1) 548 | ** 549 | ** the optimum has been computed for the half spaces 550 | ** 0 , next[0], next[next[0]] , ... , prev[0] 551 | ** the next plane that needs to be tested is 0 552 | ** 553 | ** m is the index of the first plane that is NOT on the list 554 | ** i.e. m is the terminal marker for the linked list. 555 | ** 556 | ** the objective function is dot(x,nvec)/dot(x,dvec) 557 | ** if you want the program to solve standard d dimensional linear programming 558 | ** problems then n_vec = ( x0, x1, x2, ..., xd-1, 0) 559 | ** and d_vec = ( 0, 0, 0, ..., 0, 1) 560 | ** and halves[0] = (0, 0, ... , 1) 561 | ** 562 | ** work points to (max_size+3)*(d+2)*(d-1)/2 double space 563 | */ 564 | { 565 | int status, imax; 566 | double *new_opt, *new_n_vec, *new_d_vec, *new_halves, *new_work; 567 | const double *plane_i; 568 | 569 | double val = 0.0; 570 | for (int j = 0; j <= d; j++) 571 | { 572 | val += d_vec[j] * d_vec[j]; 573 | } 574 | const bool d_vec_zero = (val < (d + 1) * eps * eps); 575 | 576 | /* find the unconstrained minimum */ 577 | status = lp_no_con(n_vec, d_vec, opt); 578 | if (m <= 0) 579 | { 580 | return status; 581 | } 582 | 583 | /* allocate memory for next level of recursion */ 584 | new_opt = work; 585 | new_n_vec = new_opt + d; 586 | new_d_vec = new_n_vec + d; 587 | new_halves = new_d_vec + d; 588 | new_work = new_halves + max_size * d; 589 | for (int i = 0; i != m; i = next[i]) 590 | { 591 | /* if the optimum is not in half space i then project the problem onto that plane */ 592 | plane_i = halves + i * (d + 1); 593 | /* determine if the optimum is on the correct side of plane_i */ 594 | val = 0.0; 595 | for (int j = 0; j <= d; j++) 596 | { 597 | val += opt[j] * plane_i[j]; 598 | } 599 | if (val < -(d + 1) * eps) 600 | { 601 | /* find the largest of the coefficients to eliminate */ 602 | findimax(plane_i, &imax); 603 | /* eliminate that variable */ 604 | if (i != 0) 605 | { 606 | const double fac = 1.0 / plane_i[imax]; 607 | for (int j = 0; j != i; j = next[j]) 608 | { 609 | const double *old_plane = halves + j * (d + 1); 610 | const double crit = old_plane[imax] * fac; 611 | double *new_plane = new_halves + j * d; 612 | for (int k = 0; k <= d; k++) 613 | { 614 | const int l = k < imax ? k : k - 1; 615 | new_plane[l] = k != imax ? old_plane[k] - plane_i[k] * crit : new_plane[l]; 616 | } 617 | } 618 | } 619 | /* project the objective function to lower dimension */ 620 | if (d_vec_zero) 621 | { 622 | vector_down(plane_i, imax, n_vec, new_n_vec); 623 | for (int j = 0; j < d; j++) 624 | { 625 | new_d_vec[j] = 0.0; 626 | } 627 | } 628 | else 629 | { 630 | plane_down(plane_i, imax, n_vec, new_n_vec); 631 | plane_down(plane_i, imax, d_vec, new_d_vec); 632 | } 633 | /* solve sub problem */ 634 | status = linfracprog(new_halves, max_size, i, new_n_vec, 635 | new_d_vec, new_opt, new_work, next, prev); 636 | /* back substitution */ 637 | if (status != INFEASIBLE) 638 | { 639 | vector_up(plane_i, imax, new_opt, opt); 640 | 641 | /* inline code for unit */ 642 | double mag = 0.0; 643 | for (int j = 0; j <= d; j++) 644 | { 645 | mag += opt[j] * opt[j]; 646 | } 647 | mag = 1.0 / sqrt(mag); 648 | for (int j = 0; j <= d; j++) 649 | { 650 | opt[j] *= mag; 651 | } 652 | } 653 | else 654 | { 655 | return status; 656 | } 657 | /* place this offensive plane in second place */ 658 | i = move_to_front(i, next, prev); 659 | } 660 | } 661 | return status; 662 | } 663 | 664 | template <> 665 | inline int linfracprog<1>(const double *halves, 666 | const int max_size, 667 | const int m, 668 | const double *n_vec, 669 | const double *d_vec, 670 | double *opt, 671 | double *work, 672 | int *next, 673 | int *prev) 674 | { 675 | if (m > 0) 676 | { 677 | return lp_base_case((const double(*)[2])halves, m, 678 | n_vec, d_vec, opt, next, prev); 679 | } 680 | else 681 | { 682 | return lp_no_con<1>(n_vec, d_vec, opt); 683 | } 684 | } 685 | 686 | inline void rand_permutation(const int n, 687 | int *p) 688 | { 689 | typedef std::uniform_int_distribution rand_int; 690 | typedef rand_int::param_type rand_range; 691 | static std::mt19937_64 gen; 692 | static rand_int rdi(0, 1); 693 | int j, k; 694 | for (int i = 0; i < n; i++) 695 | { 696 | p[i] = i; 697 | } 698 | for (int i = 0; i < n; i++) 699 | { 700 | rdi.param(rand_range(0, n - i - 1)); 701 | j = rdi(gen) + i; 702 | k = p[j]; 703 | p[j] = p[i]; 704 | p[i] = k; 705 | } 706 | } 707 | 708 | template 709 | inline double linprog(const Eigen::Matrix &c, 710 | const Eigen::Matrix &A, 711 | const Eigen::Matrix &b, 712 | Eigen::Matrix &x) 713 | /* 714 | ** min cTx, s.t. Ax<=b 715 | ** dim(x) << dim(b) 716 | */ 717 | { 718 | int m = b.size() + 1; 719 | x.setZero(); 720 | if (m <= 1) 721 | { 722 | return c.cwiseAbs().maxCoeff() > 0.0 ? -INFINITY : 0.0; 723 | } 724 | 725 | Eigen::VectorXi perm(m - 1); 726 | Eigen::VectorXi next(m); 727 | /* original allocated size is m, here changed to m + 1 for legal tail accessing */ 728 | Eigen::VectorXi prev(m + 1); 729 | Eigen::Matrix n_vec; 730 | Eigen::Matrix d_vec; 731 | Eigen::Matrix opt; 732 | Eigen::Matrix halves(d + 1, m); 733 | Eigen::VectorXd work((m + 3) * (d + 2) * (d - 1) / 2); 734 | 735 | halves.col(0).setZero(); 736 | halves(d, 0) = 1.0; 737 | halves.topRightCorner(d, m - 1) = -A.transpose(); 738 | halves.bottomRightCorner(1, m - 1) = b.transpose(); 739 | /* normalize all halves as required in linfracprog */ 740 | halves.colwise().normalize(); 741 | n_vec.head(d) = c; 742 | n_vec(d) = 0.0; 743 | d_vec.setZero(); 744 | d_vec(d) = 1.0; 745 | 746 | /* randomize the input planes */ 747 | rand_permutation(m - 1, perm.data()); 748 | /* previous to 0 is actually never used */ 749 | prev(0) = 0; 750 | /* link the zero position in at the beginning */ 751 | next(0) = perm(0) + 1; 752 | prev(perm(0) + 1) = 0; 753 | /* link the other planes */ 754 | for (int i = 0; i < m - 2; i++) 755 | { 756 | next(perm(i) + 1) = perm(i + 1) + 1; 757 | prev(perm(i + 1) + 1) = perm(i) + 1; 758 | } 759 | /* flag the last plane */ 760 | next(perm(m - 2) + 1) = m; 761 | 762 | int status = sdlp::linfracprog(halves.data(), m, m, 763 | n_vec.data(), d_vec.data(), 764 | opt.data(), work.data(), 765 | next.data(), prev.data()); 766 | 767 | /* handle states for linprog whose definitions differ from linfracprog */ 768 | double minimum = INFINITY; 769 | if (status != sdlp::INFEASIBLE) 770 | { 771 | if (opt(d) != 0.0 && status != sdlp::UNBOUNDED) 772 | { 773 | x = opt.head(d) / opt(d); 774 | minimum = c.dot(x); 775 | } 776 | 777 | if (opt(d) == 0.0 || status == sdlp::UNBOUNDED) 778 | { 779 | x = opt.head(d); 780 | minimum = -INFINITY; 781 | } 782 | } 783 | 784 | return minimum; 785 | } 786 | 787 | } // namespace sdlp 788 | 789 | #endif 790 | --------------------------------------------------------------------------------