CGAL 6.0.1 - Linear and Quadratic Programming Solver
Loading...
Searching...
No Matches
QP_solver/first_nonnegative_lp_from_iterators.cpp
// example: construct a nonnegative linear program from given iterators
// the LP below is the first nonnegative linear program example
// in the user manual
#include <iostream>
#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**, // for A
int*, // for b
int*> // for c
Program;
int main() {
int Ax[] = {1, -1}; // column for x
int Ay[] = {1, 2}; // column for y
int* A[] = {Ax, Ay}; // A comes columnwise
int b[] = {7, 4}; // right-hand side
r( CGAL::SMALLER); // constraints are "<="
int c[] = {0, -32};
// now construct the linear program; the first two parameters are
// the number of variables and the number of constraints (rows of A)
Program lp (2, 2, A, b, r, c); // constant term defaults to 0
// solve the program, using ET as the exact type
// output solution
std:: cout << s;
return 0;
}
An object of class Nonnegative_linear_program_from_iterators describes a linear program of the form.
Definition: QP_models.h:314
An object of class Quadratic_program_solution represents the solution of a linear or convex quadratic...
Definition: QP_solution.h:65
Quadratic_program_solution< ET > solve_nonnegative_linear_program(const NonnegativeLinearProgram &lp, const ET &, const Quadratic_program_options &options=Quadratic_program_options())
This function solves a nonnegative linear program, using some exact Integral Domain ET for its comput...