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