CGAL 6.1 - Linear and Quadratic Programming Solver
Loading...
Searching...
No Matches
QP_solver/first_qp_from_iterators.cpp
// example: construct a quadratic program from given iterators
// the QP below is the first quadratic 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
bool*, // for fl
int*, // for l
bool*, // for fu
int*, // for u
int**, // for D
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 "<="
bool fl[] = {true, true}; // both x, y are lower-bounded
int l[] = {0, 0};
bool fu[] = {false, true}; // only y is upper-bounded
int u[] = {0, 4}; // x's u-entry is ignored
int D1[] = {2}; // 2D_{1,1}
int D2[] = {0, 8}; // 2D_{2,1}, 2D_{2,2}
int* D[] = {D1, D2}; // D-entries on/below diagonal
int c[] = {0, -32};
int c0 = 64; // constant term
// now construct the quadratic program; the first two parameters are
// the number of variables and the number of constraints (rows of A)
Program qp (2, 2, A, b, r, fl, l, fu, u, D, c, c0);
// solve the program, using ET as the exact type
Solution s = CGAL::solve_quadratic_program(qp, ET());
// output solution
std::cout << s;
return 0;
}
An object of class Quadratic_program_from_iterators describes a convex quadratic program of the form.
Definition: QP_models.h:486
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_quadratic_program(const QuadraticProgram &qp, const ET &, const Quadratic_program_options &options=Quadratic_program_options())
This function solves a quadratic program, using some exact Integral Domain ET for its computations.