gonum.org/v1/gonum@v0.14.0/optimize/convex/lp/simplex.go (about) 1 // Copyright ©2016 The Gonum Authors. All rights reserved. 2 // Use of this source code is governed by a BSD-style 3 // license that can be found in the LICENSE file. 4 5 // package lp implements routines for solving linear programs. 6 package lp 7 8 import ( 9 "errors" 10 "fmt" 11 "math" 12 13 "gonum.org/v1/gonum/floats" 14 "gonum.org/v1/gonum/mat" 15 ) 16 17 // TODO(btracey): Could have a solver structure with an abstract factorizer. With 18 // this transformation the same high-level code could handle both Dense and Sparse. 19 // TODO(btracey): Need to improve error handling. Only want to panic if condition number inf. 20 // TODO(btracey): Performance enhancements. There are currently lots of linear 21 // solves that can be improved by doing rank-one updates. For example, the swap 22 // step is just a rank-one update. 23 // TODO(btracey): Better handling on the linear solve errors. If the condition 24 // number is not inf and the equation solved "well", should keep moving. 25 26 var ( 27 ErrBland = errors.New("lp: bland: all replacements are negative or cause ill-conditioned ab") 28 ErrInfeasible = errors.New("lp: problem is infeasible") 29 ErrLinSolve = errors.New("lp: linear solve failure") 30 ErrUnbounded = errors.New("lp: problem is unbounded") 31 ErrSingular = errors.New("lp: A is singular") 32 ErrZeroColumn = errors.New("lp: A has a column of all zeros") 33 ErrZeroRow = errors.New("lp: A has a row of all zeros") 34 ) 35 36 const badShape = "lp: size mismatch" 37 38 // TODO(btracey): Should these tolerances be part of a settings struct? 39 40 const ( 41 // initPosTol is the tolerance on the initial condition being feasible. Strictly, 42 // the x should be positive, but instead it must be greater than -initPosTol. 43 initPosTol = 1e-13 44 // blandNegTol is the tolerance on the value being greater than 0 in the bland test. 45 blandNegTol = 1e-14 46 // rRoundTol is the tolerance for rounding values to zero when testing if 47 // constraints are met. 48 rRoundTol = 1e-13 49 // dRoundTol is the tolerance for testing if values are zero for the problem 50 // being unbounded. 51 dRoundTol = 1e-13 52 // phaseIZeroTol tests if the Phase I problem returned a feasible solution. 53 phaseIZeroTol = 1e-12 54 // blandZeroTol is the tolerance on testing if the bland solution can move. 55 blandZeroTol = 1e-12 56 ) 57 58 // Simplex solves a linear program in standard form using Danzig's Simplex 59 // algorithm. The standard form of a linear program is: 60 // 61 // minimize cᵀ x 62 // s.t. A*x = b 63 // x >= 0 . 64 // 65 // The input tol sets how close to the optimal solution is found (specifically, 66 // when the maximal reduced cost is below tol). An error will be returned if the 67 // problem is infeasible or unbounded. In rare cases, numeric errors can cause 68 // the Simplex to fail. In this case, an error will be returned along with the 69 // most recently found feasible solution. 70 // 71 // The Convert function can be used to transform a general LP into standard form. 72 // 73 // The input matrix A must have at least as many columns as rows, len(c) must 74 // equal the number of columns of A, and len(b) must equal the number of rows of 75 // A or Simplex will panic. A must also have full row rank and may not contain any 76 // columns with all zeros, or Simplex will return an error. 77 // 78 // initialBasic can be used to set the initial set of indices for a feasible 79 // solution to the LP. If an initial feasible solution is not known, initialBasic 80 // may be nil. If initialBasic is non-nil, len(initialBasic) must equal the number 81 // of rows of A and must be an actual feasible solution to the LP, otherwise 82 // Simplex will panic. 83 // 84 // A description of the Simplex algorithm can be found in Ch. 8 of 85 // 86 // Strang, Gilbert. "Linear Algebra and Applications." Academic, New York (1976). 87 // 88 // For a detailed video introduction, see lectures 11-13 of UC Math 352 89 // 90 // https://www.youtube.com/watch?v=ESzYPFkY3og&index=11&list=PLh464gFUoJWOmBYla3zbZbc4nv2AXez6X. 91 func Simplex(c []float64, A mat.Matrix, b []float64, tol float64, initialBasic []int) (optF float64, optX []float64, err error) { 92 ans, x, _, err := simplex(initialBasic, c, A, b, tol) 93 return ans, x, err 94 } 95 96 func simplex(initialBasic []int, c []float64, A mat.Matrix, b []float64, tol float64) (float64, []float64, []int, error) { 97 err := verifyInputs(initialBasic, c, A, b) 98 if err != nil { 99 if err == ErrUnbounded { 100 return math.Inf(-1), nil, nil, ErrUnbounded 101 } 102 return math.NaN(), nil, nil, err 103 } 104 m, n := A.Dims() 105 106 if m == n { 107 // Problem is exactly constrained, perform a linear solve. 108 bVec := mat.NewVecDense(len(b), b) 109 x := make([]float64, n) 110 xVec := mat.NewVecDense(n, x) 111 err := xVec.SolveVec(A, bVec) 112 if err != nil { 113 return math.NaN(), nil, nil, ErrSingular 114 } 115 for _, v := range x { 116 if v < 0 { 117 return math.NaN(), nil, nil, ErrInfeasible 118 } 119 } 120 f := floats.Dot(x, c) 121 return f, x, nil, nil 122 } 123 124 // There is at least one optimal solution to the LP which is at the intersection 125 // to a set of constraint boundaries. For a standard form LP with m variables 126 // and n equality constraints, at least m-n elements of x must equal zero 127 // at optimality. The Simplex algorithm solves the standard-form LP by starting 128 // at an initial constraint vertex and successively moving to adjacent constraint 129 // vertices. At every vertex, the set of non-zero x values is the "basic 130 // feasible solution". The list of non-zero x's are maintained in basicIdxs, 131 // the respective columns of A are in ab, and the actual non-zero values of 132 // x are in xb. 133 // 134 // The LP is equality constrained such that A * x = b. This can be expanded 135 // to 136 // ab * xb + an * xn = b 137 // where ab are the columns of a in the basic set, and an are all of the 138 // other columns. Since each element of xn is zero by definition, this means 139 // that for all feasible solutions xb = ab^-1 * b. 140 // 141 // Before the simplex algorithm can start, an initial feasible solution must 142 // be found. If initialBasic is non-nil a feasible solution has been supplied. 143 // Otherwise the "Phase I" problem must be solved to find an initial feasible 144 // solution. 145 146 var basicIdxs []int // The indices of the non-zero x values. 147 var ab *mat.Dense // The subset of columns of A listed in basicIdxs. 148 var xb []float64 // The non-zero elements of x. xb = ab^-1 b 149 150 if initialBasic != nil { 151 // InitialBasic supplied. Panic if incorrect length or infeasible. 152 if len(initialBasic) != m { 153 panic("lp: incorrect number of initial vectors") 154 } 155 ab = mat.NewDense(m, len(initialBasic), nil) 156 extractColumns(ab, A, initialBasic) 157 xb = make([]float64, m) 158 err = initializeFromBasic(xb, ab, b) 159 if err != nil { 160 panic(err) 161 } 162 basicIdxs = make([]int, len(initialBasic)) 163 copy(basicIdxs, initialBasic) 164 } else { 165 // No initial basis supplied. Solve the PhaseI problem. 166 basicIdxs, ab, xb, err = findInitialBasic(A, b) 167 if err != nil { 168 return math.NaN(), nil, nil, err 169 } 170 } 171 172 // basicIdxs contains the indexes for an initial feasible solution, 173 // ab contains the extracted columns of A, and xb contains the feasible 174 // solution. All x not in the basic set are 0 by construction. 175 176 // nonBasicIdx is the set of nonbasic variables. 177 nonBasicIdx := make([]int, 0, n-m) 178 inBasic := make(map[int]struct{}) 179 for _, v := range basicIdxs { 180 inBasic[v] = struct{}{} 181 } 182 for i := 0; i < n; i++ { 183 _, ok := inBasic[i] 184 if !ok { 185 nonBasicIdx = append(nonBasicIdx, i) 186 } 187 } 188 189 // cb is the subset of c for the basic variables. an and cn 190 // are the equivalents to ab and cb but for the nonbasic variables. 191 cb := make([]float64, len(basicIdxs)) 192 for i, idx := range basicIdxs { 193 cb[i] = c[idx] 194 } 195 cn := make([]float64, len(nonBasicIdx)) 196 for i, idx := range nonBasicIdx { 197 cn[i] = c[idx] 198 } 199 an := mat.NewDense(m, len(nonBasicIdx), nil) 200 extractColumns(an, A, nonBasicIdx) 201 202 bVec := mat.NewVecDense(len(b), b) 203 cbVec := mat.NewVecDense(len(cb), cb) 204 205 // Temporary data needed each iteration. (Described later) 206 r := make([]float64, n-m) 207 move := make([]float64, m) 208 209 // Solve the linear program starting from the initial feasible set. This is 210 // the "Phase 2" problem. 211 // 212 // Algorithm: 213 // 1) Compute the "reduced costs" for the non-basic variables. The reduced 214 // costs are the lagrange multipliers of the constraints. 215 // r = cn - anᵀ * ab¯ᵀ * cb 216 // 2) If all of the reduced costs are positive, no improvement is possible, 217 // and the solution is optimal (xn can only increase because of 218 // non-negativity constraints). Otherwise, the solution can be improved and 219 // one element will be exchanged in the basic set. 220 // 3) Choose the x_n with the most negative value of r. Call this value xe. 221 // This variable will be swapped into the basic set. 222 // 4) Increase xe until the next constraint boundary is met. This will happen 223 // when the first element in xb becomes 0. The distance xe can increase before 224 // a given element in xb becomes negative can be found from 225 // xb = Ab^-1 b - Ab^-1 An xn 226 // = Ab^-1 b - Ab^-1 Ae xe 227 // = bhat + d x_e 228 // xe = bhat_i / - d_i 229 // where Ae is the column of A corresponding to xe. 230 // The constraining basic index is the first index for which this is true, 231 // so remove the element which is min_i (bhat_i / -d_i), assuming d_i is negative. 232 // If no d_i is less than 0, then the problem is unbounded. 233 // 5) If the new xe is 0 (that is, bhat_i == 0), then this location is at 234 // the intersection of several constraints. Use the Bland rule instead 235 // of the rule in step 4 to avoid cycling. 236 for { 237 // Compute reduced costs -- r = cn - anᵀ ab¯ᵀ cb 238 var tmp mat.VecDense 239 err = tmp.SolveVec(ab.T(), cbVec) 240 if err != nil { 241 break 242 } 243 data := make([]float64, n-m) 244 tmp2 := mat.NewVecDense(n-m, data) 245 tmp2.MulVec(an.T(), &tmp) 246 floats.SubTo(r, cn, data) 247 248 // Replace the most negative element in the simplex. If there are no 249 // negative entries then the optimal solution has been found. 250 minIdx := floats.MinIdx(r) 251 if r[minIdx] >= -tol { 252 break 253 } 254 255 for i, v := range r { 256 if math.Abs(v) < rRoundTol { 257 r[i] = 0 258 } 259 } 260 261 // Compute the moving distance. 262 err = computeMove(move, minIdx, A, ab, xb, nonBasicIdx) 263 if err != nil { 264 if err == ErrUnbounded { 265 return math.Inf(-1), nil, nil, ErrUnbounded 266 } 267 break 268 } 269 270 // Replace the basic index along the tightest constraint. 271 replace := floats.MinIdx(move) 272 if move[replace] <= 0 { 273 replace, minIdx, err = replaceBland(A, ab, xb, basicIdxs, nonBasicIdx, r, move) 274 if err != nil { 275 if err == ErrUnbounded { 276 return math.Inf(-1), nil, nil, ErrUnbounded 277 } 278 break 279 } 280 } 281 282 // Replace the constrained basicIdx with the newIdx. 283 basicIdxs[replace], nonBasicIdx[minIdx] = nonBasicIdx[minIdx], basicIdxs[replace] 284 cb[replace], cn[minIdx] = cn[minIdx], cb[replace] 285 tmpCol1 := mat.Col(nil, replace, ab) 286 tmpCol2 := mat.Col(nil, minIdx, an) 287 ab.SetCol(replace, tmpCol2) 288 an.SetCol(minIdx, tmpCol1) 289 290 // Compute the new xb. 291 xbVec := mat.NewVecDense(len(xb), xb) 292 err = xbVec.SolveVec(ab, bVec) 293 if err != nil { 294 break 295 } 296 } 297 // Found the optimum successfully or died trying. The basic variables get 298 // their values, and the non-basic variables are all zero. 299 opt := floats.Dot(cb, xb) 300 xopt := make([]float64, n) 301 for i, v := range basicIdxs { 302 xopt[v] = xb[i] 303 } 304 return opt, xopt, basicIdxs, err 305 } 306 307 // computeMove computes how far can be moved replacing each index. The results 308 // are stored into move. 309 func computeMove(move []float64, minIdx int, A mat.Matrix, ab *mat.Dense, xb []float64, nonBasicIdx []int) error { 310 // Find ae. 311 col := mat.Col(nil, nonBasicIdx[minIdx], A) 312 aCol := mat.NewVecDense(len(col), col) 313 314 // d = - Ab^-1 Ae 315 nb, _ := ab.Dims() 316 d := make([]float64, nb) 317 dVec := mat.NewVecDense(nb, d) 318 err := dVec.SolveVec(ab, aCol) 319 if err != nil { 320 return ErrLinSolve 321 } 322 floats.Scale(-1, d) 323 324 for i, v := range d { 325 if math.Abs(v) < dRoundTol { 326 d[i] = 0 327 } 328 } 329 330 // If no di < 0, then problem is unbounded. 331 if floats.Min(d) >= 0 { 332 return ErrUnbounded 333 } 334 335 // move = bhat_i / - d_i, assuming d is negative. 336 bHat := xb // ab^-1 b 337 for i, v := range d { 338 if v >= 0 { 339 move[i] = math.Inf(1) 340 } else { 341 move[i] = bHat[i] / math.Abs(v) 342 } 343 } 344 return nil 345 } 346 347 // replaceBland uses the Bland rule to find the indices to swap if the minimum 348 // move is 0. The indices to be swapped are replace and minIdx (following the 349 // nomenclature in the main routine). 350 func replaceBland(A mat.Matrix, ab *mat.Dense, xb []float64, basicIdxs, nonBasicIdx []int, r, move []float64) (replace, minIdx int, err error) { 351 m, _ := A.Dims() 352 // Use the traditional bland rule, except don't replace a constraint which 353 // causes the new ab to be singular. 354 for i, v := range r { 355 if v > -blandNegTol { 356 continue 357 } 358 minIdx = i 359 err = computeMove(move, minIdx, A, ab, xb, nonBasicIdx) 360 if err != nil { 361 // Either unbounded or something went wrong. 362 return -1, -1, err 363 } 364 replace = floats.MinIdx(move) 365 if math.Abs(move[replace]) > blandZeroTol { 366 // Large enough that it shouldn't be a problem 367 return replace, minIdx, nil 368 } 369 // Find a zero index where replacement is non-singular. 370 biCopy := make([]int, len(basicIdxs)) 371 for replace, v := range move { 372 if v > blandZeroTol { 373 continue 374 } 375 copy(biCopy, basicIdxs) 376 biCopy[replace] = nonBasicIdx[minIdx] 377 abTmp := mat.NewDense(m, len(biCopy), nil) 378 extractColumns(abTmp, A, biCopy) 379 // If the condition number is reasonable, use this index. 380 if mat.Cond(abTmp, 1) < 1e16 { 381 return replace, minIdx, nil 382 } 383 } 384 } 385 return -1, -1, ErrBland 386 } 387 388 func verifyInputs(initialBasic []int, c []float64, A mat.Matrix, b []float64) error { 389 m, n := A.Dims() 390 if m > n { 391 panic("lp: more equality constraints than variables") 392 } 393 if len(c) != n { 394 panic("lp: c vector incorrect length") 395 } 396 if len(b) != m { 397 panic("lp: b vector incorrect length") 398 } 399 if len(c) != n { 400 panic("lp: c vector incorrect length") 401 } 402 if len(initialBasic) != 0 && len(initialBasic) != m { 403 panic("lp: initialBasic incorrect length") 404 } 405 406 // Do some sanity checks so that ab does not become singular during the 407 // simplex solution. If the ZeroRow checks are removed then the code for 408 // finding a set of linearly independent columns must be improved. 409 410 // Check that if a row of A only has zero elements that corresponding 411 // element in b is zero, otherwise the problem is infeasible. 412 // Otherwise return ErrZeroRow. 413 for i := 0; i < m; i++ { 414 isZero := true 415 for j := 0; j < n; j++ { 416 if A.At(i, j) != 0 { 417 isZero = false 418 break 419 } 420 } 421 if isZero && b[i] != 0 { 422 // Infeasible 423 return ErrInfeasible 424 } else if isZero { 425 return ErrZeroRow 426 } 427 } 428 // Check that if a column only has zero elements that the respective C vector 429 // is positive (otherwise unbounded). Otherwise return ErrZeroColumn. 430 for j := 0; j < n; j++ { 431 isZero := true 432 for i := 0; i < m; i++ { 433 if A.At(i, j) != 0 { 434 isZero = false 435 break 436 } 437 } 438 if isZero && c[j] < 0 { 439 return ErrUnbounded 440 } else if isZero { 441 return ErrZeroColumn 442 } 443 } 444 return nil 445 } 446 447 // initializeFromBasic initializes the basic feasible solution given a set of 448 // basic indices. It extracts the columns of A specified by basicIdxs and finds 449 // the x values at that location. These are stored into xb. 450 // 451 // If the columns of A are not linearly independent or if the initial set is not 452 // feasible, an error is returned. 453 func initializeFromBasic(xb []float64, ab *mat.Dense, b []float64) error { 454 m, _ := ab.Dims() 455 if len(xb) != m { 456 panic("simplex: bad xb length") 457 } 458 xbMat := mat.NewVecDense(m, xb) 459 460 err := xbMat.SolveVec(ab, mat.NewVecDense(m, b)) 461 if err != nil { 462 return errors.New("lp: subcolumns of A for supplied initial basic singular") 463 } 464 // The solve ensures that the equality constraints are met (ab * xb = b). 465 // Thus, the solution is feasible if and only if all of the x's are positive. 466 allPos := true 467 for _, v := range xb { 468 if v < -initPosTol { 469 allPos = false 470 break 471 } 472 } 473 if !allPos { 474 return errors.New("lp: supplied subcolumns not a feasible solution") 475 } 476 return nil 477 } 478 479 // extractColumns copies the columns specified by cols into the columns of dst. 480 func extractColumns(dst *mat.Dense, A mat.Matrix, cols []int) { 481 r, c := dst.Dims() 482 ra, _ := A.Dims() 483 if ra != r { 484 panic("simplex: row mismatch") 485 } 486 if c != len(cols) { 487 panic("simplex: column mismatch") 488 } 489 col := make([]float64, r) 490 for j, idx := range cols { 491 mat.Col(col, idx, A) 492 dst.SetCol(j, col) 493 } 494 } 495 496 // findInitialBasic finds an initial basic solution, and returns the basic 497 // indices, ab, and xb. 498 func findInitialBasic(A mat.Matrix, b []float64) ([]int, *mat.Dense, []float64, error) { 499 m, n := A.Dims() 500 basicIdxs := findLinearlyIndependent(A) 501 if len(basicIdxs) != m { 502 return nil, nil, nil, ErrSingular 503 } 504 505 // It may be that this linearly independent basis is also a feasible set. If 506 // so, the Phase I problem can be avoided. 507 ab := mat.NewDense(m, len(basicIdxs), nil) 508 extractColumns(ab, A, basicIdxs) 509 xb := make([]float64, m) 510 err := initializeFromBasic(xb, ab, b) 511 if err == nil { 512 return basicIdxs, ab, xb, nil 513 } 514 515 // This set was not feasible. Instead the "Phase I" problem must be solved 516 // to find an initial feasible set of basis. 517 // 518 // Method: Construct an LP whose optimal solution is a feasible solution 519 // to the original LP. 520 // 1) Introduce an artificial variable x_{n+1}. 521 // 2) Let x_j be the most negative element of x_b (largest constraint violation). 522 // 3) Add the artificial variable to A with: 523 // a_{n+1} = b - \sum_{i in basicIdxs} a_i + a_j 524 // swap j with n+1 in the basicIdxs. 525 // 4) Define a new LP: 526 // minimize x_{n+1} 527 // subject to [A A_{n+1}][x_1 ... x_{n+1}] = b 528 // x, x_{n+1} >= 0 529 // 5) Solve this LP. If x_{n+1} != 0, then the problem is infeasible, otherwise 530 // the found basis can be used as an initial basis for phase II. 531 // 532 // The extra column in Step 3 is defined such that the vector of 1s is an 533 // initial feasible solution. 534 535 // Find the largest constraint violator. 536 // Compute a_{n+1} = b - \sum{i in basicIdxs}a_i + a_j. j is in basicIDx, so 537 // instead just subtract the basicIdx columns that are not minIDx. 538 minIdx := floats.MinIdx(xb) 539 aX1 := make([]float64, m) 540 copy(aX1, b) 541 col := make([]float64, m) 542 for i, v := range basicIdxs { 543 if i == minIdx { 544 continue 545 } 546 mat.Col(col, v, A) 547 floats.Sub(aX1, col) 548 } 549 550 // Construct the new LP. 551 // aNew = [A, a_{n+1}] 552 // bNew = b 553 // cNew = 1 for x_{n+1} 554 aNew := mat.NewDense(m, n+1, nil) 555 aNew.Copy(A) 556 aNew.SetCol(n, aX1) 557 basicIdxs[minIdx] = n // swap minIdx with n in the basic set. 558 c := make([]float64, n+1) 559 c[n] = 1 560 561 // Solve the Phase I linear program. 562 _, xOpt, newBasic, err := simplex(basicIdxs, c, aNew, b, 1e-10) 563 if err != nil { 564 return nil, nil, nil, fmt.Errorf("lp: error finding feasible basis: %s", err) 565 } 566 567 // The original LP is infeasible if the added variable has non-zero value 568 // in the optimal solution to the Phase I problem. 569 if math.Abs(xOpt[n]) > phaseIZeroTol { 570 return nil, nil, nil, ErrInfeasible 571 } 572 573 // The basis found in Phase I is a feasible solution to the original LP if 574 // the added variable is not in the basis. 575 addedIdx := -1 576 for i, v := range newBasic { 577 if v == n { 578 addedIdx = i 579 } 580 xb[i] = xOpt[v] 581 } 582 if addedIdx == -1 { 583 extractColumns(ab, A, newBasic) 584 return newBasic, ab, xb, nil 585 } 586 587 // The value of the added variable is in the basis, but it has a zero value. 588 // See if exchanging another variable into the basic set finds a feasible 589 // solution. 590 basicMap := make(map[int]struct{}) 591 for _, v := range newBasic { 592 basicMap[v] = struct{}{} 593 } 594 var set bool 595 for i := range xOpt { 596 if _, inBasic := basicMap[i]; inBasic { 597 continue 598 } 599 newBasic[addedIdx] = i 600 if set { 601 mat.Col(col, i, A) 602 ab.SetCol(addedIdx, col) 603 } else { 604 extractColumns(ab, A, newBasic) 605 set = true 606 } 607 err := initializeFromBasic(xb, ab, b) 608 if err == nil { 609 return newBasic, ab, xb, nil 610 } 611 } 612 return nil, nil, nil, ErrInfeasible 613 } 614 615 // findLinearlyIndependent finds a set of linearly independent columns of A, and 616 // returns the column indexes of the linearly independent columns. 617 func findLinearlyIndependent(A mat.Matrix) []int { 618 m, n := A.Dims() 619 idxs := make([]int, 0, m) 620 columns := mat.NewDense(m, m, nil) 621 newCol := make([]float64, m) 622 // Walk in reverse order because slack variables are typically the last columns 623 // of A. 624 for i := n - 1; i >= 0; i-- { 625 if len(idxs) == m { 626 break 627 } 628 mat.Col(newCol, i, A) 629 columns.SetCol(len(idxs), newCol) 630 if len(idxs) == 0 { 631 // A column is linearly independent from the null set. 632 // If all-zero column of A are allowed, this code needs to be adjusted. 633 idxs = append(idxs, i) 634 continue 635 } 636 if mat.Cond(columns.Slice(0, m, 0, len(idxs)+1), 1) > 1e12 { 637 // Not linearly independent. 638 continue 639 } 640 idxs = append(idxs, i) 641 } 642 return idxs 643 }