CGAL 6.1 - Linear and Quadratic Programming Solver
Loading...
Searching...
No Matches
QP_solver/first_lp.cpp
// example: construct a linear program from data
// the LP below is the first linear program example in the user manual
#include <iostream>
#include <cassert>
#include <CGAL/QP_models.h>
#include <CGAL/QP_functions.h>
// choose exact integral type
#ifdef CGAL_USE_GMP
#include <CGAL/Gmpz.h>
typedef CGAL::Gmpz ET;
#else
#include <CGAL/MP_Float.h>
typedef CGAL::MP_Float ET;
#endif
// program and solution types
int main() {
// by default, we have a nonnegative LP with Ax <= b
Program lp (CGAL::SMALLER, true, 0, false, 0);
// now set the non-default entries
const int X = 0;
const int Y = 1;
lp.set_a(X, 0, 1); lp.set_a(Y, 0, 1); lp.set_b(0, 7); // x + y <= 7
lp.set_a(X, 1, -1); lp.set_a(Y, 1, 2); lp.set_b(1, 4); // -x + 2y <= 4
lp.set_u(Y, true, 4); // y <= 4
lp.set_c(Y, -32); // -32y
lp.set_c0(64); // +64
// solve the program, using ET as the exact type
Solution s = CGAL::solve_linear_program(lp, ET());
assert (s.solves_linear_program(lp));
// output solution
std::cout << s;
return 0;
}
An object of class Quadratic_program_solution represents the solution of a linear or convex quadratic...
Definition: QP_solution.h:65
An object of class Quadratic_program describes a convex quadratic program of the form.
Definition: QP_models.h:788
Quadratic_program_solution< ET > solve_linear_program(const LinearProgram &lp, const ET &, const Quadratic_program_options &options=Quadratic_program_options())
This function solves a linear program, using some exact Integral Domain ET for its computations.