Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

1

2

3

4

5

6

7

8

9

10

11

12

13

14

15

16

17

18

19

20

21

22

23

24

25

26

27

28

29

30

31

32

33

34

35

36

37

38

39

40

41

42

43

44

45

46

47

48

49

50

51

52

53

54

55

56

57

58

59

60

61

62

63

64

65

66

67

68

69

70

71

72

73

74

75

76

77

78

79

80

81

82

83

84

85

86

87

88

89

90

91

92

93

94

95

96

97

98

99

100

101

102

103

104

105

106

107

108

109

110

111

112

113

114

115

116

117

118

119

120

121

122

123

124

125

126

127

128

129

130

131

132

133

134

135

136

137

138

139

140

141

142

143

144

145

146

147

148

149

150

151

152

153

154

155

156

157

158

159

160

161

162

163

164

165

166

167

168

169

170

171

172

173

174

175

176

177

178

179

180

181

182

183

184

185

186

187

188

189

190

191

192

193

194

195

196

197

198

199

200

201

202

203

204

205

206

207

208

209

210

211

212

213

214

215

216

217

218

219

220

221

222

223

224

225

226

227

228

229

230

231

232

233

234

235

236

237

238

239

240

241

242

243

244

245

246

247

248

249

250

251

252

253

254

255

256

257

258

259

260

261

262

263

264

265

266

267

268

269

270

271

272

273

274

275

276

277

278

279

280

281

282

283

284

285

286

287

288

289

290

291

292

293

294

295

296

297

298

299

300

301

302

303

304

305

306

307

308

309

310

311

312

313

314

315

316

317

318

319

320

321

322

323

324

325

326

327

328

329

330

331

332

333

334

335

336

337

338

339

340

341

342

343

344

345

346

347

348

349

350

351

352

353

354

355

356

357

358

359

360

361

362

363

364

365

366

367

368

369

370

371

372

373

374

375

376

377

378

379

380

381

382

383

384

385

386

387

388

389

390

391

392

393

394

395

396

397

398

399

400

401

402

403

404

405

406

407

408

409

410

411

412

413

414

415

416

417

418

419

420

421

422

423

424

425

426

427

428

429

430

431

432

433

434

435

436

437

438

439

440

441

442

443

444

445

446

447

448

449

450

451

452

453

454

455

456

457

458

459

460

461

462

463

464

465

466

467

468

469

470

471

472

473

474

475

476

477

478

479

480

481

482

483

484

485

486

487

488

489

490

491

492

493

494

495

496

497

498

499

500

501

502

503

504

505

506

507

508

509

510

511

512

513

514

515

516

517

518

519

520

521

522

523

524

525

526

527

528

529

530

531

532

533

534

535

536

537

538

539

540

541

542

543

544

545

546

547

548

549

550

551

552

553

554

555

556

557

558

559

560

561

562

563

564

565

566

567

568

569

570

571

572

573

574

575

576

577

578

579

580

581

582

583

584

585

586

587

588

589

590

591

592

593

594

595

596

597

598

599

600

601

602

603

604

605

606

607

608

609

610

611

612

613

614

615

616

617

618

619

620

621

622

623

624

625

626

627

628

629

630

631

632

633

634

635

636

637

638

639

640

641

642

643

644

645

646

647

648

649

650

651

652

653

654

655

656

657

658

659

660

661

662

663

664

665

666

667

668

669

670

671

672

673

674

675

676

677

678

679

680

681

682

683

684

685

686

687

688

689

690

691

692

693

694

695

696

697

698

699

700

701

702

703

704

705

706

707

708

709

710

711

712

713

714

715

716

717

718

719

720

721

722

723

724

725

726

727

728

729

730

731

732

733

734

735

736

737

738

739

740

741

742

743

744

745

746

747

748

749

750

751

752

753

754

755

756

757

758

759

760

761

762

763

764

765

766

767

768

769

770

771

772

773

774

775

776

777

778

779

780

781

782

783

784

785

786

787

788

789

790

791

792

793

794

795

796

797

798

799

800

801

802

803

804

805

806

807

808

809

810

811

812

813

814

815

816

817

818

819

820

821

822

823

824

825

826

827

828

829

830

831

832

833

834

835

836

837

838

839

840

841

842

843

844

845

846

847

848

849

850

851

852

853

854

855

856

857

858

859

860

861

862

863

864

865

866

867

868

869

870

871

872

873

874

875

876

877

878

879

880

881

882

883

884

885

886

887

888

889

890

891

892

893

894

895

896

897

898

899

900

901

902

903

904

905

906

907

908

909

910

911

912

913

914

915

916

917

918

919

920

921

922

923

924

925

926

927

928

929

930

931

932

933

934

935

936

937

938

939

940

941

942

943

944

945

946

947

948

949

950

951

952

953

954

955

956

957

958

959

960

961

962

963

964

965

966

967

968

969

970

971

972

973

974

975

976

977

978

979

980

981

982

983

984

985

986

987

988

989

990

991

992

993

994

995

996

997

998

999

1000

1001

1002

1003

1004

1005

1006

1007

1008

1009

1010

1011

1012

1013

1014

1015

1016

1017

1018

1019

1020

1021

1022

1023

1024

1025

1026

1027

1028

1029

1030

1031

1032

1033

1034

1035

1036

1037

1038

1039

1040

1041

1042

1043

1044

1045

1046

1047

1048

1049

1050

1051

1052

1053

1054

1055

1056

1057

1058

1059

1060

1061

1062

1063

1064

1065

1066

1067

1068

1069

1070

1071

1072

1073

1074

1075

1076

1077

1078

1079

1080

1081

1082

1083

1084

1085

1086

1087

1088

1089

1090

1091

1092

1093

1094

1095

1096

1097

1098

1099

1100

1101

1102

1103

1104

1105

1106

1107

1108

1109

1110

1111

1112

1113

1114

1115

1116

1117

1118

1119

1120

1121

1122

1123

1124

1125

1126

1127

1128

1129

1130

1131

1132

1133

1134

1135

1136

1137

1138

1139

1140

1141

1142

1143

1144

1145

1146

1147

1148

1149

1150

1151

1152

1153

1154

1155

1156

1157

1158

1159

1160

1161

1162

1163

1164

1165

1166

1167

1168

1169

1170

1171

1172

1173

1174

1175

1176

1177

1178

1179

1180

1181

1182

1183

1184

1185

1186

1187

1188

1189

1190

1191

1192

1193

1194

1195

1196

1197

1198

1199

1200

1201

1202

1203

1204

1205

1206

1207

1208

1209

1210

1211

1212

1213

1214

1215

1216

1217

1218

1219

1220

1221

1222

1223

1224

1225

1226

1227

1228

1229

1230

1231

1232

1233

1234

1235

1236

1237

1238

1239

1240

1241

1242

1243

1244

1245

1246

1247

1248

1249

1250

1251

1252

1253

1254

1255

1256

1257

1258

1259

1260

1261

1262

1263

1264

1265

1266

1267

1268

1269

1270

1271

1272

1273

1274

1275

1276

1277

1278

1279

1280

1281

1282

1283

1284

1285

1286

1287

1288

1289

1290

1291

1292

1293

1294

1295

1296

1297

1298

1299

1300

1301

1302

1303

1304

1305

1306

1307

1308

1309

1310

1311

1312

1313

1314

1315

1316

1317

1318

1319

1320

1321

1322

1323

1324

1325

1326

1327

1328

1329

1330

1331

1332

1333

1334

1335

1336

1337

1338

1339

1340

1341

1342

1343

1344

1345

1346

1347

1348

1349

1350

1351

1352

1353

1354

1355

1356

1357

1358

1359

1360

1361

1362

1363

1364

1365

1366

1367

1368

1369

1370

1371

1372

1373

1374

1375

1376

1377

1378

1379

1380

1381

1382

1383

1384

1385

1386

1387

1388

1389

1390

1391

1392

1393

1394

1395

1396

1397

1398

1399

1400

1401

1402

1403

1404

1405

1406

1407

1408

1409

1410

1411

1412

1413

1414

1415

1416

1417

1418

1419

1420

1421

1422

1423

1424

1425

1426

1427

1428

1429

1430

1431

1432

1433

1434

1435

1436

1437

1438

1439

1440

1441

1442

1443

1444

1445

1446

1447

1448

1449

1450

1451

1452

1453

1454

1455

1456

1457

1458

1459

1460

1461

1462

1463

1464

1465

1466

1467

1468

1469

1470

1471

1472

1473

1474

1475

1476

1477

1478

1479

1480

1481

1482

1483

1484

1485

1486

1487

1488

1489

1490

1491

1492

1493

1494

1495

1496

1497

1498

1499

1500

1501

1502

1503

1504

1505

1506

1507

1508

1509

1510

1511

1512

1513

1514

1515

1516

1517

1518

1519

1520

1521

1522

1523

1524

1525

1526

1527

1528

1529

1530

1531

1532

1533

1534

1535

1536

1537

1538

1539

1540

1541

1542

1543

1544

1545

1546

1547

1548

1549

1550

1551

1552

1553

1554

1555

1556

1557

1558

1559

1560

1561

1562

1563

1564

1565

1566

1567

1568

1569

1570

1571

1572

1573

1574

1575

1576

1577

1578

1579

1580

1581

1582

1583

1584

1585

1586

1587

1588

1589

1590

1591

1592

1593

1594

1595

1596

1597

1598

1599

1600

1601

1602

1603

1604

1605

1606

1607

1608

1609

1610

1611

1612

1613

1614

1615

1616

1617

1618

1619

1620

1621

1622

1623

1624

1625

1626

1627

1628

1629

1630

1631

1632

1633

1634

1635

1636

1637

1638

1639

1640

1641

1642

1643

1644

1645

1646

1647

1648

1649

1650

1651

1652

1653

1654

1655

1656

1657

1658

1659

1660

1661

1662

1663

1664

1665

1666

1667

1668

1669

1670

1671

1672

1673

1674

1675

1676

1677

1678

1679

1680

1681

1682

1683

1684

1685

1686

1687

1688

1689

1690

1691

1692

1693

1694

1695

1696

1697

1698

1699

1700

1701

1702

1703

1704

1705

1706

1707

1708

1709

1710

1711

1712

1713

1714

1715

1716

1717

1718

1719

1720

1721

1722

1723

1724

1725

1726

1727

1728

1729

1730

1731

1732

1733

1734

1735

1736

1737

1738

1739

1740

1741

1742

1743

1744

1745

1746

1747

1748

1749

1750

1751

1752

1753

1754

1755

1756

1757

1758

1759

1760

1761

1762

1763

1764

1765

1766

1767

1768

1769

1770

1771

1772

1773

1774

1775

1776

1777

1778

1779

1780

1781

1782

1783

1784

1785

1786

1787

1788

1789

1790

1791

1792

1793

1794

1795

1796

1797

1798

1799

1800

1801

1802

1803

1804

1805

1806

1807

1808

1809

1810

1811

1812

1813

1814

1815

1816

1817

1818

1819

1820

1821

1822

1823

1824

1825

1826

1827

1828

1829

1830

1831

1832

1833

1834

1835

1836

1837

1838

1839

1840

1841

1842

1843

1844

1845

1846

1847

1848

1849

1850

1851

1852

1853

1854

1855

1856

1857

1858

1859

1860

1861

1862

1863

1864

1865

1866

1867

1868

1869

1870

1871

1872

1873

1874

1875

1876

1877

1878

1879

1880

1881

1882

1883

1884

1885

1886

1887

1888

1889

1890

1891

1892

1893

1894

1895

1896

1897

1898

1899

1900

1901

1902

1903

1904

1905

1906

1907

1908

1909

1910

1911

1912

1913

1914

1915

1916

1917

1918

1919

1920

1921

1922

1923

1924

1925

1926

1927

1928

1929

1930

1931

1932

1933

1934

1935

1936

1937

1938

1939

1940

1941

1942

1943

1944

1945

1946

1947

1948

1949

1950

1951

1952

1953

1954

1955

1956

1957

1958

1959

1960

1961

1962

1963

1964

1965

1966

1967

1968

1969

1970

1971

1972

1973

1974

1975

1976

1977

1978

1979

1980

1981

1982

1983

1984

1985

1986

1987

1988

1989

1990

1991

1992

1993

1994

1995

1996

1997

1998

1999

2000

2001

2002

2003

2004

2005

2006

2007

2008

2009

2010

2011

2012

2013

2014

2015

2016

2017

2018

2019

2020

2021

2022

2023

2024

2025

2026

2027

2028

2029

2030

2031

2032

2033

2034

2035

2036

2037

2038

2039

2040

2041

2042

2043

2044

2045

2046

2047

2048

2049

2050

2051

2052

2053

2054

2055

2056

2057

2058

2059

2060

2061

2062

2063

2064

2065

2066

2067

2068

2069

2070

2071

2072

2073

2074

2075

2076

2077

2078

2079

2080

2081

2082

2083

2084

2085

2086

2087

2088

2089

2090

2091

2092

2093

2094

2095

2096

2097

2098

2099

2100

2101

2102

2103

2104

2105

2106

2107

2108

2109

2110

2111

2112

2113

2114

2115

2116

2117

2118

2119

2120

2121

2122

2123

2124

2125

2126

2127

2128

2129

2130

2131

2132

2133

2134

2135

2136

2137

2138

2139

2140

2141

2142

2143

2144

2145

2146

2147

2148

2149

2150

2151

2152

2153

2154

2155

2156

2157

2158

2159

2160

2161

2162

2163

2164

2165

2166

2167

2168

2169

2170

2171

2172

2173

2174

2175

2176

2177

2178

2179

2180

2181

2182

2183

2184

2185

2186

2187

2188

2189

2190

2191

2192

2193

2194

2195

2196

2197

2198

2199

2200

2201

2202

2203

2204

2205

2206

2207

2208

2209

2210

2211

2212

2213

2214

2215

2216

2217

2218

2219

2220

2221

2222

2223

2224

2225

2226

2227

2228

2229

2230

2231

2232

2233

2234

2235

2236

2237

2238

2239

2240

2241

2242

2243

2244

2245

2246

2247

2248

2249

2250

2251

2252

2253

2254

2255

2256

2257

2258

2259

2260

2261

2262

2263

2264

2265

2266

2267

2268

2269

2270

2271

2272

2273

2274

2275

2276

2277

2278

2279

2280

2281

2282

2283

2284

2285

2286

2287

2288

2289

2290

2291

2292

2293

2294

2295

2296

2297

2298

2299

2300

2301

2302

2303

2304

2305

2306

2307

2308

2309

2310

2311

2312

2313

2314

2315

2316

2317

2318

2319

2320

2321

2322

2323

2324

2325

2326

2327

2328

2329

2330

2331

2332

2333

2334

2335

2336

2337

2338

2339

2340

2341

2342

2343

2344

2345

2346

2347

2348

2349

2350

2351

2352

2353

2354

2355

2356

2357

2358

2359

2360

2361

2362

2363

2364

2365

2366

2367

2368

2369

2370

2371

2372

2373

2374

2375

2376

2377

2378

2379

2380

2381

2382

2383

2384

2385

2386

2387

2388

2389

2390

2391

2392

2393

2394

2395

2396

2397

2398

2399

2400

2401

2402

2403

2404

2405

2406

2407

2408

2409

2410

2411

2412

2413

2414

2415

2416

2417

2418

2419

2420

2421

2422

2423

2424

2425

2426

2427

2428

2429

2430

2431

2432

2433

2434

2435

2436

2437

2438

2439

2440

2441

2442

2443

2444

2445

2446

2447

2448

2449

2450

2451

2452

2453

2454

2455

2456

2457

2458

2459

2460

2461

2462

2463

2464

2465

2466

2467

2468

2469

2470

2471

2472

2473

2474

2475

2476

2477

2478

2479

2480

2481

2482

2483

2484

2485

2486

2487

2488

2489

2490

2491

2492

2493

2494

2495

2496

2497

2498

2499

2500

2501

2502

2503

2504

2505

2506

2507

2508

2509

2510

2511

2512

2513

2514

2515

2516

2517

2518

2519

2520

2521

2522

2523

2524

2525

2526

2527

2528

2529

2530

2531

2532

2533

2534

2535

2536

2537

2538

2539

2540

2541

2542

2543

2544

2545

2546

2547

2548

2549

2550

2551

2552

2553

2554

2555

2556

2557

2558

2559

2560

2561

2562

2563

2564

2565

2566

2567

2568

2569

2570

2571

2572

2573

2574

2575

2576

2577

2578

2579

2580

2581

2582

2583

2584

2585

2586

2587

2588

2589

2590

2591

2592

2593

2594

2595

2596

2597

2598

2599

2600

2601

2602

2603

2604

2605

2606

2607

2608

2609

2610

2611

2612

2613

2614

2615

2616

2617

2618

2619

2620

2621

2622

2623

2624

2625

2626

2627

2628

2629

2630

2631

2632

2633

2634

2635

2636

2637

2638

2639

2640

2641

2642

2643

2644

2645

2646

2647

2648

2649

2650

2651

2652

2653

2654

2655

2656

2657

2658

2659

2660

2661

2662

2663

2664

2665

2666

2667

2668

2669

2670

2671

2672

2673

2674

2675

2676

2677

2678

2679

2680

2681

2682

2683

2684

2685

2686

2687

2688

2689

2690

2691

2692

2693

2694

2695

2696

2697

2698

2699

2700

2701

2702

2703

2704

2705

2706

2707

2708

2709

2710

2711

2712

2713

2714

2715

2716

2717

2718

2719

2720

2721

2722

2723

2724

2725

2726

2727

2728

2729

2730

2731

2732

2733

2734

2735

2736

2737

2738

2739

2740

2741

2742

2743

2744

2745

2746

2747

2748

2749

2750

2751

2752

2753

2754

2755

2756

2757

2758

2759

2760

2761

2762

2763

2764

2765

2766

2767

2768

2769

2770

2771

2772

2773

2774

2775

2776

2777

2778

2779

2780

2781

2782

2783

2784

2785

2786

2787

2788

2789

2790

2791

2792

2793

2794

2795

2796

2797

2798

2799

2800

2801

2802

2803

2804

2805

2806

2807

2808

2809

2810

2811

2812

2813

2814

2815

2816

2817

2818

2819

2820

2821

2822

2823

2824

2825

2826

2827

2828

2829

2830

2831

2832

2833

2834

2835

2836

2837

2838

2839

2840

2841

2842

2843

2844

2845

2846

2847

2848

2849

2850

2851

2852

2853

2854

2855

2856

2857

2858

2859

2860

2861

2862

2863

2864

2865

2866

2867

2868

2869

2870

2871

2872

2873

2874

2875

2876

2877

2878

2879

2880

2881

2882

2883

2884

2885

2886

2887

2888

2889

2890

2891

2892

2893

2894

2895

2896

2897

2898

2899

2900

2901

2902

2903

2904

2905

2906

2907

2908

2909

2910

2911

2912

2913

2914

2915

2916

2917

2918

2919

2920

2921

2922

2923

2924

2925

2926

2927

2928

2929

2930

2931

2932

2933

2934

2935

2936

2937

2938

2939

2940

2941

2942

2943

2944

2945

2946

2947

2948

2949

2950

2951

2952

2953

2954

2955

2956

2957

2958

2959

2960

2961

2962

2963

2964

2965

2966

2967

2968

2969

2970

2971

2972

2973

2974

2975

2976

2977

2978

2979

2980

2981

2982

2983

2984

2985

2986

2987

2988

2989

2990

2991

2992

2993

2994

2995

2996

2997

2998

2999

3000

3001

3002

3003

3004

3005

3006

3007

3008

3009

3010

3011

3012

3013

3014

3015

3016

3017

3018

3019

3020

3021

3022

3023

3024

3025

3026

3027

3028

3029

3030

3031

3032

3033

3034

3035

3036

3037

3038

3039

3040

3041

3042

3043

3044

3045

3046

3047

3048

3049

3050

3051

3052

3053

3054

3055

3056

3057

3058

3059

3060

3061

3062

3063

3064

""" 

GLPK Backend 

  

AUTHORS: 

  

- Nathann Cohen (2010-10): initial implementation 

- John Perry (2012-01): glp_simplex preprocessing 

- John Perry and Raniere Gaia Silva (2012-03): solver parameters 

- Christian Kuper (2012-10): Additions for sensitivity analysis 

""" 

  

#***************************************************************************** 

# Copyright (C) 2010 Nathann Cohen <nathann.cohen@gmail.com> 

# 

# This program is free software: you can redistribute it and/or modify 

# it under the terms of the GNU General Public License as published by 

# the Free Software Foundation, either version 2 of the License, or 

# (at your option) any later version. 

# http://www.gnu.org/licenses/ 

#***************************************************************************** 

from __future__ import print_function 

  

from cysignals.memory cimport sig_malloc, sig_free 

from cysignals.signals cimport sig_on, sig_off 

  

from sage.ext.memory_allocator cimport MemoryAllocator 

from sage.numerical.mip import MIPSolverException 

from sage.libs.glpk.error import GLPKError 

from sage.libs.glpk.constants cimport * 

from sage.libs.glpk.lp cimport * 

from libc.float cimport DBL_MAX 

from libc.limits cimport INT_MAX 

  

  

cdef class GLPKBackend(GenericBackend): 

  

""" 

MIP Backend that uses the GLPK solver. 

  

TESTS: 

  

General backend testsuite:: 

  

sage: p = MixedIntegerLinearProgram(solver="GLPK") 

sage: TestSuite(p.get_backend()).run(skip="_test_pickling") 

""" 

  

def __cinit__(self, maximization = True): 

""" 

Constructor 

  

EXAMPLES:: 

  

sage: p = MixedIntegerLinearProgram(solver="GLPK") 

""" 

self.lp = glp_create_prob() 

self.simplex_or_intopt = glp_intopt_only 

self.smcp = <glp_smcp* > sig_malloc(sizeof(glp_smcp)) 

glp_init_smcp(self.smcp) 

self.iocp = <glp_iocp* > sig_malloc(sizeof(glp_iocp)) 

glp_init_iocp(self.iocp) 

  

self.iocp.cb_func = glp_callback # callback function 

self.iocp.cb_info = <void *> &(self.search_tree_data) # callback data 

self.iocp.presolve = GLP_ON 

self.set_verbosity(0) 

self.obj_constant_term = 0.0 

  

if maximization: 

self.set_sense(+1) 

else: 

self.set_sense(-1) 

  

cpdef int add_variable(self, lower_bound=0.0, upper_bound=None, binary=False, continuous=False, integer=False, obj=0.0, name=None) except -1: 

""" 

Add a variable. 

  

This amounts to adding a new column to the matrix. By default, 

the variable is both positive, real and the coefficient in the 

objective function is 0.0. 

  

INPUT: 

  

- ``lower_bound`` - the lower bound of the variable (default: 0) 

  

- ``upper_bound`` - the upper bound of the variable (default: ``None``) 

  

- ``binary`` - ``True`` if the variable is binary (default: ``False``). 

  

- ``continuous`` - ``True`` if the variable is binary (default: ``True``). 

  

- ``integer`` - ``True`` if the variable is binary (default: ``False``). 

  

- ``obj`` - (optional) coefficient of this variable in the objective function (default: 0.0) 

  

- ``name`` - an optional name for the newly added variable (default: ``None``). 

  

OUTPUT: The index of the newly created variable 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.ncols() 

0 

sage: p.add_variable() 

0 

sage: p.ncols() 

1 

sage: p.add_variable(binary=True) 

1 

sage: p.add_variable(lower_bound=-2.0, integer=True) 

2 

sage: p.add_variable(continuous=True, integer=True) 

Traceback (most recent call last): 

... 

ValueError: ... 

sage: p.add_variable(name='x', obj=1.0) 

3 

sage: p.col_name(3) 

'x' 

sage: p.objective_coefficient(3) 

1.0 

""" 

cdef int vtype = int(bool(binary)) + int(bool(continuous)) + int(bool(integer)) 

if vtype == 0: 

continuous = True 

elif vtype != 1: 

raise ValueError("Exactly one parameter of 'binary', 'integer' and 'continuous' must be 'True'.") 

  

glp_add_cols(self.lp, 1) 

cdef int n_var = glp_get_num_cols(self.lp) 

  

self.variable_lower_bound(n_var - 1, lower_bound) 

self.variable_upper_bound(n_var - 1, upper_bound) 

  

if continuous: 

glp_set_col_kind(self.lp, n_var, GLP_CV) 

elif binary: 

glp_set_col_kind(self.lp, n_var, GLP_BV) 

elif integer: 

glp_set_col_kind(self.lp, n_var, GLP_IV) 

  

if name is not None: 

glp_set_col_name(self.lp, n_var, name) 

  

if obj: 

self.objective_coefficient(n_var - 1, obj) 

  

return n_var - 1 

  

cpdef int add_variables(self, int number, lower_bound=0.0, upper_bound=None, binary=False, continuous=False, integer=False, obj=0.0, names=None) except -1: 

""" 

Add ``number`` new variables. 

  

This amounts to adding new columns to the matrix. By default, 

the variables are both positive, real and their coefficient in 

the objective function is 0.0. 

  

INPUT: 

  

- ``n`` - the number of new variables (must be > 0) 

  

- ``lower_bound`` - the lower bound of the variable (default: 0) 

  

- ``upper_bound`` - the upper bound of the variable (default: ``None``) 

  

- ``binary`` - ``True`` if the variable is binary (default: ``False``). 

  

- ``continuous`` - ``True`` if the variable is binary (default: ``True``). 

  

- ``integer`` - ``True`` if the variable is binary (default: ``False``). 

  

- ``obj`` - (optional) coefficient of all variables in the objective function (default: 0.0) 

  

- ``names`` - optional list of names (default: ``None``) 

  

OUTPUT: The index of the variable created last. 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.ncols() 

0 

sage: p.add_variables(5) 

4 

sage: p.ncols() 

5 

sage: p.add_variables(2, lower_bound=-2.0, integer=True, obj=42.0, names=['a','b']) 

6 

  

TESTS: 

  

Check that arguments are used:: 

  

sage: p.col_bounds(5) # tol 1e-8 

(-2.0, None) 

sage: p.is_variable_integer(5) 

True 

sage: p.col_name(5) 

'a' 

sage: p.objective_coefficient(5) 

42.0 

""" 

cdef int vtype = int(bool(binary)) + int(bool(continuous)) + int(bool(integer)) 

if vtype == 0: 

continuous = True 

elif vtype != 1: 

raise ValueError("Exactly one parameter of 'binary', 'integer' and 'continuous' must be 'True'.") 

  

glp_add_cols(self.lp, number) 

  

cdef int n_var 

n_var = glp_get_num_cols(self.lp) 

  

cdef int i 

  

for 0<= i < number: 

self.variable_lower_bound(n_var - i - 1, lower_bound) 

self.variable_upper_bound(n_var - i - 1, upper_bound) 

if continuous: 

glp_set_col_kind(self.lp, n_var - i, GLP_CV) 

elif binary: 

glp_set_col_kind(self.lp, n_var - i, GLP_BV) 

elif integer: 

glp_set_col_kind(self.lp, n_var - i, GLP_IV) 

  

if obj: 

self.objective_coefficient(n_var - i - 1, obj) 

  

if names is not None: 

glp_set_col_name(self.lp, n_var - i, names[number - i - 1]) 

  

return n_var - 1 

  

cpdef set_variable_type(self, int variable, int vtype): 

""" 

Set the type of a variable 

  

INPUT: 

  

- ``variable`` (integer) -- the variable's id 

  

- ``vtype`` (integer) : 

  

* 1 Integer 

* 0 Binary 

* -1 Real 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.ncols() 

0 

sage: p.add_variable() 

0 

sage: p.set_variable_type(0,1) 

sage: p.is_variable_integer(0) 

True 

""" 

  

if vtype==1: 

glp_set_col_kind(self.lp, variable+1, GLP_IV) 

  

elif vtype==0: 

glp_set_col_kind(self.lp, variable+1, GLP_BV) 

  

else: 

glp_set_col_kind(self.lp, variable+1, GLP_CV) 

  

cpdef set_sense(self, int sense): 

""" 

Set the direction (maximization/minimization). 

  

INPUT: 

  

- ``sense`` (integer) : 

  

* +1 => Maximization 

* -1 => Minimization 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.is_maximization() 

True 

sage: p.set_sense(-1) 

sage: p.is_maximization() 

False 

""" 

if sense == 1: 

glp_set_obj_dir(self.lp, GLP_MAX) 

else: 

glp_set_obj_dir(self.lp, GLP_MIN) 

  

cpdef objective_coefficient(self, int variable, coeff=None): 

""" 

Set or get the coefficient of a variable in the objective function 

  

INPUT: 

  

- ``variable`` (integer) -- the variable's id 

  

- ``coeff`` (double) -- its coefficient or ``None`` for 

reading (default: ``None``) 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.add_variable() 

0 

sage: p.objective_coefficient(0) 

0.0 

sage: p.objective_coefficient(0,2) 

sage: p.objective_coefficient(0) 

2.0 

""" 

if coeff is None: 

return glp_get_obj_coef(self.lp, variable + 1) 

else: 

glp_set_obj_coef(self.lp, variable + 1, coeff) 

  

cpdef problem_name(self, char * name = NULL): 

""" 

Return or define the problem's name 

  

INPUT: 

  

- ``name`` (``char *``) -- the problem's name. When set to 

``NULL`` (default), the method returns the problem's name. 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.problem_name("There once was a french fry") 

sage: print(p.problem_name()) 

There once was a french fry 

""" 

cdef char * n 

  

if name == NULL: 

n = <char *> glp_get_prob_name(self.lp) 

if n == NULL: 

return "" 

else: 

return n 

  

else: 

if len(name) > 255: 

raise ValueError("Problem name for GLPK must not be longer than 255 characters.") 

glp_set_prob_name(self.lp, name) 

  

cpdef set_objective(self, list coeff, d = 0.0): 

""" 

Set the objective function. 

  

INPUT: 

  

- ``coeff`` - a list of real values, whose ith element is the 

coefficient of the ith variable in the objective function. 

  

- ``d`` (double) -- the constant term in the linear function (set to `0` by default) 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.add_variables(5) 

4 

sage: p.set_objective([1, 1, 2, 1, 3]) 

sage: [p.objective_coefficient(x) for x in range(5)] 

[1.0, 1.0, 2.0, 1.0, 3.0] 

""" 

cdef int i 

  

for i,v in enumerate(coeff): 

glp_set_obj_coef(self.lp, i+1, v) 

  

glp_set_obj_coef(self.lp, 0, d) 

  

self.obj_constant_term = d 

  

cpdef set_verbosity(self, int level): 

""" 

Set the verbosity level 

  

INPUT: 

  

- ``level`` (integer) -- From 0 (no verbosity) to 3. 

  

EXAMPLES:: 

  

sage: p.<x> = MixedIntegerLinearProgram(solver="GLPK") 

sage: p.add_constraint(10 * x[0] <= 1) 

sage: p.add_constraint(5 * x[1] <= 1) 

sage: p.set_objective(x[0] + x[1]) 

sage: p.solve() 

0.30000000000000004 

sage: p.get_backend().set_verbosity(3) 

sage: p.solve() 

GLPK Integer Optimizer, v4.63 

2 rows, 2 columns, 2 non-zeros 

0 integer variables, none of which are binary 

Preprocessing... 

Objective value = 3.000000000e-01 

INTEGER OPTIMAL SOLUTION FOUND BY MIP PREPROCESSOR 

0.30000000000000004 

  

:: 

  

sage: p.<x> = MixedIntegerLinearProgram(solver="GLPK/exact") 

sage: p.add_constraint(10 * x[0] <= 1) 

sage: p.add_constraint(5 * x[1] <= 1) 

sage: p.set_objective(x[0] + x[1]) 

sage: p.solve() 

0.3 

sage: p.get_backend().set_verbosity(2) 

sage: p.solve() 

* 2: objval = 0.3 (0) 

* 2: objval = 0.3 (0) 

0.3 

sage: p.get_backend().set_verbosity(3) 

sage: p.solve() 

glp_exact: 2 rows, 2 columns, 2 non-zeros 

GNU MP bignum library is being used 

* 2: objval = 0.3 (0) 

* 2: objval = 0.3 (0) 

OPTIMAL SOLUTION FOUND 

0.3 

""" 

if level == 0: 

self.iocp.msg_lev = GLP_MSG_OFF 

self.smcp.msg_lev = GLP_MSG_OFF 

elif level == 1: 

self.iocp.msg_lev = GLP_MSG_ERR 

self.smcp.msg_lev = GLP_MSG_ERR 

elif level == 2: 

self.iocp.msg_lev = GLP_MSG_ON 

self.smcp.msg_lev = GLP_MSG_ON 

else: 

self.iocp.msg_lev = GLP_MSG_ALL 

self.smcp.msg_lev = GLP_MSG_ALL 

  

cpdef remove_constraint(self, int i): 

r""" 

Remove a constraint from self. 

  

INPUT: 

  

- ``i`` -- index of the constraint to remove 

  

EXAMPLES:: 

  

sage: p = MixedIntegerLinearProgram(solver='GLPK') 

sage: x, y = p['x'], p['y'] 

sage: p.add_constraint(2*x + 3*y <= 6) 

sage: p.add_constraint(3*x + 2*y <= 6) 

sage: p.add_constraint(x >= 0) 

sage: p.set_objective(x + y + 7) 

sage: p.set_integer(x); p.set_integer(y) 

sage: p.solve() 

9.0 

sage: p.remove_constraint(0) 

sage: p.solve() 

10.0 

  

Removing fancy constraints does not make Sage crash:: 

  

sage: MixedIntegerLinearProgram(solver = "GLPK").remove_constraint(-2) 

Traceback (most recent call last): 

... 

ValueError: The constraint's index i must satisfy 0 <= i < number_of_constraints 

""" 

cdef int rows[2] 

  

if i < 0 or i >= glp_get_num_rows(self.lp): 

raise ValueError("The constraint's index i must satisfy 0 <= i < number_of_constraints") 

  

rows[1] = i + 1 

glp_del_rows(self.lp, 1, rows) 

glp_std_basis(self.lp) 

  

cpdef remove_constraints(self, constraints): 

r""" 

Remove several constraints. 

  

INPUT: 

  

- ``constraints`` -- an iterable containing the indices of the rows to remove. 

  

EXAMPLES:: 

  

sage: p = MixedIntegerLinearProgram(solver='GLPK') 

sage: x, y = p['x'], p['y'] 

sage: p.add_constraint(2*x + 3*y <= 6) 

sage: p.add_constraint(3*x + 2*y <= 6) 

sage: p.add_constraint(x >= 0) 

sage: p.set_objective(x + y + 7) 

sage: p.set_integer(x); p.set_integer(y) 

sage: p.solve() 

9.0 

sage: p.remove_constraints([0]) 

sage: p.solve() 

10.0 

sage: p.get_values([x,y]) 

[0.0, 3.0] 

  

TESTS: 

  

Removing fancy constraints does not make Sage crash:: 

  

sage: MixedIntegerLinearProgram(solver= "GLPK").remove_constraints([0, -2]) 

Traceback (most recent call last): 

... 

ValueError: The constraint's index i must satisfy 0 <= i < number_of_constraints 

""" 

cdef int i, c 

cdef int m = len(constraints) 

cdef int * rows = <int *>sig_malloc((m + 1) * sizeof(int *)) 

cdef int nrows = glp_get_num_rows(self.lp) 

  

for i in xrange(m): 

  

c = constraints[i] 

if c < 0 or c >= nrows: 

sig_free(rows) 

raise ValueError("The constraint's index i must satisfy 0 <= i < number_of_constraints") 

  

rows[i+1] = c + 1 

  

glp_del_rows(self.lp, m, rows) 

sig_free(rows) 

glp_std_basis(self.lp) 

  

cpdef add_linear_constraint(self, coefficients, lower_bound, upper_bound, name=None): 

""" 

Add a linear constraint. 

  

INPUT: 

  

- ``coefficients`` an iterable with ``(c,v)`` pairs where ``c`` 

is a variable index (integer) and ``v`` is a value (real 

value). 

  

- ``lower_bound`` - a lower bound, either a real value or ``None`` 

  

- ``upper_bound`` - an upper bound, either a real value or ``None`` 

  

- ``name`` - an optional name for this row (default: ``None``) 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.add_variables(5) 

4 

sage: p.add_linear_constraint(list(zip(range(5), range(5))), 2.0, 2.0) 

sage: p.row(0) 

([4, 3, 2, 1], [4.0, 3.0, 2.0, 1.0]) 

sage: p.row_bounds(0) 

(2.0, 2.0) 

sage: p.add_linear_constraint(list(zip(range(5), range(5))), 1.0, 1.0, name='foo') 

sage: p.row_name(1) 

'foo' 

  

TESTS: 

  

This used to crash Sage, but was fixed in :trac:`19525`:: 

  

sage: p = MixedIntegerLinearProgram(solver="glpk") 

sage: q = MixedIntegerLinearProgram(solver="glpk") 

sage: q.add_constraint(p.new_variable()[0] <= 1) 

Traceback (most recent call last): 

... 

GLPKError: glp_set_mat_row: i = 1; len = 1; invalid row length 

Error detected in file api/prob1.c at line ... 

""" 

if lower_bound is None and upper_bound is None: 

raise ValueError("At least one of 'upper_bound' or 'lower_bound' must be set.") 

  

glp_add_rows(self.lp, 1) 

cdef int n = glp_get_num_rows(self.lp) 

  

cdef MemoryAllocator mem = MemoryAllocator() 

cdef int * row_i 

cdef double * row_values 

  

row_i = <int*>mem.allocarray(len(coefficients)+1, sizeof(int)) 

row_values = <double*>mem.allocarray(len(coefficients)+1, sizeof(double)) 

  

cdef Py_ssize_t i = 1 

for c,v in coefficients: 

row_i[i] = c+1 

row_values[i] = v 

i += 1 

  

sig_on() 

glp_set_mat_row(self.lp, n, len(coefficients), row_i, row_values) 

sig_off() 

  

if upper_bound is not None and lower_bound is None: 

glp_set_row_bnds(self.lp, n, GLP_UP, upper_bound, upper_bound) 

elif lower_bound is not None and upper_bound is None: 

glp_set_row_bnds(self.lp, n, GLP_LO, lower_bound, lower_bound) 

elif upper_bound is not None and lower_bound is not None: 

if lower_bound == upper_bound: 

glp_set_row_bnds(self.lp, n, GLP_FX, lower_bound, upper_bound) 

else: 

glp_set_row_bnds(self.lp, n, GLP_DB, lower_bound, upper_bound) 

  

if name is not None: 

glp_set_row_name(self.lp, n, name) 

  

cpdef add_linear_constraints(self, int number, lower_bound, upper_bound, names=None): 

""" 

Add ``'number`` linear constraints. 

  

INPUT: 

  

- ``number`` (integer) -- the number of constraints to add. 

  

- ``lower_bound`` - a lower bound, either a real value or ``None`` 

  

- ``upper_bound`` - an upper bound, either a real value or ``None`` 

  

- ``names`` - an optional list of names (default: ``None``) 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.add_variables(5) 

4 

sage: p.add_linear_constraints(5, None, 2) 

sage: p.row(4) 

([], []) 

sage: p.row_bounds(4) 

(None, 2.0) 

sage: p.add_linear_constraints(2, None, 2, names=['foo','bar']) 

""" 

if lower_bound is None and upper_bound is None: 

raise ValueError("At least one of 'upper_bound' or 'lower_bound' must be set.") 

  

glp_add_rows(self.lp, number) 

cdef int n = glp_get_num_rows(self.lp) 

  

cdef int i 

for 0<= i < number: 

if upper_bound is not None and lower_bound is None: 

glp_set_row_bnds(self.lp, n-i, GLP_UP, upper_bound, upper_bound) 

elif lower_bound is not None and upper_bound is None: 

glp_set_row_bnds(self.lp, n-i, GLP_LO, lower_bound, lower_bound) 

elif upper_bound is not None and lower_bound is not None: 

if lower_bound == upper_bound: 

glp_set_row_bnds(self.lp, n-i, GLP_FX, lower_bound, upper_bound) 

else: 

glp_set_row_bnds(self.lp, n-i, GLP_DB, lower_bound, upper_bound) 

if names is not None: 

glp_set_row_name(self.lp, n-i, names[number-i-1]) 

  

cpdef row(self, int index): 

r""" 

Return a row 

  

INPUT: 

  

- ``index`` (integer) -- the constraint's id. 

  

OUTPUT: 

  

A pair ``(indices, coeffs)`` where ``indices`` lists the 

entries whose coefficient is nonzero, and to which ``coeffs`` 

associates their coefficient on the model of the 

``add_linear_constraint`` method. 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.add_variables(5) 

4 

sage: p.add_linear_constraint(list(zip(range(5), range(5))), 2, 2) 

sage: p.row(0) 

([4, 3, 2, 1], [4.0, 3.0, 2.0, 1.0]) 

sage: p.row_bounds(0) 

(2.0, 2.0) 

""" 

cdef int n = glp_get_num_cols(self.lp) 

cdef MemoryAllocator mem = MemoryAllocator() 

cdef int * c_indices = <int*>mem.allocarray(n+1, sizeof(int)) 

cdef double * c_values = <double*>mem.allocarray(n+1, sizeof(double)) 

cdef list indices = [] 

cdef list values = [] 

cdef int i,j 

  

i = glp_get_mat_row(self.lp, index + 1, c_indices, c_values) 

for 0 < j <= i: 

indices.append(c_indices[j]-1) 

values.append(c_values[j]) 

  

return (indices, values) 

  

cpdef row_bounds(self, int index): 

""" 

Return the bounds of a specific constraint. 

  

INPUT: 

  

- ``index`` (integer) -- the constraint's id. 

  

OUTPUT: 

  

A pair ``(lower_bound, upper_bound)``. Each of them can be set 

to ``None`` if the constraint is not bounded in the 

corresponding direction, and is a real value otherwise. 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.add_variables(5) 

4 

sage: p.add_linear_constraint(list(zip(range(5), range(5))), 2, 2) 

sage: p.row(0) 

([4, 3, 2, 1], [4.0, 3.0, 2.0, 1.0]) 

sage: p.row_bounds(0) 

(2.0, 2.0) 

""" 

cdef double ub 

cdef double lb 

  

ub = glp_get_row_ub(self.lp, index + 1) 

lb = glp_get_row_lb(self.lp, index +1) 

  

return ( 

(lb if lb != -DBL_MAX else None), 

(ub if ub != +DBL_MAX else None) 

) 

  

cpdef col_bounds(self, int index): 

""" 

Return the bounds of a specific variable. 

  

INPUT: 

  

- ``index`` (integer) -- the variable's id. 

  

OUTPUT: 

  

A pair ``(lower_bound, upper_bound)``. Each of them can be set 

to ``None`` if the variable is not bounded in the 

corresponding direction, and is a real value otherwise. 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.add_variable() 

0 

sage: p.col_bounds(0) 

(0.0, None) 

sage: p.variable_upper_bound(0, 5) 

sage: p.col_bounds(0) 

(0.0, 5.0) 

""" 

  

cdef double ub 

cdef double lb 

  

ub = glp_get_col_ub(self.lp, index +1) 

lb = glp_get_col_lb(self.lp, index +1) 

  

return ( 

(lb if lb != -DBL_MAX else None), 

(ub if ub != +DBL_MAX else None) 

) 

  

cpdef add_col(self, list indices, list coeffs): 

""" 

Add a column. 

  

INPUT: 

  

- ``indices`` (list of integers) -- this list contains the 

indices of the constraints in which the variable's 

coefficient is nonzero 

  

- ``coeffs`` (list of real values) -- associates a coefficient 

to the variable in each of the constraints in which it 

appears. Namely, the ith entry of ``coeffs`` corresponds to 

the coefficient of the variable in the constraint 

represented by the ith entry in ``indices``. 

  

.. NOTE:: 

  

``indices`` and ``coeffs`` are expected to be of the same 

length. 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.ncols() 

0 

sage: p.nrows() 

0 

sage: p.add_linear_constraints(5, 0, None) 

sage: p.add_col(range(5), range(5)) 

sage: p.nrows() 

5 

""" 

  

glp_add_cols(self.lp, 1) 

cdef int n = glp_get_num_cols(self.lp) 

  

cdef int * col_i 

cdef double * col_values 

  

col_i = <int *> sig_malloc((len(indices)+1) * sizeof(int)) 

col_values = <double *> sig_malloc((len(indices)+1) * sizeof(double)) 

  

for i,v in enumerate(indices): 

col_i[i+1] = v+1 

for i,v in enumerate(coeffs): 

col_values[i+1] = v 

  

glp_set_mat_col(self.lp, n, len(indices), col_i, col_values) 

glp_set_col_bnds(self.lp, n, GLP_LO, 0,0) 

  

  

cpdef int solve(self) except -1: 

""" 

Solve the problem. 

  

Sage uses GLPK's implementation of the branch-and-cut algorithm 

(``glp_intopt``) to solve the mixed-integer linear program. 

This algorithm can be requested explicitly by setting the solver 

parameter "simplex_or_intopt" to "intopt_only". 

(If all variables are continuous, the algorithm reduces to solving the 

linear program by the simplex method.) 

  

EXAMPLES:: 

  

sage: lp = MixedIntegerLinearProgram(solver = 'GLPK', maximization = False) 

sage: x, y = lp[0], lp[1] 

sage: lp.add_constraint(-2*x + y <= 1) 

sage: lp.add_constraint(x - y <= 1) 

sage: lp.add_constraint(x + y >= 2) 

sage: lp.set_objective(x + y) 

sage: lp.set_integer(x) 

sage: lp.set_integer(y) 

sage: lp.solve() 

2.0 

sage: lp.get_values([x, y]) 

[1.0, 1.0] 

  

.. NOTE:: 

  

This method raises ``MIPSolverException`` exceptions when 

the solution can not be computed for any reason (none 

exists, or the LP solver was not able to find it, etc...) 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.add_linear_constraints(5, 0, None) 

sage: p.add_col(range(5), range(5)) 

sage: p.solve() 

0 

sage: p.objective_coefficient(0,1) 

sage: p.solve() 

Traceback (most recent call last): 

... 

MIPSolverException: ... 

  

.. WARNING:: 

  

Sage uses GLPK's ``glp_intopt`` to find solutions. 

This routine sometimes FAILS CATASTROPHICALLY 

when given a system it cannot solve. (:trac:`12309`.) 

Here, "catastrophic" can mean either "infinite loop" or 

segmentation fault. Upstream considers this behavior 

"essentially innate" to their design, and suggests 

preprocessing it with ``glp_simplex`` first. 

Thus, if you suspect that your system is infeasible, 

set the ``preprocessing`` option first. 

  

EXAMPLES:: 

  

sage: lp = MixedIntegerLinearProgram(solver = "GLPK") 

sage: v = lp.new_variable(nonnegative=True) 

sage: lp.add_constraint(v[1] +v[2] -2.0 *v[3], max=-1.0) 

sage: lp.add_constraint(v[0] -4.0/3 *v[1] +1.0/3 *v[2], max=-1.0/3) 

sage: lp.add_constraint(v[0] +0.5 *v[1] -0.5 *v[2] +0.25 *v[3], max=-0.25) 

sage: lp.solve() 

0.0 

sage: lp.add_constraint(v[0] +4.0 *v[1] -v[2] +v[3], max=-1.0) 

sage: lp.solve() 

Traceback (most recent call last): 

... 

GLPKError: Assertion failed: ... 

sage: lp.solver_parameter("simplex_or_intopt", "simplex_then_intopt") 

sage: lp.solve() 

Traceback (most recent call last): 

... 

MIPSolverException: GLPK: Problem has no feasible solution 

  

If we switch to "simplex_only", the integrality constraints are ignored, 

and we get an optimal solution to the continuous relaxation. 

  

EXAMPLES:: 

  

sage: lp = MixedIntegerLinearProgram(solver = 'GLPK', maximization = False) 

sage: x, y = lp[0], lp[1] 

sage: lp.add_constraint(-2*x + y <= 1) 

sage: lp.add_constraint(x - y <= 1) 

sage: lp.add_constraint(x + y >= 2) 

sage: lp.set_objective(x + y) 

sage: lp.set_integer(x) 

sage: lp.set_integer(y) 

sage: lp.solver_parameter("simplex_or_intopt", "simplex_only") # use simplex only 

sage: lp.solve() 

2.0 

sage: lp.get_values([x, y]) 

[1.5, 0.5] 

  

If one solves a linear program and wishes to access dual information 

(`get_col_dual` etc.) or tableau data (`get_row_stat` etc.), 

one needs to switch to "simplex_only" before solving. 

  

GLPK also has an exact rational simplex solver. The only 

access to data is via double-precision floats, however. It 

reconstructs rationals from doubles and also provides results 

as doubles. 

  

EXAMPLES:: 

  

sage: lp.solver_parameter("simplex_or_intopt", "exact_simplex_only") # use exact simplex only 

sage: lp.solve() 

2.0 

sage: lp.get_values([x, y]) 

[1.5, 0.5] 

  

If you need the rational solution, you need to retrieve the 

basis information via ``get_col_stat`` and ``get_row_stat`` 

and calculate the corresponding basic solution. Below we only 

test that the basis information is indeed available. 

Calculating the corresponding basic solution is left as an 

exercise. 

  

EXAMPLES:: 

  

sage: lp.get_backend().get_row_stat(0) 

1 

sage: lp.get_backend().get_col_stat(0) 

1 

  

Below we test that integers that can be exactly represented by 

IEEE 754 double-precision floating point numbers survive the 

rational reconstruction done by ``glp_exact`` and the subsequent 

conversion to double-precision floating point numbers. 

  

EXAMPLES:: 

  

sage: lp = MixedIntegerLinearProgram(solver = 'GLPK', maximization = True) 

sage: test = 2^53 - 43 

sage: lp.solver_parameter("simplex_or_intopt", "exact_simplex_only") # use exact simplex only 

sage: x = lp[0] 

sage: lp.add_constraint(x <= test) 

sage: lp.set_objective(x) 

sage: lp.solve() == test # yes, we want an exact comparison here 

True 

sage: lp.get_values(x) == test # yes, we want an exact comparison here 

True 

  

Below we test that GLPK backend can detect unboundedness in 

"simplex_only" mode (:trac:`18838`). 

  

EXAMPLES:: 

  

sage: lp = MixedIntegerLinearProgram(maximization=True, solver = "GLPK") 

sage: lp.set_objective(lp[0]) 

sage: lp.solver_parameter("simplex_or_intopt", "simplex_only") 

sage: lp.solve() 

Traceback (most recent call last): 

... 

MIPSolverException: GLPK: Problem has unbounded solution 

sage: lp.set_objective(lp[1]) 

sage: lp.solver_parameter("primal_v_dual", "GLP_DUAL") 

sage: lp.solve() 

Traceback (most recent call last): 

... 

MIPSolverException: GLPK: Problem has unbounded solution 

sage: lp.solver_parameter("simplex_or_intopt", "simplex_then_intopt") 

sage: lp.solve() 

Traceback (most recent call last): 

... 

MIPSolverException: GLPK: The LP (relaxation) problem has no dual feasible solution 

sage: lp.solver_parameter("simplex_or_intopt", "intopt_only") 

sage: lp.solve() 

Traceback (most recent call last): 

... 

MIPSolverException: GLPK: The LP (relaxation) problem has no dual feasible solution 

sage: lp.set_max(lp[1],5) 

sage: lp.solve() 

5.0 

  

Solving a LP within the acceptable gap. No exception is raised, even if 

the result is not optimal. To do this, we try to compute the maximum 

number of disjoint balls (of diameter 1) in a hypercube:: 

  

sage: g = graphs.CubeGraph(9) 

sage: p = MixedIntegerLinearProgram(solver="GLPK") 

sage: p.solver_parameter("mip_gap_tolerance",100) 

sage: b = p.new_variable(binary=True) 

sage: p.set_objective(p.sum(b[v] for v in g)) 

sage: for v in g: 

....: p.add_constraint(b[v]+p.sum(b[u] for u in g.neighbors(v)) <= 1) 

sage: p.add_constraint(b[v] == 1) # Force an easy non-0 solution 

sage: p.solve() # rel tol 100 

1 

  

Same, now with a time limit:: 

  

sage: p.solver_parameter("mip_gap_tolerance",1) 

sage: p.solver_parameter("timelimit",3.0) 

sage: p.solve() # rel tol 100 

1 

""" 

  

cdef int solve_status 

cdef int solution_status 

global solve_status_msg 

global solution_status_msg 

  

if (self.simplex_or_intopt == glp_simplex_only 

or self.simplex_or_intopt == glp_simplex_then_intopt 

or self.simplex_or_intopt == glp_exact_simplex_only): 

if self.simplex_or_intopt == glp_exact_simplex_only: 

solve_status = glp_exact(self.lp, self.smcp) 

else: 

solve_status = glp_simplex(self.lp, self.smcp) 

solution_status = glp_get_status(self.lp) 

  

if ((self.simplex_or_intopt == glp_intopt_only) 

or (self.simplex_or_intopt == glp_simplex_then_intopt) and (solution_status != GLP_UNDEF) and (solution_status != GLP_NOFEAS)): 

sig_on() 

solve_status = glp_intopt(self.lp, self.iocp) 

solution_status = glp_mip_status(self.lp) 

sig_off() 

  

if solution_status == GLP_OPT: 

pass 

elif (solution_status == GLP_FEAS) and (solve_status == GLP_ETMLIM or solve_status == GLP_EITLIM \ 

or solve_status == GLP_EMIPGAP or solve_status == GLP_EOBJLL or solve_status == GLP_EOBJUL): 

# no exception when time limit or iteration limit or mip gap tolerances or objective limits reached. 

pass 

elif solution_status == GLP_UNDEF: 

raise MIPSolverException("GLPK: "+solve_status_msg.get(solve_status, "unknown error during call to GLPK : "+str(solve_status))) 

else: 

raise MIPSolverException("GLPK: "+solution_status_msg.get(solution_status, "unknown error during call to GLPK : "+str(solution_status))) 

return 0 

  

cpdef get_objective_value(self): 

""" 

Returns the value of the objective function. 

  

.. NOTE:: 

  

Behaviour is undefined unless ``solve`` has been called before. 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.add_variables(2) 

1 

sage: p.add_linear_constraint([[0, 1], [1, 2]], None, 3) 

sage: p.set_objective([2, 5]) 

sage: p.solve() 

0 

sage: p.get_objective_value() 

7.5 

sage: p.get_variable_value(0) # abs tol 1e-15 

0.0 

sage: p.get_variable_value(1) 

1.5 

""" 

if (self.simplex_or_intopt != glp_simplex_only 

and self.simplex_or_intopt != glp_exact_simplex_only): 

return glp_mip_obj_val(self.lp) 

else: 

return glp_get_obj_val(self.lp) 

  

cpdef best_known_objective_bound(self): 

r""" 

Return the value of the currently best known bound. 

  

This method returns the current best upper (resp. lower) bound on the 

optimal value of the objective function in a maximization 

(resp. minimization) problem. It is equal to the output of 

:meth:`get_objective_value` if the MILP found an optimal solution, but 

it can differ if it was interrupted manually or after a time limit (cf 

:meth:`solver_parameter`). 

  

.. NOTE:: 

  

Has no meaning unless ``solve`` has been called before. 

  

EXAMPLES:: 

  

sage: g = graphs.CubeGraph(9) 

sage: p = MixedIntegerLinearProgram(solver="GLPK") 

sage: p.solver_parameter("mip_gap_tolerance",100) 

sage: b = p.new_variable(binary=True) 

sage: p.set_objective(p.sum(b[v] for v in g)) 

sage: for v in g: 

....: p.add_constraint(b[v]+p.sum(b[u] for u in g.neighbors(v)) <= 1) 

sage: p.add_constraint(b[v] == 1) # Force an easy non-0 solution 

sage: p.solve() # rel tol 100 

1.0 

sage: backend = p.get_backend() 

sage: backend.best_known_objective_bound() # random 

48.0 

""" 

return self.search_tree_data.best_bound 

  

cpdef get_relative_objective_gap(self): 

r""" 

Return the relative objective gap of the best known solution. 

  

For a minimization problem, this value is computed by 

`(\texttt{bestinteger} - \texttt{bestobjective}) / (1e-10 + 

|\texttt{bestobjective}|)`, where ``bestinteger`` is the value returned 

by :meth:`get_objective_value` and ``bestobjective`` is the value 

returned by :meth:`best_known_objective_bound`. For a maximization 

problem, the value is computed by `(\texttt{bestobjective} - 

\texttt{bestinteger}) / (1e-10 + |\texttt{bestobjective}|)`. 

  

.. NOTE:: 

  

Has no meaning unless ``solve`` has been called before. 

  

EXAMPLES:: 

  

sage: g = graphs.CubeGraph(9) 

sage: p = MixedIntegerLinearProgram(solver="GLPK") 

sage: p.solver_parameter("mip_gap_tolerance",100) 

sage: b = p.new_variable(binary=True) 

sage: p.set_objective(p.sum(b[v] for v in g)) 

sage: for v in g: 

....: p.add_constraint(b[v]+p.sum(b[u] for u in g.neighbors(v)) <= 1) 

sage: p.add_constraint(b[v] == 1) # Force an easy non-0 solution 

sage: p.solve() # rel tol 100 

1.0 

sage: backend = p.get_backend() 

sage: backend.get_relative_objective_gap() # random 

46.99999999999999 

  

TESTS: 

  

Just make sure that the variable *has* been defined, and is not just 

undefined:: 

  

sage: backend.get_relative_objective_gap() > 1 

True 

""" 

return self.search_tree_data.mip_gap 

  

cpdef get_variable_value(self, int variable): 

""" 

Returns the value of a variable given by the solver. 

  

.. NOTE:: 

  

Behaviour is undefined unless ``solve`` has been called before. 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.add_variables(2) 

1 

sage: p.add_linear_constraint([[0, 1], [1, 2]], None, 3) 

sage: p.set_objective([2, 5]) 

sage: p.solve() 

0 

sage: p.get_objective_value() 

7.5 

sage: p.get_variable_value(0) # abs tol 1e-15 

0.0 

sage: p.get_variable_value(1) 

1.5 

""" 

if (self.simplex_or_intopt != glp_simplex_only 

and self.simplex_or_intopt != glp_exact_simplex_only): 

return glp_mip_col_val(self.lp, variable+1) 

else: 

return glp_get_col_prim(self.lp, variable+1) 

  

cpdef get_row_prim(self, int i): 

r""" 

Returns the value of the auxiliary variable associated with i-th row. 

  

.. NOTE:: 

  

Behaviour is undefined unless ``solve`` has been called before. 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: lp = get_solver(solver = "GLPK") 

sage: lp.add_variables(3) 

2 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [8, 6, 1])), None, 48) 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [4, 2, 1.5])), None, 20) 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [2, 1.5, 0.5])), None, 8) 

sage: lp.set_objective([60, 30, 20]) 

sage: import sage.numerical.backends.glpk_backend as backend 

sage: lp.solver_parameter(backend.glp_simplex_or_intopt, backend.glp_simplex_only) 

sage: lp.solve() 

0 

sage: lp.get_objective_value() 

280.0 

sage: lp.get_row_prim(0) 

24.0 

sage: lp.get_row_prim(1) 

20.0 

sage: lp.get_row_prim(2) 

8.0 

""" 

return glp_get_row_prim(self.lp, i+1) 

  

cpdef int ncols(self): 

""" 

Return the number of columns/variables. 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.ncols() 

0 

sage: p.add_variables(2) 

1 

sage: p.ncols() 

2 

""" 

return glp_get_num_cols(self.lp) 

  

cpdef int nrows(self): 

""" 

Return the number of rows/constraints. 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.nrows() 

0 

sage: p.add_linear_constraints(2, 2, None) 

sage: p.nrows() 

2 

""" 

  

return glp_get_num_rows(self.lp) 

  

cpdef col_name(self, int index): 

""" 

Return the ``index`` th col name 

  

INPUT: 

  

- ``index`` (integer) -- the col's id 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.add_variable(name='I am a variable') 

0 

sage: p.col_name(0) 

'I am a variable' 

""" 

cdef char * s 

  

glp_create_index(self.lp) 

s = <char*> glp_get_col_name(self.lp, index + 1) 

  

if s != NULL: 

return s 

else: 

return "" 

  

cpdef row_name(self, int index): 

""" 

Return the ``index`` th row name 

  

INPUT: 

  

- ``index`` (integer) -- the row's id 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.add_linear_constraints(1, 2, None, names=['Empty constraint 1']) 

sage: p.row_name(0) 

'Empty constraint 1' 

""" 

cdef char * s 

  

glp_create_index(self.lp) 

s = <char*> glp_get_row_name(self.lp, index + 1) 

  

if s != NULL: 

return s 

else: 

return "" 

  

cpdef bint is_variable_binary(self, int index): 

""" 

Test whether the given variable is of binary type. 

  

INPUT: 

  

- ``index`` (integer) -- the variable's id 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.ncols() 

0 

sage: p.add_variable() 

0 

sage: p.set_variable_type(0,0) 

sage: p.is_variable_binary(0) 

True 

  

""" 

return glp_get_col_kind(self.lp, index + 1) == GLP_BV 

  

cpdef bint is_variable_integer(self, int index): 

""" 

Test whether the given variable is of integer type. 

  

INPUT: 

  

- ``index`` (integer) -- the variable's id 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.ncols() 

0 

sage: p.add_variable() 

0 

sage: p.set_variable_type(0,1) 

sage: p.is_variable_integer(0) 

True 

""" 

return glp_get_col_kind(self.lp, index + 1) == GLP_IV 

  

cpdef bint is_variable_continuous(self, int index): 

""" 

Test whether the given variable is of continuous/real type. 

  

INPUT: 

  

- ``index`` (integer) -- the variable's id 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.ncols() 

0 

sage: p.add_variable() 

0 

sage: p.is_variable_continuous(0) 

True 

sage: p.set_variable_type(0,1) 

sage: p.is_variable_continuous(0) 

False 

  

""" 

return glp_get_col_kind(self.lp, index + 1) == GLP_CV 

  

cpdef bint is_maximization(self): 

""" 

Test whether the problem is a maximization 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.is_maximization() 

True 

sage: p.set_sense(-1) 

sage: p.is_maximization() 

False 

""" 

  

return glp_get_obj_dir(self.lp) == GLP_MAX 

  

cpdef variable_upper_bound(self, int index, value = False): 

""" 

Return or define the upper bound on a variable 

  

INPUT: 

  

- ``index`` (integer) -- the variable's id 

  

- ``value`` -- real value, or ``None`` to mean that the 

variable has not upper bound. When set to ``False`` 

(default), the method returns the current value. 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.add_variable() 

0 

sage: p.col_bounds(0) 

(0.0, None) 

sage: p.variable_upper_bound(0, 5) 

sage: p.col_bounds(0) 

(0.0, 5.0) 

  

TESTS: 

  

:trac:`14581`:: 

  

sage: P = MixedIntegerLinearProgram(solver="GLPK") 

sage: x = P["x"] 

sage: P.set_max(x, 0) 

sage: P.get_max(x) 

0.0 

  

Check that :trac:`10232` is fixed:: 

  

sage: p = get_solver(solver="GLPK") 

sage: p.variable_upper_bound(2) 

Traceback (most recent call last): 

... 

GLPKError: ... 

sage: p.variable_upper_bound(3, 5) 

Traceback (most recent call last): 

... 

GLPKError: ... 

  

sage: p.add_variable() 

0 

sage: p.variable_upper_bound(0, 'hey!') 

Traceback (most recent call last): 

... 

TypeError: a float is required 

""" 

cdef double x 

cdef double min 

cdef double dvalue 

  

if value is False: 

sig_on() 

x = glp_get_col_ub(self.lp, index +1) 

sig_off() 

if x == DBL_MAX: 

return None 

else: 

return x 

else: 

sig_on() 

min = glp_get_col_lb(self.lp, index + 1) 

sig_off() 

  

if value is None: 

sig_on() 

if min == -DBL_MAX: 

glp_set_col_bnds(self.lp, index + 1, GLP_FR, 0, 0) 

else: 

glp_set_col_bnds(self.lp, index + 1, GLP_LO, min, 0) 

sig_off() 

else: 

dvalue = <double?> value 

  

sig_on() 

if min == -DBL_MAX: 

glp_set_col_bnds(self.lp, index + 1, GLP_UP, 0, dvalue) 

  

elif min == dvalue: 

glp_set_col_bnds(self.lp, index + 1, GLP_FX, dvalue, dvalue) 

else: 

glp_set_col_bnds(self.lp, index + 1, GLP_DB, min, dvalue) 

sig_off() 

  

cpdef variable_lower_bound(self, int index, value = False): 

""" 

Return or define the lower bound on a variable 

  

INPUT: 

  

- ``index`` (integer) -- the variable's id 

  

- ``value`` -- real value, or ``None`` to mean that the 

variable has not lower bound. When set to ``False`` 

(default), the method returns the current value. 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.add_variable() 

0 

sage: p.col_bounds(0) 

(0.0, None) 

sage: p.variable_lower_bound(0, 5) 

sage: p.col_bounds(0) 

(5.0, None) 

  

TESTS: 

  

:trac:`14581`:: 

  

sage: P = MixedIntegerLinearProgram(solver="GLPK") 

sage: x = P["x"] 

sage: P.set_min(x, 5) 

sage: P.set_min(x, 0) 

sage: P.get_min(x) 

0.0 

  

Check that :trac:`10232` is fixed:: 

  

sage: p = get_solver(solver="GLPK") 

sage: p.variable_lower_bound(2) 

Traceback (most recent call last): 

... 

GLPKError: ... 

sage: p.variable_lower_bound(3, 5) 

Traceback (most recent call last): 

... 

GLPKError: ... 

  

sage: p.add_variable() 

0 

sage: p.variable_lower_bound(0, 'hey!') 

Traceback (most recent call last): 

... 

TypeError: a float is required 

""" 

cdef double x 

cdef double max 

cdef double dvalue 

  

if value is False: 

sig_on() 

x = glp_get_col_lb(self.lp, index +1) 

sig_off() 

if x == -DBL_MAX: 

return None 

else: 

return x 

else: 

sig_on() 

max = glp_get_col_ub(self.lp, index + 1) 

sig_off() 

  

if value is None: 

sig_on() 

if max == DBL_MAX: 

glp_set_col_bnds(self.lp, index + 1, GLP_FR, 0.0, 0.0) 

else: 

glp_set_col_bnds(self.lp, index + 1, GLP_UP, 0.0, max) 

sig_off() 

  

else: 

dvalue = <double?> value 

  

sig_on() 

if max == DBL_MAX: 

glp_set_col_bnds(self.lp, index + 1, GLP_LO, value, 0.0) 

elif max == value: 

glp_set_col_bnds(self.lp, index + 1, GLP_FX, value, value) 

else: 

glp_set_col_bnds(self.lp, index + 1, GLP_DB, value, max) 

sig_off() 

  

cpdef write_lp(self, char * filename): 

""" 

Write the problem to a .lp file 

  

INPUT: 

  

- ``filename`` (string) 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.add_variables(2) 

1 

sage: p.add_linear_constraint([[0, 1], [1, 2]], None, 3) 

sage: p.set_objective([2, 5]) 

sage: p.write_lp(os.path.join(SAGE_TMP, "lp_problem.lp")) 

Writing problem data to ... 

9 lines were written 

""" 

glp_write_lp(self.lp, NULL, filename) 

  

cpdef write_mps(self, char * filename, int modern): 

""" 

Write the problem to a .mps file 

  

INPUT: 

  

- ``filename`` (string) 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.add_variables(2) 

1 

sage: p.add_linear_constraint([[0, 1], [1, 2]], None, 3) 

sage: p.set_objective([2, 5]) 

sage: p.write_mps(os.path.join(SAGE_TMP, "lp_problem.mps"), 2) 

Writing problem data to ... 

17 records were written 

""" 

glp_write_mps(self.lp, modern, NULL, filename) 

  

cpdef __copy__(self): 

""" 

Returns a copy of self. 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = MixedIntegerLinearProgram(solver = "GLPK") 

sage: b = p.new_variable() 

sage: p.add_constraint(b[1] + b[2] <= 6) 

sage: p.set_objective(b[1] + b[2]) 

sage: copy(p).solve() 

6.0 

""" 

cdef GLPKBackend p = type(self)(maximization = (1 if self.is_maximization() else -1)) 

p.simplex_or_intopt = self.simplex_or_intopt 

p.iocp.tm_lim = self.iocp.tm_lim 

glp_copy_prob(p.lp, self.lp, 1) 

return p 

  

  

cpdef solver_parameter(self, name, value = None): 

""" 

Return or define a solver parameter 

  

INPUT: 

  

- ``name`` (string) -- the parameter 

  

- ``value`` -- the parameter's value if it is to be defined, 

or ``None`` (default) to obtain its current value. 

  

You can supply the name of a parameter and its value using either a 

string or a ``glp_`` constant (which are defined as Cython variables of 

this module). 

  

In most cases, you can use the same name for a parameter as that given 

in the GLPK documentation, which is available by downloading GLPK from 

<http://www.gnu.org/software/glpk/>. The exceptions relate to parameters 

common to both methods; these require you to append ``_simplex`` or 

``_intopt`` to the name to resolve ambiguity, since the interface allows 

access to both. 

  

We have also provided more meaningful names, to assist readability. 

  

Parameter **names** are specified in lower case. 

To use a constant instead of a string, prepend ``glp_`` to the name. 

For example, both ``glp_gmi_cuts`` or ``"gmi_cuts"`` control whether 

to solve using Gomory cuts. 

  

Parameter **values** are specified as strings in upper case, 

or as constants in lower case. For example, both ``glp_on`` and ``"GLP_ON"`` 

specify the same thing. 

  

Naturally, you can use ``True`` and ``False`` in cases where ``glp_on`` and ``glp_off`` 

would be used. 

  

A list of parameter names, with their possible values: 

  

**General-purpose parameters:** 

  

.. list-table:: 

:widths: 30 70 

  

* - ``timelimit`` 

  

- specify the time limit IN SECONDS. This affects both simplex 

and intopt. 

  

* - ``timelimit_simplex`` and ``timelimit_intopt`` 

  

- specify the time limit IN MILLISECONDS. (This is glpk's 

default.) 

  

* - ``simplex_or_intopt`` 

  

- specify which of ``simplex``, ``exact`` and ``intopt`` routines 

in GLPK to use. 

This is controlled by setting ``simplex_or_intopt`` to 

``glp_simplex_only``, ``glp_exact_simplex_only``, 

``glp_intopt_only`` and ``glp_simplex_then_intopt``, respectively. 

The latter is useful to deal with a problem in GLPK where 

problems with no solution hang when using integer optimization; 

if you specify ``glp_simplex_then_intopt``, 

sage will try simplex first, then perform integer optimization 

only if a solution of the LP relaxation exists. 

  

* - ``verbosity_intopt`` and ``verbosity_simplex`` 

  

- one of ``GLP_MSG_OFF``, ``GLP_MSG_ERR``, ``GLP_MSG_ON``, or 

``GLP_MSG_ALL``. The default is ``GLP_MSG_OFF``. 

  

* - ``output_frequency_intopt`` and ``output_frequency_simplex`` 

  

- the output frequency, in milliseconds. Default is 5000. 

  

* - ``output_delay_intopt`` and ``output_delay_simplex`` 

  

- the output delay, in milliseconds, regarding the use of the 

simplex method on the LP relaxation. Default is 10000. 

  

  

**intopt-specific parameters:** 

  

.. list-table:: 

:widths: 30 70 

  

* - ``branching`` 

- - ``GLP_BR_FFV`` first fractional variable 

- ``GLP_BR_LFV`` last fractional variable 

- ``GLP_BR_MFV`` most fractional variable 

- ``GLP_BR_DTH`` Driebeck-Tomlin heuristic (default) 

- ``GLP_BR_PCH`` hybrid pseudocost heuristic 

  

* - ``backtracking`` 

- - ``GLP_BT_DFS`` depth first search 

- ``GLP_BT_BFS`` breadth first search 

- ``GLP_BT_BLB`` best local bound (default) 

- ``GLP_BT_BPH`` best projection heuristic 

  

* - ``preprocessing`` 

- - ``GLP_PP_NONE`` 

- ``GLP_PP_ROOT`` preprocessing only at root level 

- ``GLP_PP_ALL`` (default) 

  

  

* - ``feasibility_pump`` 

  

- ``GLP_ON`` or ``GLP_OFF`` (default) 

  

* - ``gomory_cuts`` 

  

- ``GLP_ON`` or ``GLP_OFF`` (default) 

  

* - ``mixed_int_rounding_cuts`` 

  

- ``GLP_ON`` or ``GLP_OFF`` (default) 

  

* - ``mixed_cover_cuts`` 

  

- ``GLP_ON`` or ``GLP_OFF`` (default) 

  

* - ``clique_cuts`` 

  

- ``GLP_ON`` or ``GLP_OFF`` (default) 

  

* - ``absolute_tolerance`` 

  

- (double) used to check if optimal solution to LP relaxation is 

integer feasible. GLPK manual advises, "do not change... without 

detailed understanding of its purpose." 

  

* - ``relative_tolerance`` 

  

- (double) used to check if objective value in LP relaxation is not 

better than best known integer solution. GLPK manual advises, "do 

not change... without detailed understanding of its purpose." 

  

* - ``mip_gap_tolerance`` 

  

- (double) relative mip gap tolerance. Default is 0.0. 

  

* - ``presolve_intopt`` 

  

- ``GLP_ON`` (default) or ``GLP_OFF``. 

  

* - ``binarize`` 

  

- ``GLP_ON`` or ``GLP_OFF`` (default) 

  

  

**simplex-specific parameters:** 

  

.. list-table:: 

:widths: 30 70 

  

* - ``primal_v_dual`` 

  

- - ``GLP_PRIMAL`` (default) 

- ``GLP_DUAL`` 

- ``GLP_DUALP`` 

  

* - ``pricing`` 

  

- - ``GLP_PT_STD`` standard (textbook) 

- ``GLP_PT_PSE`` projected steepest edge (default) 

  

* - ``ratio_test`` 

  

- - ``GLP_RT_STD`` standard (textbook) 

- ``GLP_RT_HAR`` Harris' two-pass ratio test (default) 

  

* - ``tolerance_primal`` 

  

- (double) tolerance used to check if basic solution is primal 

feasible. GLPK manual advises, "do not change... without 

detailed understanding of its purpose." 

  

* - ``tolerance_dual`` 

  

- (double) tolerance used to check if basic solution is dual 

feasible. GLPK manual advises, "do not change... without 

detailed understanding of its purpose." 

  

* - ``tolerance_pivot`` 

  

- (double) tolerance used to choose pivot. GLPK manual advises, 

"do not change... without detailed understanding of its 

purpose." 

  

* - ``obj_lower_limit`` 

  

- (double) lower limit of the objective function. The default is 

``-DBL_MAX``. 

  

* - ``obj_upper_limit`` 

  

- (double) upper limit of the objective function. The default is 

``DBL_MAX``. 

  

* - ``iteration_limit`` 

  

- (int) iteration limit of the simplex algorithm. The default is 

``INT_MAX``. 

  

* - ``presolve_simplex`` 

  

- ``GLP_ON`` or ``GLP_OFF`` (default). 

  

.. NOTE:: 

  

The coverage for GLPK's control parameters for simplex and integer optimization 

is nearly complete. The only thing lacking is a wrapper for callback routines. 

  

To date, no attempt has been made to expose the interior point methods. 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.solver_parameter("timelimit", 60) 

sage: p.solver_parameter("timelimit") 

60.0 

  

- Don't forget the difference between ``timelimit`` and ``timelimit_intopt`` 

  

:: 

  

sage: p.solver_parameter("timelimit_intopt") 

60000 

  

If you don't care for an integer answer, you can ask for an LP 

relaxation instead. The default solver performs integer optimization, 

but you can switch to the standard simplex algorithm through the 

``glp_simplex_or_intopt`` parameter. 

  

EXAMPLES:: 

  

sage: lp = MixedIntegerLinearProgram(solver = 'GLPK', maximization = False) 

sage: x, y = lp[0], lp[1] 

sage: lp.add_constraint(-2*x + y <= 1) 

sage: lp.add_constraint(x - y <= 1) 

sage: lp.add_constraint(x + y >= 2) 

sage: lp.set_integer(x); lp.set_integer(y) 

sage: lp.set_objective(x + y) 

sage: lp.solve() 

2.0 

sage: lp.get_values([x,y]) 

[1.0, 1.0] 

sage: import sage.numerical.backends.glpk_backend as backend 

sage: lp.solver_parameter(backend.glp_simplex_or_intopt, backend.glp_simplex_only) 

sage: lp.solve() 

2.0 

sage: lp.get_values([x,y]) 

[1.5, 0.5] 

  

You can get GLPK to spout all sorts of information at you. 

The default is to turn this off, but sometimes (debugging) it's very useful:: 

  

sage: lp.solver_parameter(backend.glp_simplex_or_intopt, backend.glp_simplex_then_intopt) 

sage: lp.solver_parameter(backend.glp_mir_cuts, backend.glp_on) 

sage: lp.solver_parameter(backend.glp_msg_lev_intopt, backend.glp_msg_all) 

sage: lp.solver_parameter(backend.glp_mir_cuts) 

1 

  

If you actually try to solve ``lp``, you will get a lot of detailed information. 

""" 

  

if type(name) == str: 

if name == "out_frq" or name == "out_dly" or name == "tm_lim" or name == "msg_lev": 

raise ValueError("To set parameter " + name + " you must specify the solver. Append either _simplex or _intopt.") 

name = solver_parameter_names_dict[name] 

  

if type(value) == str: value = solver_parameter_values[value] 

  

if name == timelimit_intopt: 

if value is None: return self.iocp.tm_lim 

else: self.iocp.tm_lim = value 

  

if name == timelimit_seconds: 

if value is None: return self.iocp.tm_lim / 1000.0 

else: 

self.iocp.tm_lim = value * 1000.0 

self.smcp.tm_lim = value * 1000.0 

  

elif name == timelimit_simplex: 

if value is None: return self.smcp.tm_lim 

else: self.smcp.tm_lim = value 

  

elif name == simplex_or_intopt: 

if value is None: return self.simplex_or_intopt 

if not value in (simplex_only,intopt_only,simplex_then_intopt,exact_simplex_only): 

raise MIPSolverException("GLPK: invalid value for simplex_or_intopt; see documentation") 

self.simplex_or_intopt = value 

  

elif name == msg_lev_simplex: 

if value is None: return self.smcp.msg_lev 

else: self.smcp.msg_lev = value 

  

elif name == msg_lev_intopt: 

if value is None: return self.iocp.msg_lev 

else: self.iocp.msg_lev = value 

  

elif name == br_tech: 

if value is None: return self.iocp.br_tech 

else: self.iocp.br_tech = value 

  

elif name == bt_tech: 

if value is None: return self.iocp.bt_tech 

else: self.iocp.bt_tech = value 

  

elif name == pp_tech: 

if value is None: return self.iocp.pp_tech 

else: self.iocp.pp_tech = value 

  

elif name == fp_heur: 

if value is None: return self.iocp.fp_heur 

else: 

if value: value = GLP_ON 

else: value = GLP_OFF 

self.iocp.fp_heur = value 

  

elif name == gmi_cuts: 

if value is None: return self.iocp.gmi_cuts 

else: 

if value: value = GLP_ON 

else: value = GLP_OFF 

self.iocp.gmi_cuts = value 

  

elif name == mir_cuts: 

if value is None: return self.iocp.mir_cuts 

else: 

if value: value = GLP_ON 

else: value = GLP_OFF 

self.iocp.mir_cuts = value 

  

elif name == cov_cuts: 

if value is None: return self.iocp.cov_cuts 

else: 

if value: value = GLP_ON 

else: value = GLP_OFF 

self.iocp.cov_cuts = value 

  

elif name == clq_cuts: 

if value is None: return self.iocp.clq_cuts 

else: 

if value: value = GLP_ON 

else: value = GLP_OFF 

self.iocp.clq_cuts = value 

  

elif name == tol_int: 

if value is None: return self.iocp.tol_int 

else: self.iocp.tol_int = value 

  

elif name == tol_obj: 

if value is None: return self.iocp.tol_obj 

else: self.iocp.tol_obj = value 

  

elif name == mip_gap: 

if value is None: return self.iocp.mip_gap 

else: self.iocp.mip_gap = value 

  

elif name == tm_lim_intopt: 

if value is None: return self.iocp.tm_lim 

else: self.iocp.tm_lim = value 

  

elif name == out_frq_intopt: 

if value is None: return self.iocp.out_frq 

else: self.iocp.out_frq = value 

  

elif name == out_dly_intopt: 

if value is None: return self.iocp.out_dly 

else: self.iocp.out_dly = value 

  

elif name == presolve_intopt: 

if value is None: return self.iocp.presolve 

else: 

if value: value = GLP_ON 

else: value = GLP_OFF 

self.iocp.presolve = value 

  

elif name == binarize: 

if value is None: return self.iocp.binarize 

else: 

if value: value = GLP_ON 

else: value = GLP_OFF 

self.iocp.binarize = value 

  

elif name == msg_lev_simplex: 

if value is None: return self.smcp.msg_lev 

else: self.smcp.msg_lev = value 

  

elif name == meth: 

if value is None: return self.smcp.meth 

else: self.smcp.meth = value 

  

elif name == pricing: 

if value is None: return self.smcp.pricing 

else: self.smcp.pricing = value 

  

elif name == r_test: 

if value is None: return self.smcp.r_test 

else: self.smcp.r_test = value 

  

elif name == tol_bnd: 

if value is None: return self.smcp.tol_bnd 

else: self.smcp.tol_bnd = value 

  

elif name == tol_dj: 

if value is None: return self.smcp.tol_dj 

else: self.smcp.tol_dj = value 

  

elif name == tol_piv: 

if value is None: return self.smcp.tol_piv 

else: self.smcp.tol_piv = value 

  

elif name == obj_ll: 

if value is None: return self.smcp.obj_ll 

else: self.smcp.obj_ll = value 

  

elif name == obj_ul: 

if value is None: return self.smcp.obj_ul 

else: self.smcp.obj_ul = value 

  

elif name == it_lim: 

if value is None: return self.smcp.it_lim 

else: self.smcp.it_lim = value 

  

elif name == tm_lim_simplex: 

if value is None: return self.smcp.tm_lim 

else: self.smcp.tm_lim = value 

  

elif name == out_frq_simplex: 

if value is None: return self.smcp.out_frq 

else: self.smcp.out_frq = value 

  

elif name == out_dly_simplex: 

if value is None: return self.smcp.out_dly 

else: self.smcp.out_dly = value 

  

elif name == presolve_simplex: 

if value is None: return self.smcp.presolve 

else: 

if value: value = GLP_ON 

else: value = GLP_OFF 

self.smcp.presolve = value 

  

else: 

raise ValueError("This parameter is not available.") 

  

  

cpdef bint is_variable_basic(self, int index): 

""" 

Test whether the given variable is basic. 

  

This assumes that the problem has been solved with the simplex method 

and a basis is available. Otherwise an exception will be raised. 

  

INPUT: 

  

- ``index`` (integer) -- the variable's id 

  

EXAMPLES:: 

  

sage: p = MixedIntegerLinearProgram(maximization=True,\ 

solver="GLPK") 

sage: x = p.new_variable(nonnegative=True) 

sage: p.add_constraint(-x[0] + x[1] <= 2) 

sage: p.add_constraint(8 * x[0] + 2 * x[1] <= 17) 

sage: p.set_objective(5.5 * x[0] - 3 * x[1]) 

sage: b = p.get_backend() 

sage: import sage.numerical.backends.glpk_backend as backend 

sage: b.solver_parameter(backend.glp_simplex_or_intopt, backend.glp_simplex_only) 

sage: b.solve() 

0 

sage: b.is_variable_basic(0) 

True 

sage: b.is_variable_basic(1) 

False 

""" 

return self.get_col_stat(index) == GLP_BS 

  

cpdef bint is_variable_nonbasic_at_lower_bound(self, int index): 

""" 

Test whether the given variable is nonbasic at lower bound. 

This assumes that the problem has been solved with the simplex method 

and a basis is available. Otherwise an exception will be raised. 

  

INPUT: 

  

- ``index`` (integer) -- the variable's id 

  

EXAMPLES:: 

  

sage: p = MixedIntegerLinearProgram(maximization=True,\ 

solver="GLPK") 

sage: x = p.new_variable(nonnegative=True) 

sage: p.add_constraint(-x[0] + x[1] <= 2) 

sage: p.add_constraint(8 * x[0] + 2 * x[1] <= 17) 

sage: p.set_objective(5.5 * x[0] - 3 * x[1]) 

sage: b = p.get_backend() 

sage: import sage.numerical.backends.glpk_backend as backend 

sage: b.solver_parameter(backend.glp_simplex_or_intopt, backend.glp_simplex_only) 

sage: b.solve() 

0 

sage: b.is_variable_nonbasic_at_lower_bound(0) 

False 

sage: b.is_variable_nonbasic_at_lower_bound(1) 

True 

""" 

return self.get_col_stat(index) == GLP_NL 

  

cpdef bint is_slack_variable_basic(self, int index): 

""" 

Test whether the slack variable of the given row is basic. 

  

This assumes that the problem has been solved with the simplex method 

and a basis is available. Otherwise an exception will be raised. 

  

INPUT: 

  

- ``index`` (integer) -- the variable's id 

  

EXAMPLES:: 

  

sage: p = MixedIntegerLinearProgram(maximization=True,\ 

solver="GLPK") 

sage: x = p.new_variable(nonnegative=True) 

sage: p.add_constraint(-x[0] + x[1] <= 2) 

sage: p.add_constraint(8 * x[0] + 2 * x[1] <= 17) 

sage: p.set_objective(5.5 * x[0] - 3 * x[1]) 

sage: b = p.get_backend() 

sage: import sage.numerical.backends.glpk_backend as backend 

sage: b.solver_parameter(backend.glp_simplex_or_intopt, backend.glp_simplex_only) 

sage: b.solve() 

0 

sage: b.is_slack_variable_basic(0) 

True 

sage: b.is_slack_variable_basic(1) 

False 

""" 

return self.get_row_stat(index) == GLP_BS 

  

cpdef bint is_slack_variable_nonbasic_at_lower_bound(self, int index): 

""" 

Test whether the slack variable of the given row is nonbasic at lower bound. 

  

This assumes that the problem has been solved with the simplex method 

and a basis is available. Otherwise an exception will be raised. 

  

INPUT: 

  

- ``index`` (integer) -- the variable's id 

  

EXAMPLES:: 

  

sage: p = MixedIntegerLinearProgram(maximization=True,\ 

solver="GLPK") 

sage: x = p.new_variable(nonnegative=True) 

sage: p.add_constraint(-x[0] + x[1] <= 2) 

sage: p.add_constraint(8 * x[0] + 2 * x[1] <= 17) 

sage: p.set_objective(5.5 * x[0] - 3 * x[1]) 

sage: b = p.get_backend() 

sage: import sage.numerical.backends.glpk_backend as backend 

sage: b.solver_parameter(backend.glp_simplex_or_intopt, backend.glp_simplex_only) 

sage: b.solve() 

0 

sage: b.is_slack_variable_nonbasic_at_lower_bound(0) 

False 

sage: b.is_slack_variable_nonbasic_at_lower_bound(1) 

True 

""" 

return self.get_row_stat(index) == GLP_NU 

  

cpdef int print_ranges(self, char * filename = NULL) except -1: 

r""" 

Print results of a sensitivity analysis 

  

If no filename is given as an input the results of the 

sensitivity analysis are displayed on the screen. If a 

filename is given they are written to a file. 

  

INPUT: 

  

- ``filename`` -- (optional) name of the file 

  

OUTPUT: 

  

Zero if the operations was successful otherwise nonzero. 

  

.. NOTE:: 

  

This method is only effective if an optimal solution has been found 

for the lp using the simplex algorithm. In all other cases an error 

message is printed. 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.add_variables(2) 

1 

sage: p.add_linear_constraint(list(zip([0, 1], [1, 2])), None, 3) 

sage: p.set_objective([2, 5]) 

sage: import sage.numerical.backends.glpk_backend as backend 

sage: p.solver_parameter(backend.glp_simplex_or_intopt, backend.glp_simplex_only) 

sage: p.print_ranges() 

glp_print_ranges: optimal basic solution required 

1 

sage: p.solve() 

0 

sage: p.print_ranges() 

Write sensitivity analysis report to ... 

GLPK ... - SENSITIVITY ANALYSIS REPORT Page 1 

<BLANKLINE> 

Problem: 

Objective: 7.5 (MAXimum) 

<BLANKLINE> 

No. Row name St Activity Slack Lower bound Activity Obj coef Obj value at Limiting 

Marginal Upper bound range range break point variable 

------ ------------ -- ------------- ------------- ------------- ------------- ------------- ------------- ------------ 

1 NU 3.00000 . -Inf . -2.50000 . 

2.50000 3.00000 +Inf +Inf +Inf 

<BLANKLINE> 

GLPK ... - SENSITIVITY ANALYSIS REPORT Page 2 

<BLANKLINE> 

Problem: 

Objective: 7.5 (MAXimum) 

<BLANKLINE> 

No. Column name St Activity Obj coef Lower bound Activity Obj coef Obj value at Limiting 

Marginal Upper bound range range break point variable 

------ ------------ -- ------------- ------------- ------------- ------------- ------------- ------------- ------------ 

1 NL . 2.00000 . -Inf -Inf +Inf 

-.50000 +Inf 3.00000 2.50000 6.00000 

<BLANKLINE> 

2 BS 1.50000 5.00000 . -Inf 4.00000 6.00000 

. +Inf 1.50000 +Inf +Inf 

<BLANKLINE> 

End of report 

<BLANKLINE> 

0 

""" 

  

from sage.misc.all import SAGE_TMP 

  

if filename == NULL: 

fname = SAGE_TMP+"/ranges.tmp" 

else: 

fname = filename 

  

res = glp_print_ranges(self.lp, 0, 0, 0, fname) 

  

if filename == NULL: 

if res == 0: 

with open(fname) as f: 

for line in f: 

print(line, end=" ") 

print("\n") 

  

return res 

  

cpdef double get_row_dual(self, int variable): 

r""" 

Returns the dual value of a constraint. 

  

The dual value of the ith row is also the value of the ith variable 

of the dual problem. 

  

The dual value of a constraint is the shadow price of the constraint. 

The shadow price is the amount by which the objective value will change 

if the constraints bounds change by one unit under the precondition 

that the basis remains the same. 

  

INPUT: 

  

- ``variable`` -- The number of the constraint 

  

.. NOTE:: 

  

Behaviour is undefined unless ``solve`` has been called before. 

If the simplex algorithm has not been used for solving 0.0 will 

be returned. 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: lp = get_solver(solver = "GLPK") 

sage: lp.add_variables(3) 

2 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [8, 6, 1])), None, 48) 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [4, 2, 1.5])), None, 20) 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [2, 1.5, 0.5])), None, 8) 

sage: lp.set_objective([60, 30, 20]) 

sage: import sage.numerical.backends.glpk_backend as backend 

sage: lp.solver_parameter(backend.glp_simplex_or_intopt, backend.glp_simplex_only) 

sage: lp.solve() 

0 

sage: lp.get_row_dual(0) # tolerance 0.00001 

0.0 

sage: lp.get_row_dual(1) # tolerance 0.00001 

10.0 

  

  

""" 

  

if self.simplex_or_intopt == simplex_only: 

return glp_get_row_dual(self.lp, variable+1) 

else: 

return 0.0 

  

cpdef double get_col_dual(self, int variable): 

""" 

Returns the dual value (reduced cost) of a variable 

  

The dual value is the reduced cost of a variable. 

The reduced cost is the amount by which the objective coefficient 

of a non basic variable has to change to become a basic variable. 

  

INPUT: 

  

- ``variable`` -- The number of the variable 

  

.. NOTE:: 

  

Behaviour is undefined unless ``solve`` has been called before. 

If the simplex algorithm has not been used for solving just a 

0.0 will be returned. 

  

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: p = get_solver(solver = "GLPK") 

sage: p.add_variables(3) 

2 

sage: p.add_linear_constraint(list(zip([0, 1, 2], [8, 6, 1])), None, 48) 

sage: p.add_linear_constraint(list(zip([0, 1, 2], [4, 2, 1.5])), None, 20) 

sage: p.add_linear_constraint(list(zip([0, 1, 2], [2, 1.5, 0.5])), None, 8) 

sage: p.set_objective([60, 30, 20]) 

sage: import sage.numerical.backends.glpk_backend as backend 

sage: p.solver_parameter(backend.glp_simplex_or_intopt, backend.glp_simplex_only) 

sage: p.solve() 

0 

sage: p.get_col_dual(1) 

-5.0 

  

""" 

if self.simplex_or_intopt == simplex_only: 

return glp_get_col_dual(self.lp, variable+1) 

else: 

return 0.0 

  

cpdef int get_row_stat(self, int i) except? -1: 

""" 

Retrieve the status of a constraint. 

  

INPUT: 

  

- ``i`` -- The index of the constraint 

  

OUTPUT: 

  

- Returns current status assigned to the auxiliary variable associated with i-th row: 

  

* GLP_BS = 1 basic variable 

* GLP_NL = 2 non-basic variable on lower bound 

* GLP_NU = 3 non-basic variable on upper bound 

* GLP_NF = 4 non-basic free (unbounded) variable 

* GLP_NS = 5 non-basic fixed variable 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: lp = get_solver(solver = "GLPK") 

sage: lp.add_variables(3) 

2 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [8, 6, 1])), None, 48) 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [4, 2, 1.5])), None, 20) 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [2, 1.5, 0.5])), None, 8) 

sage: lp.set_objective([60, 30, 20]) 

sage: import sage.numerical.backends.glpk_backend as backend 

sage: lp.solver_parameter(backend.glp_simplex_or_intopt, backend.glp_simplex_only) 

sage: lp.solve() 

0 

sage: lp.get_row_stat(0) 

1 

sage: lp.get_row_stat(1) 

3 

sage: lp.get_row_stat(-1) 

Traceback (most recent call last): 

... 

ValueError: The constraint's index i must satisfy 0 <= i < number_of_constraints 

""" 

if i < 0 or i >= glp_get_num_rows(self.lp): 

raise ValueError("The constraint's index i must satisfy 0 <= i < number_of_constraints") 

return glp_get_row_stat(self.lp, i+1) 

  

cpdef int get_col_stat(self, int j) except? -1: 

""" 

Retrieve the status of a variable. 

  

INPUT: 

  

- ``j`` -- The index of the variable 

  

OUTPUT: 

  

- Returns current status assigned to the structural variable associated with j-th column: 

  

* GLP_BS = 1 basic variable 

* GLP_NL = 2 non-basic variable on lower bound 

* GLP_NU = 3 non-basic variable on upper bound 

* GLP_NF = 4 non-basic free (unbounded) variable 

* GLP_NS = 5 non-basic fixed variable 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: lp = get_solver(solver = "GLPK") 

sage: lp.add_variables(3) 

2 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [8, 6, 1])), None, 48) 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [4, 2, 1.5])), None, 20) 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [2, 1.5, 0.5])), None, 8) 

sage: lp.set_objective([60, 30, 20]) 

sage: import sage.numerical.backends.glpk_backend as backend 

sage: lp.solver_parameter(backend.glp_simplex_or_intopt, backend.glp_simplex_only) 

sage: lp.solve() 

0 

sage: lp.get_col_stat(0) 

1 

sage: lp.get_col_stat(1) 

2 

sage: lp.get_col_stat(100) 

Traceback (most recent call last): 

... 

ValueError: The variable's index j must satisfy 0 <= j < number_of_variables 

""" 

if j < 0 or j >= glp_get_num_cols(self.lp): 

raise ValueError("The variable's index j must satisfy 0 <= j < number_of_variables") 

  

return glp_get_col_stat(self.lp, j+1) 

  

cpdef set_row_stat(self, int i, int stat): 

r""" 

Set the status of a constraint. 

  

INPUT: 

  

- ``i`` -- The index of the constraint 

  

- ``stat`` -- The status to set to 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: lp = get_solver(solver = "GLPK") 

sage: lp.add_variables(3) 

2 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [8, 6, 1])), None, 48) 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [4, 2, 1.5])), None, 20) 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [2, 1.5, 0.5])), None, 8) 

sage: lp.set_objective([60, 30, 20]) 

sage: import sage.numerical.backends.glpk_backend as backend 

sage: lp.solver_parameter(backend.glp_simplex_or_intopt, backend.glp_simplex_only) 

sage: lp.solve() 

0 

sage: lp.get_row_stat(0) 

1 

sage: lp.set_row_stat(0, 3) 

sage: lp.get_row_stat(0) 

3 

""" 

if i < 0 or i >= glp_get_num_rows(self.lp): 

raise ValueError("The constraint's index i must satisfy 0 <= i < number_of_constraints") 

  

glp_set_row_stat(self.lp, i+1, stat) 

  

cpdef set_col_stat(self, int j, int stat): 

r""" 

Set the status of a variable. 

  

INPUT: 

  

- ``j`` -- The index of the constraint 

  

- ``stat`` -- The status to set to 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: lp = get_solver(solver = "GLPK") 

sage: lp.add_variables(3) 

2 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [8, 6, 1])), None, 48) 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [4, 2, 1.5])), None, 20) 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [2, 1.5, 0.5])), None, 8) 

sage: lp.set_objective([60, 30, 20]) 

sage: import sage.numerical.backends.glpk_backend as backend 

sage: lp.solver_parameter(backend.glp_simplex_or_intopt, backend.glp_simplex_only) 

sage: lp.solve() 

0 

sage: lp.get_col_stat(0) 

1 

sage: lp.set_col_stat(0, 2) 

sage: lp.get_col_stat(0) 

2 

""" 

if j < 0 or j >= glp_get_num_cols(self.lp): 

raise ValueError("The variable's index j must satisfy 0 <= j < number_of_variables") 

  

glp_set_col_stat(self.lp, j+1, stat) 

 

cpdef int warm_up(self): 

r""" 

Warm up the basis using current statuses assigned to rows and cols. 

  

OUTPUT: 

  

- Returns the warming up status 

  

* 0 The operation has been successfully performed. 

* GLP_EBADB The basis matrix is invalid.  

* GLP_ESING The basis matrix is singular within the working precision. 

* GLP_ECOND The basis matrix is ill-conditioned. 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: lp = get_solver(solver = "GLPK") 

sage: lp.add_variables(3) 

2 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [8, 6, 1])), None, 48) 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [4, 2, 1.5])), None, 20) 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [2, 1.5, 0.5])), None, 8) 

sage: lp.set_objective([60, 30, 20]) 

sage: import sage.numerical.backends.glpk_backend as backend 

sage: lp.solver_parameter(backend.glp_simplex_or_intopt, backend.glp_simplex_only) 

sage: lp.solve() 

0 

sage: lp.get_objective_value() 

280.0 

sage: lp.set_row_stat(0,3) 

sage: lp.set_col_stat(1,1) 

sage: lp.warm_up() 

0 

""" 

return glp_warm_up(self.lp) 

  

cpdef eval_tab_row(self, int k): 

r""" 

Computes a row of the current simplex tableau. 

  

A row corresponds to some basic variable specified by the parameter 

``k`` as follows: 

  

- if `0 \leq k \leq m-1`, the basic variable is `k`-th auxiliary 

variable, 

  

- if `m \leq k \leq m+n-1`, the basic variable is `(k-m)`-th structual 

variable, 

  

where `m` is the number of rows and `n` is the number of columns in the 

specified problem object. 

  

.. NOTE:: 

  

The basis factorization must exist. 

Otherwise, a ``MIPSolverException`` will be raised. 

  

INPUT: 

  

- ``k`` (integer) -- the id of the basic variable. 

  

OUTPUT: 

  

A pair ``(indices, coeffs)`` where ``indices`` lists the 

entries whose coefficient is nonzero, and to which ``coeffs`` 

associates their coefficient in the computed row 

of the current simplex tableau. 

  

.. NOTE:: 

  

Elements in ``indices`` have the same sense as index ``k``. 

All these variables are non-basic by definition. 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: lp = get_solver(solver = "GLPK") 

sage: lp.add_variables(3) 

2 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [8, 6, 1])), None, 48) 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [4, 2, 1.5])), None, 20) 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [2, 1.5, 0.5])), None, 8) 

sage: lp.set_objective([60, 30, 20]) 

sage: import sage.numerical.backends.glpk_backend as backend 

sage: lp.solver_parameter(backend.glp_simplex_or_intopt, backend.glp_simplex_only) 

sage: lp.eval_tab_row(0) 

Traceback (most recent call last): 

... 

MIPSolverException: ... 

sage: lp.solve() 

0 

sage: lp.eval_tab_row(0) 

([1, 2, 4], [-2.0, 8.0, -2.0]) 

sage: lp.eval_tab_row(3) 

([1, 2, 4], [-0.5, 1.5, -1.25]) 

sage: lp.eval_tab_row(5) 

([1, 2, 4], [2.0, -4.0, 2.0]) 

sage: lp.eval_tab_row(1) 

Traceback (most recent call last): 

... 

MIPSolverException: ... 

sage: lp.eval_tab_row(-1) 

Traceback (most recent call last): 

... 

ValueError: ... 

  

""" 

cdef int n = glp_get_num_cols(self.lp) 

cdef int i,j 

  

if k < 0 or k >= n + glp_get_num_rows(self.lp): 

raise ValueError("k = %s; Variable number out of range" % k) 

  

cdef MemoryAllocator mem = MemoryAllocator() 

cdef int * c_indices = <int*>mem.allocarray(n+1, sizeof(int)) 

cdef double * c_values = <double*>mem.allocarray(n+1, sizeof(double)) 

  

try: 

sig_on() # to catch GLPKError 

i = glp_eval_tab_row(self.lp, k + 1, c_indices, c_values) 

sig_off() 

except GLPKError: 

raise MIPSolverException('GLPK: basis factorization does not exist; or variable must be basic') 

  

indices = [c_indices[j+1] - 1 for j in range(i)] 

values = [c_values[j+1] for j in range(i)] 

return (indices, values) 

  

cpdef eval_tab_col(self, int k): 

r""" 

Computes a column of the current simplex tableau. 

  

A (column) corresponds to some non-basic variable specified by the 

parameter ``k`` as follows: 

  

- if `0 \leq k \leq m-1`, the non-basic variable is `k`-th auxiliary 

variable, 

  

- if `m \leq k \leq m+n-1`, the non-basic variable is `(k-m)`-th 

structual variable, 

  

where `m` is the number of rows and `n` is the number of columns 

in the specified problem object. 

  

.. NOTE:: 

  

The basis factorization must exist. 

Otherwise a ``MIPSolverException`` will be raised. 

  

INPUT: 

  

- ``k`` (integer) -- the id of the non-basic variable. 

  

OUTPUT: 

  

A pair ``(indices, coeffs)`` where ``indices`` lists the 

entries whose coefficient is nonzero, and to which ``coeffs`` 

associates their coefficient in the computed column 

of the current simplex tableau. 

  

.. NOTE:: 

  

Elements in ``indices`` have the same sense as index `k`. 

All these variables are basic by definition. 

  

EXAMPLES:: 

  

sage: from sage.numerical.backends.generic_backend import get_solver 

sage: lp = get_solver(solver = "GLPK") 

sage: lp.add_variables(3) 

2 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [8, 6, 1])), None, 48) 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [4, 2, 1.5])), None, 20) 

sage: lp.add_linear_constraint(list(zip([0, 1, 2], [2, 1.5, 0.5])), None, 8) 

sage: lp.set_objective([60, 30, 20]) 

sage: import sage.numerical.backends.glpk_backend as backend 

sage: lp.solver_parameter(backend.glp_simplex_or_intopt, backend.glp_simplex_only) 

sage: lp.eval_tab_col(1) 

Traceback (most recent call last): 

... 

MIPSolverException: ... 

sage: lp.solve() 

0 

sage: lp.eval_tab_col(1) 

([0, 5, 3], [-2.0, 2.0, -0.5]) 

sage: lp.eval_tab_col(2) 

([0, 5, 3], [8.0, -4.0, 1.5]) 

sage: lp.eval_tab_col(4) 

([0, 5, 3], [-2.0, 2.0, -1.25]) 

sage: lp.eval_tab_col(0) 

Traceback (most recent call last): 

... 

MIPSolverException: ... 

sage: lp.eval_tab_col(-1) 

Traceback (most recent call last): 

... 

ValueError: ... 

  

""" 

cdef int m = glp_get_num_rows(self.lp) 

cdef int i,j 

  

if k < 0 or k >= m + glp_get_num_cols(self.lp): 

raise ValueError("k = %s; Variable number out of range" % k) 

  

cdef MemoryAllocator mem = MemoryAllocator() 

cdef int * c_indices = <int*>mem.allocarray(m+1, sizeof(int)) 

cdef double * c_values = <double*>mem.allocarray(m+1, sizeof(double)) 

  

try: 

sig_on() # To catch GLPKError 

i = glp_eval_tab_col(self.lp, k + 1, c_indices, c_values) 

sig_off() 

except GLPKError: 

raise MIPSolverException('GLPK: basis factorization does not exist; or variable must be non-basic') 

  

indices = [c_indices[j+1] - 1 for j in range(i)] 

values = [c_values[j+1] for j in range(i)] 

return (indices, values) 

  

def __dealloc__(self): 

""" 

Destructor 

""" 

glp_delete_prob(self.lp) 

sig_free(self.iocp) 

sig_free(self.smcp) 

  

cdef void glp_callback(glp_tree* tree, void* info): 

r""" 

A callback routine called by glp_intopt 

  

This function fills the ``search_tree_data`` structure of a ``GLPKBackend`` 

object. It is fed to ``glp_intopt`` (which calls it while the search tree is 

being built) through ``iocp.cb_func``. 

  

INPUT: 

  

- ``tree`` -- a pointer toward ``glp_tree``, which is GLPK's search tree 

  

- ``info`` -- a ``void *`` to let the function know *where* it should store 

the data we need. The value of ``info`` is equal to the one stored in 

iocp.cb_info. 

  

""" 

cdef search_tree_data_t * data = <search_tree_data_t *> info 

data.mip_gap = glp_ios_mip_gap(tree) 

  

cdef int node_id = glp_ios_best_node(tree) 

data.best_bound = glp_ios_node_bound(tree, node_id) 

  

# parameter names 

  

cdef enum solver_parameter_names: 

timelimit_seconds, timelimit_simplex, timelimit_intopt, simplex_or_intopt, msg_lev_simplex, 

msg_lev_intopt, br_tech, bt_tech, pp_tech, fp_heur, gmi_cuts, 

mir_cuts, cov_cuts, clq_cuts, tol_int, tol_obj, mip_gap, 

tm_lim_intopt, out_frq_intopt, out_dly_intopt, presolve_intopt, 

binarize, meth, pricing, r_test, tol_bnd, tol_dj, tol_piv, obj_ll, 

obj_ul, it_lim, tm_lim_simplex, out_frq_simplex, out_dly_simplex, 

presolve_simplex 

  

glp_tm_lim_simplex = timelimit_simplex 

glp_tm_lim_intopt = timelimit_intopt 

glp_simplex_or_intopt = simplex_or_intopt 

glp_msg_lev_intopt = glp_verbosity_intopt = msg_lev_intopt 

glp_msg_lev_simplex = glp_verbosity_simplex = msg_lev_simplex 

glp_br_tech = glp_branching = br_tech 

glp_bt_tech = glp_backtracking = bt_tech 

glp_pp_tech = glp_preprocessing = pp_tech 

glp_fp_heur = glp_feasibility_pump = fp_heur 

glp_gmi_cuts = glp_gomory_cuts = gmi_cuts 

glp_mir_cuts = glp_mixed_int_rounding_cuts = mir_cuts 

glp_cov_cuts = glp_mixed_cover_cuts = cov_cuts 

glp_clq_cuts = glp_clique_cuts = clq_cuts 

glp_tol_int = glp_absolute_tolerance = tol_int 

glp_tol_obj = glp_relative_tolerance = tol_obj 

glp_mip_gap = glp_mip_gap_tolerance = mip_gap 

glp_out_frq_intopt = glp_output_frequency_intopt = out_frq_intopt 

glp_out_dly_intopt = glp_output_delay_intopt = out_dly_intopt 

glp_presolve_intopt = presolve_intopt 

glp_binarize = binarize 

glp_meth = glp_primal_v_dual = meth 

glp_pricing = pricing 

glp_r_test = glp_ratio_test = r_test 

glp_tol_bnd = glp_tolerance_primal = tol_bnd 

glp_tol_dj = glp_tolerance_dual = tol_dj 

glp_tol_piv = glp_tolerance_pivot = tol_piv 

glp_obj_ll = glp_obj_lower_limit = obj_ll 

glp_obj_ul = glp_obj_upper_limit = obj_ul 

glp_it_lim = glp_iteration_limit = it_lim 

glp_out_frq_simplex = glp_output_frequency_intopt = out_frq_simplex 

glp_out_dly_simplex = glp_output_delay_simplex = out_dly_simplex 

glp_presolve_simplex = presolve_simplex 

  

solver_parameter_names_dict = { 

'timelimit': timelimit_seconds, 

'timelimit_intopt': timelimit_intopt, 

'tm_lim_intopt': timelimit_intopt, 

'timelimit_simplex': timelimit_simplex, 

'tm_lim_simplex': timelimit_simplex, 

'simplex_or_intopt': simplex_or_intopt, 

'msg_lev_simplex': msg_lev_simplex, 'verbosity_simplex': msg_lev_simplex, 

'msg_lev_intopt': msg_lev_intopt, 'verbosity_intopt': msg_lev_intopt, 

'br_tech': br_tech, 'branching': br_tech, 

'bt_tech': bt_tech, 'backtracking': bt_tech, 

'pp_tech': pp_tech, 'preprocessing': pp_tech, 

'fp_heur': fp_heur, 'feasibility_pump': fp_heur, 

'gmi_cuts': gmi_cuts, 'gomory_cuts': gmi_cuts, 

'mir_cuts': mir_cuts, 'mixed_int_rounding_cuts': mir_cuts, 

'cov_cuts': cov_cuts, 'mixed_cover_cuts': cov_cuts, 

'clq_cuts': clq_cuts, 'clique_cuts': clq_cuts, 

'tol_int': tol_int, 'absolute_tolerance': tol_int, 

'tol_obj': tol_obj, 'relative_tolerance': tol_obj, 

'mip_gap': mip_gap, 'mip_gap_tolerance': mip_gap, 

'out_frq_intopt': out_frq_intopt, 'output_frequency_intopt': out_frq_intopt, 

'out_dly_intopt': out_dly_intopt, 'output_delay_intopt': out_dly_intopt, 

'presolve_intopt': presolve_intopt, 'binarize': binarize, 

'meth': meth, 'primal_v_dual': meth, 

'pricing': pricing, 

'r_test': r_test, 'ratio_test': r_test, 

'tol_bnd': tol_bnd, 'tolerance_primal': tol_bnd, 

'tol_dj': tol_dj, 'tolerance_dual': tol_dj, 

'tol_piv': tol_piv, 'tolerance_pivot': tol_piv, 

'obj_ll': obj_ll, 'obj_lower_limit': obj_ll, 

'obj_ul': obj_ul, 'obj_upper_limit': obj_ul, 

'it_lim': it_lim, 'iteration_limit': it_lim, 

'out_frq_simplex': out_frq_simplex, 'output_frequency_intopt': out_frq_simplex, 

'out_dly_simplex': out_dly_simplex, 'output_delay_simplex': out_dly_simplex, 

'presolve_simplex': presolve_simplex 

} 

  

# parameter values 

  

glp_msg_off = GLP_MSG_OFF 

glp_msg_on = GLP_MSG_ON 

glp_msg_err = GLP_MSG_ERR 

glp_msg_all = GLP_MSG_ALL 

glp_msg_dbg = GLP_MSG_DBG 

  

glp_primal = GLP_PRIMAL 

glp_dual = GLP_DUAL 

glp_dualp = GLP_DUALP 

  

glp_pt_std = GLP_PT_STD 

glp_pt_pse = GLP_PT_PSE 

  

glp_rt_std = GLP_RT_STD 

glp_rt_har = GLP_RT_HAR 

  

dbl_max = DBL_MAX 

int_max = INT_MAX 

  

glp_on = GLP_ON 

glp_off = GLP_OFF 

  

glp_br_ffv = GLP_BR_FFV 

glp_br_lfv = GLP_BR_LFV 

glp_br_mfv = GLP_BR_MFV 

glp_br_dth = GLP_BR_DTH 

glp_br_pch = GLP_BR_PCH 

  

glp_bt_dfs = GLP_BT_DFS 

glp_bt_bfs = GLP_BT_BFS 

glp_bt_blb = GLP_BT_BLB 

glp_bt_bph = GLP_BT_BPH 

  

glp_pp_none = GLP_PP_NONE 

glp_pp_root = GLP_PP_ROOT 

glp_pp_all = GLP_PP_ALL 

  

glp_max = GLP_MAX 

glp_min = GLP_MIN 

glp_up = GLP_UP 

glp_fr = GLP_FR 

glp_db = GLP_DB 

glp_fx = GLP_FX 

glp_lo = GLP_LO 

glp_cv = GLP_CV 

glp_iv = GLP_IV 

glp_bv = GLP_BV 

glp_mps_deck = GLP_MPS_DECK 

glp_mps_file = GLP_MPS_FILE 

  

glp_undef = GLP_UNDEF 

glp_opt = GLP_OPT 

glp_feas = GLP_FEAS 

glp_nofeas = GLP_NOFEAS 

  

glp_bs = GLP_BS 

glp_nl = GLP_NL 

glp_nu = GLP_NU 

glp_nf = GLP_NF 

  

cdef enum more_parameter_values: 

simplex_only, simplex_then_intopt, intopt_only, exact_simplex_only 

  

glp_simplex_only = simplex_only 

glp_simplex_then_intopt = simplex_then_intopt 

glp_intopt_only = intopt_only 

glp_exact_simplex_only = exact_simplex_only 

  

# dictionaries for those who prefer to use strings 

  

solver_parameter_values = { 

  

'simplex_only': simplex_only, 

'simplex_then_intopt': simplex_then_intopt, 

'intopt_only': intopt_only, 

'exact_simplex_only': exact_simplex_only, 

  

'GLP_MSG_OFF' : GLP_MSG_OFF, 

'GLP_MSG_ON' : GLP_MSG_ON, 

'GLP_MSG_ERR' : GLP_MSG_ERR, 

'GLP_MSG_ALL' : GLP_MSG_ALL, 

'GLP_MSG_DBG' : GLP_MSG_DBG, 

  

'GLP_PRIMAL' : GLP_PRIMAL, 

'GLP_DUAL' : GLP_DUAL, 

'GLP_DUALP' : GLP_DUALP, 

  

'GLP_PT_STD' : GLP_PT_STD, 

'GLP_PT_PSE' : GLP_PT_PSE, 

  

'GLP_RT_STD' : GLP_RT_STD, 

'GLP_RT_HAR' : GLP_RT_HAR, 

  

'DBL_MAX' : DBL_MAX, 

'INT_MAX' : INT_MAX, 

  

'GLP_ON' : GLP_ON, 

'GLP_OFF' : GLP_OFF, 

  

'GLP_BR_FFV' : GLP_BR_FFV, 

'GLP_BR_LFV' : GLP_BR_LFV, 

'GLP_BR_MFV' : GLP_BR_MFV, 

'GLP_BR_DTH' : GLP_BR_DTH, 

'GLP_BR_PCH' : GLP_BR_PCH, 

  

'GLP_BT_DFS' : GLP_BT_DFS, 

'GLP_BT_BFS' : GLP_BT_BFS, 

'GLP_BT_BLB' : GLP_BT_BLB, 

'GLP_BT_BPH' : GLP_BT_BPH, 

  

'GLP_PP_NONE' : GLP_PP_NONE, 

'GLP_PP_ROOT' : GLP_PP_ROOT, 

'GLP_PP_ALL' : GLP_PP_ALL, 

  

'GLP_MAX' : GLP_MAX, 

'GLP_MIN' : GLP_MIN, 

'GLP_UP' : GLP_UP, 

'GLP_FR' : GLP_FR, 

'GLP_DB' : GLP_DB, 

'GLP_FX' : GLP_FX, 

'GLP_LO' : GLP_LO, 

'GLP_CV' : GLP_CV, 

'GLP_IV' : GLP_IV, 

'GLP_BV' : GLP_BV, 

'GLP_MPS_DECK' : GLP_MPS_DECK, 

'GLP_MPS_FILE' : GLP_MPS_FILE, 

  

'GLP_UNDEF' : GLP_UNDEF, 

'GLP_OPT' : GLP_OPT, 

'GLP_FEAS' : GLP_FEAS, 

'GLP_NOFEAS' : GLP_NOFEAS 

  

} 

  

cdef dict solve_status_msg = { 

GLP_EBADB: "The initial basis specified in the problem object is invalid", 

GLP_ESING: "The basis matrix corresponding to the initial basis is singular within the working precision", 

GLP_ECOND: "The basis matrix corresponding to the initial basis is ill-conditioned, i.e. its condition number is too large", 

GLP_EBOUND: "Some variables (auxiliary or structural) have incorrect bounds", 

GLP_EFAIL: "Solver failure", 

GLP_EOBJLL: "The objective lower limit has been reached", 

GLP_EOBJUL: "The objective upper limit has been reached", 

GLP_EITLIM: "The iteration limit has been exceeded", 

GLP_ETMLIM: "The time limit has been exceeded", 

GLP_ENOPFS: "The LP (relaxation) problem has no primal feasible solution", 

GLP_ENODFS: "The LP (relaxation) problem has no dual feasible solution", 

GLP_EROOT: "Optimal basis for initial LP relaxation is not provided", 

GLP_ESTOP: "The search was prematurely terminated by application", 

GLP_EMIPGAP: "The relative mip gap tolerance has been reached", 

} 

  

cdef dict solution_status_msg = { 

GLP_UNDEF: "Solution is undefined", 

GLP_FEAS: "Feasible solution found, while optimality has not been proven", 

GLP_INFEAS: "Solution is infeasible", 

GLP_NOFEAS: "Problem has no feasible solution", 

GLP_OPT: "Solution is optimal", 

GLP_UNBND: "Problem has unbounded solution", 

}