CGAL 5.3 - Linear and Quadratic Programming Solver
|
#include <CGAL/QP_models.h>
An object of class Nonnegative_linear_program_from_iterators
describes a linear program of the form.
\begin{eqnarray*} \mbox{(QP)}& \mbox{minimize} &\qpc^{T}\qpx+c_0 \\ &\mbox{subject to} & A\qpx\qprel \qpb, \\ & & \qpx \geq 0 \end{eqnarray*}
in \( n\) real variables \( \qpx=(x_0,\ldots,x_{n-1})\).
Here,
\( \qprel\) is an \( m\)-dimensional vector of relations from \( \{\leq, =, \geq\}\),
\( c_0\) is a constant.
This class is simply a wrapper for existing iterators, and it does not copy the program data.
It frequently happens that all values in one of the vectors from above are the same, for example if the system \( Ax\qprel b\) is actually a system of equations \( Ax=b\). To get an iterator over such a vector, it is not necessary to store multiple copies of the value in some container; an instance of the class Const_oneset_iterator<T>
, constructed from the value in question, does the job more efficiently.
Example
QP_solver/first_nonnegative_lp_from_iterators.cpp
QP_solver/solve_convex_hull_containment_lp.h
QP_solver/convex_hull_containment.cpp
NonnegativeLinearProgram
Quadratic_program<NT>
Quadratic_program_from_mps<NT>
Creation | |
Nonnegative_linear_program_from_iterators (int n, int m, const A_it &a, const B_it &b, const R_it &r, const C_it &c, const std::iterator_traits< C_it >value_type &c0=0) | |
constructs lp from given random-access iterators and the constant c0 . More... | |
CGAL::Nonnegative_linear_program_from_iterators< A_it, B_it, R_it, C_it >::Nonnegative_linear_program_from_iterators | ( | int | n, |
int | m, | ||
const A_it & | a, | ||
const B_it & | b, | ||
const R_it & | r, | ||
const C_it & | c, | ||
const std::iterator_traits< C_it >value_type & | c0 = 0 |
||
) |
constructs lp
from given random-access iterators and the constant c0
.
The passed iterators are merely stored, no copying of the program data takes place. How these iterators are supposed to encode the nonnegative linear program is described in NonnegativeLinearProgram
.