CGAL 4.8.1 - Linear and Quadratic Programming Solver
|
This module provides makers to construct a program, as well as functions to solve and print programs.
In case you want to construct a program from complicated iterators (whose types you don't know, or simply don't want to bother with), we provide four makers.
CGAL::make_quadratic_program_from_iterators()
CGAL::make_linear_program_from_iterators()
CGAL::make_nonnegative_quadratic_program_from_iterators()
CGAL::make_nonnegative_linear_program_from_iterators()
There are four functions to solve a program, one for each program concept.
CGAL::solve_quadratic_program()
CGAL::solve_linear_program()
CGAL::solve_nonnegative_quadratic_program()
CGAL::solve_nonnegative_linear_program()
The solution process can customized by passing an object of the class
Programs can be written to an output stream in MPSFormat, using one of the following four functions.
CGAL::print_quadratic_program()
CGAL::print_linear_program()
CGAL::print_nonnegative_quadratic_program()
CGAL::print_nonnegative_linear_program()
Classes | |
class | CGAL::Quadratic_program_options |
This is a class used for passing options to the linear and quadratic programming solvers. More... | |
Enumerations | |
enum | CGAL::Quadratic_program_pricing_strategy { CGAL::QP_CHOOSE_DEFAULT, CGAL::QP_PARTIAL_DANTZIG, CGAL::QP_DANTZIG, CGAL::QP_PARTIAL_FILTERED_DANTZIG, CGAL::QP_FILTERED_DANTZIG, CGAL::QP_BLAND } |
This is an enumeration type containing the values QP_CHOOSE_DEFAULT , QP_DANTZIG , QP_PARTIAL_DANTZIG , QP_FILTERED_DANTZIG , QP_PARTIAL_FILTERED_DANTZIG , andQP_BLAND . More... | |
Functions | |
template<LinearProgram > | |
void | CGAL::print_linear_program (std::ostream &out, const LinearProgram &lp, const std::string &problem_name=std::string("MY_MPS")) |
This function writes a linear program to an output stream (in MPSFormat ). More... | |
template<NonnegativeLinearProgram > | |
void | CGAL::print_nonnegative_linear_program (std::ostream &out, const NonnegativeLinearProgram &lp, const std::string &problem_name=std::string("MY_MPS")) |
This function writes a nonnegative linear program to an output stream (in MPSFormat ). More... | |
template<NonnegativeQuadraticProgram > | |
void | CGAL::print_nonnegative_quadratic_program (std::ostream &out, const NonnegativeQuadraticProgram &qp, const std::string &problem_name=std::string("MY_MPS")) |
This function writes a nonnegative quadratic program to an output stream (in MPSFormat ). More... | |
template<QuadraticProgram > | |
void | CGAL::print_quadratic_program (std::ostream &out, const QuadraticProgram &qp, const std::string &problem_name=std::string("MY_MPS")) |
This function writes a quadratic program to an output stream (in MPSFormat ). More... | |
template<LinearProgram , ET > | |
Quadratic_program_solution< ET > | CGAL::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. More... | |
template<NonnegativeLinearProgram , ET > | |
Quadratic_program_solution< ET > | CGAL::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 computations. More... | |
template<NonnegativeQuadraticProgram , ET > | |
Quadratic_program_solution< ET > | CGAL::solve_nonnegative_quadratic_program (const NonnegativeQuadraticProgram &qp, const ET &, const Quadratic_program_options &options=Quadratic_program_options()) |
This function solves a nonnegative quadratic program, using some exact Integral Domain ET for its computations. More... | |
template<QuadraticProgram , ET > | |
Quadratic_program_solution< ET > | CGAL::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. More... | |
template<typename A_it , typename B_it , typename R_it , typename FL_it , typename L_it , typename FU_it , typename U_it , typename C_it > | |
Linear_program_from_iterators < A_it, B_it, R_it, FL_it, L_it, FU_it, U_it, C_it > | CGAL::make_linear_program_from_iterators (int n, int m, const A_it &a, const B_it &b, const R_it &r, const FL_it &fl, const L_it &l, const FU_it &fu, const U_it &u, const C_it &c, std::iterator_traits< C_it >::value_type c0=std::iterator_traits< C_it >::value_type(0)) |
This template function creates an instance of Linear_program_from_iterators<A_it, B_it, R_it, FL_it, L_it, FU_it, U_it, C_it> from given iterators. More... | |
template<A_it , B_it , R_it , C_it > | |
Nonnegative_linear_program_from_iterators < A_it, B_it, R_it, C_it > | CGAL::make_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, std::iterator_traits< C_it >::value_type c0=std::iterator_traits< C_it >::value_type(0)) |
This template function creates an instance of Nonnegative_linear_program_from_iterators<A_it, B_it, R_it, C_it> from given iterators. More... | |
template<A_it , B_it , R_it , D_it , C_it > | |
Nonnegative_quadratic_program_from_iterators < A_it, B_it, R_it, D_it, C_it > | CGAL::make_nonnegative_quadratic_program_from_iterators (int n, int m, const A_it &a, const B_it &b, const R_it &r, const D_it &d, const C_it &c, std::iterator_traits< C_it >::value_type c0=std::iterator_traits< C_it >::value_type(0)) |
This template function creates an instance of Nonnegative_quadratic_program_from_iterators<A_it, B_it, R_it, D_it, C_it> from given iterators. More... | |
template<typename A_it , typename B_it , typename R_it , typename FL_it , typename L_it , typename FU_it , typename U_it , typename D_it , typename C_it > | |
Quadratic_program_from_iterators < A_it, B_it, R_it, FL_it, L_it, FU_it, U_it, D_it, C_it > | CGAL::make_quadratic_program_from_iterators (int n, int m, const A_it &a, const B_it &b, const R_it &r, const FL_it &fl, const L_it &l, const FU_it &fu, const U_it &u, const D_it &d, const C_it &c, std::iterator_traits< C_it >::value_type c0=std::iterator_traits< C_it >::value_type(0)) |
This template function creates an instance of Quadratic_program_from_iterators<A_it, B_it, R_it, FL_it, L_it, FU_it, U_it, D_it, C_it> from given iterators. More... | |
This is an enumeration type containing the values QP_CHOOSE_DEFAULT
, QP_DANTZIG
, QP_PARTIAL_DANTZIG
, QP_FILTERED_DANTZIG
, QP_PARTIAL_FILTERED_DANTZIG
, andQP_BLAND
.
It indicates the pricing strategy to be used in solving a linear or quadratic program. This strategy determines how the solver gets from one intermediate solution to the next during any of its iterations.
Here we briefly describe when to choose which strategy.
Quadratic_program_options
#include <CGAL/QP_options.h>
Enumerator | |
---|---|
QP_CHOOSE_DEFAULT |
This is the default value of the pricing strategy in There are only few reasons to deviate from this default, but you are free to experiment, of course. |
QP_PARTIAL_DANTZIG |
If the input type is not |
QP_DANTZIG |
If the input type is not |
QP_PARTIAL_FILTERED_DANTZIG |
If the input type is If the input type is not Note: filtered strategies may in rare cases fail due to double exponent overflows, see Section Exponent Overflow in Double Using Floating-Point Filters. In this case, the slower fallback option is the non-filtered variant |
QP_FILTERED_DANTZIG |
If the input type is If the input type is not
|
QP_BLAND |
This is hardly ever the most efficient choice, but it is guaranteed to avoid internal cycling of the solution algorithm, see Section The Solver Internally Cycles. |
Linear_program_from_iterators<A_it, B_it, R_it, FL_it, L_it, FU_it, U_it, C_it> CGAL::make_linear_program_from_iterators | ( | int | n, |
int | m, | ||
const A_it & | a, | ||
const B_it & | b, | ||
const R_it & | r, | ||
const FL_it & | fl, | ||
const L_it & | l, | ||
const FU_it & | fu, | ||
const U_it & | u, | ||
const C_it & | c, | ||
std::iterator_traits< C_it >::value_type | c0 = std::iterator_traits< C_it >::value_type(0) |
||
) |
This template function creates an instance of Linear_program_from_iterators<A_it, B_it, R_it, FL_it, L_it, FU_it, U_it, C_it>
from given iterators.
This function can be useful if the types of these iterators are too complicated (or of too little interest for you) to write them down explicitly.
Linear_program_from_iterators<A_it, B_it, R_it, FL_it, L_it, FU_it, U_it, C_it>
, constructed from the given iterators.The following example demonstrates the typical usage of makers with the simpler function make_nonnegative_linear_program_from_iterators()
.
QP_solver/solve_convex_hull_containment_lp2.h
QP_solver/convex_hull_containment2.cpp
Linear_program_from_iterators<A_it, B_it, R_it, FL_it, L_it, FU_it, U_it, C_it>
#include <CGAL/QP_models.h>
Nonnegative_linear_program_from_iterators<A_it, B_it, R_it, C_it> CGAL::make_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, | ||
std::iterator_traits< C_it >::value_type | c0 = std::iterator_traits< C_it >::value_type(0) |
||
) |
This template function creates an instance of Nonnegative_linear_program_from_iterators<A_it, B_it, R_it, C_it>
from given iterators.
This function can be useful if the types of these iterators are too complicated (or of too little interest for you) to write them down explicitly.
Nonnegative_linear_program_from_iterators<A_it, B_it, R_it, C_it>
, constructed from the given iterators.QP_solver/solve_convex_hull_containment_lp2.h
QP_solver/convex_hull_containment2.cpp
Nonnegative_linear_program_from_iterators<A_it, B_it, R_it, C_it>
#include <CGAL/QP_models.h>
Nonnegative_quadratic_program_from_iterators<A_it, B_it, R_it, D_it, C_it> CGAL::make_nonnegative_quadratic_program_from_iterators | ( | int | n, |
int | m, | ||
const A_it & | a, | ||
const B_it & | b, | ||
const R_it & | r, | ||
const D_it & | d, | ||
const C_it & | c, | ||
std::iterator_traits< C_it >::value_type | c0 = std::iterator_traits< C_it >::value_type(0) |
||
) |
This template function creates an instance of Nonnegative_quadratic_program_from_iterators<A_it, B_it, R_it, D_it, C_it>
from given iterators.
This function can be useful if the types of these iterators are too complicated (or of too little interest for you) to write them down explicitly.
Nonnegative_quadratic_program_from_iterators<A_it, B_it, R_it, D_it,C_it>
, constructed from the given iterators.The following example demonstrates the typical usage of makers with the simpler function make_nonnegative_linear_program_from_iterators()
.
QP_solver/solve_convex_hull_containment_lp2.h
QP_solver/convex_hull_containment2.cpp
Nonnegative_quadratic_program_from_iterators<A_it, B_it, R_it, D_it, C_it>
#include <CGAL/QP_models.h>
Quadratic_program_from_iterators<A_it, B_it, R_it, FL_it, L_it, FU_it, U_it, D_it, C_it> CGAL::make_quadratic_program_from_iterators | ( | int | n, |
int | m, | ||
const A_it & | a, | ||
const B_it & | b, | ||
const R_it & | r, | ||
const FL_it & | fl, | ||
const L_it & | l, | ||
const FU_it & | fu, | ||
const U_it & | u, | ||
const D_it & | d, | ||
const C_it & | c, | ||
std::iterator_traits< C_it >::value_type | c0 = std::iterator_traits< C_it >::value_type(0) |
||
) |
This template function creates an instance of Quadratic_program_from_iterators<A_it, B_it, R_it, FL_it, L_it, FU_it, U_it, D_it, C_it>
from given iterators.
This function can be useful if the types of these iterators are too complicated (or of too little interest for you) to write them down explicitly.
Quadratic_program_from_iterators<A_it, B_it, R_it, FL_it, L_it, FU_it, U_it, D_it, C_it>
, constructed from the given iterators.The following example demonstrates the typical usage of makers with the simpler function make_nonnegative_linear_program_from_iterators()
.
QP_solver/solve_convex_hull_containment_lp2.h
QP_solver/convex_hull_containment2.cpp
Quadratic_program_from_iterators<A_it, B_it, R_it, FL_it, L_it, FU_it, U_it, D_it, C_it>
#include <CGAL/QP_models.h>
void CGAL::print_linear_program | ( | std::ostream & | out, |
const LinearProgram & | lp, | ||
const std::string & | problem_name = std::string("MY_MPS") |
||
) |
This function writes a linear program to an output stream (in MPSFormat
).
The time complexity is \( \Theta (mn)\), even if the program is very sparse.
It writes the linear program lp
to out
in MPSFormat
. The name of the program will be the one provided by problem_name
.
Output operators are defined for all entry types of lp
.
LinearProgram
#include <CGAL/QP_functions.h>
void CGAL::print_nonnegative_linear_program | ( | std::ostream & | out, |
const NonnegativeLinearProgram & | lp, | ||
const std::string & | problem_name = std::string("MY_MPS") |
||
) |
This function writes a nonnegative linear program to an output stream (in MPSFormat
).
The time complexity is \( \Theta (mn)\), even if the program is very sparse.
Writes the nonnegative linear program lp
to out
in MPSFormat
. The name of the program will be the one provided by problem_name
.
Output operators are defined for all entry types of lp
.
QP_solver/print_first_nonnegative_lp.cpp
NonnegativeLinearProgram
#include <CGAL/QP_functions.h>
void CGAL::print_nonnegative_quadratic_program | ( | std::ostream & | out, |
const NonnegativeQuadraticProgram & | qp, | ||
const std::string & | problem_name = std::string("MY_MPS") |
||
) |
This function writes a nonnegative quadratic program to an output stream (in MPSFormat
).
The time complexity is \( \Theta (n^2 + mn)\), even if the program is very sparse.
Writes the nonnegative quadratic program qp
to out
in MPSFormat
. The name of the program will be the one provided by problem_name
.
Output operators are defined for all entry types of qp
.
QP_solver/print_first_nonnegative_qp.cpp
NonnegativeQuadraticProgram
#include <CGAL/QP_functions.h>
void CGAL::print_quadratic_program | ( | std::ostream & | out, |
const QuadraticProgram & | qp, | ||
const std::string & | problem_name = std::string("MY_MPS") |
||
) |
This function writes a quadratic program to an output stream (in MPSFormat
).
The time complexity is \( \Theta (n^2 + mn)\), even if the program is very sparse.
Writes the quadratic program qp
to out
in MPSFormat
. The name of the program will be the one provided by problem_name
.
Output operators are defined for all entry types of qp
.
QuadraticProgram
#include <CGAL/QP_functions.h>
Quadratic_program_solution<ET> CGAL::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.
Various options may be provided, see Quadratic_program_options
.
ET
is a model of the concepts IntegralDomain
and RealEmbeddable
; it must be an exact type, and all entries of qp
are convertible to ET
.
Here are some recommended combinations of input type (the type of the qp
entries) and ET
.
input type | ET |
---|---|
double | MP_Float , Gmpzf , or Gmpq |
int | MP_Float , or Gmpz |
any exact type NT | NT |
CGAL_QP_NO_ASSERTIONS
or NDEBUG
.lp
, solved with exact number type ET
.Quadratic_program<NT>
Quadratic_program_from_mps<NT>
Linear_program_from_iterators<A_it, B_it, R_it, FL_it, L_it, FU_it, U_it, C_it>
#include <CGAL/QP_functions.h>
Quadratic_program_solution<ET> CGAL::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 computations.
Various options may be provided, see Quadratic_program_options
.
ET
is a model of the concepts IntegralDomain
and RealEmbeddable
; it must be an exact type, and all entries of qp
are convertible to ET
.
Here are some recommended combinations of input type (the type of the qp
entries) and ET
.
input type | ET |
---|---|
double | MP_Float , Gmpzf , or Gmpq |
int | MP_Float , or Gmpz |
any exact type NT | NT |
CGAL_QP_NO_ASSERTIONS
or NDEBUG
.lp
, solved with exact number type ET
.QP_solver/first_nonnegative_lp.cpp
The models of NonnegativeLinearProgram\:
Quadratic_program<NT>
Quadratic_program_from_mps<NT>
Nonnegative_linear_program_from_iterators<A_it, B_it, R_it, C_it>
#include <CGAL/QP_functions.h>
Quadratic_program_solution<ET> CGAL::solve_nonnegative_quadratic_program | ( | const NonnegativeQuadraticProgram & | qp, |
const ET & | , | ||
const Quadratic_program_options & | options = Quadratic_program_options() |
||
) |
This function solves a nonnegative quadratic program, using some exact Integral Domain ET
for its computations.
Various options may be provided, see Quadratic_program_options
.
ET
is a model of the concepts IntegralDomain
and RealEmbeddable
; it must be an exact type, and all entries of qp
are convertible to ET
.
Here are some recommended combinations of input type (the type of the qp
entries) and ET
.
input type | ET |
---|---|
double | MP_Float , Gmpzf , or Gmpq |
int | MP_Float , or Gmpz |
any exact type NT | NT |
CGAL_QP_NO_ASSERTIONS
or NDEBUG
.qp
, solved with exact number type ET
.QP_solver/first_nonnegative_qp.cpp
The models of NonnegativeQuadraticProgram\:
Quadratic_program<NT>
Quadratic_program_from_mps<NT>
Nonnegative_quadratic_program_from_iterators<A_it, B_it, R_it, D_it, C_it>
#include <CGAL/QP_functions.h>
Quadratic_program_solution<ET> CGAL::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.
Various options may be provided, see Quadratic_program_options
.
ET
is a model of the concepts IntegralDomain
and RealEmbeddable
; it must be an exact type, and all entries of qp
are convertible to ET
.
Here are some recommended combinations of input type (the type of the qp
entries) and ET
.
input type | ET |
---|---|
double | MP_Float , Gmpzf , or Gmpq |
int | MP_Float , or Gmpz |
any exact type NT | NT |
CGAL_QP_NO_ASSERTIONS
or NDEBUG
.qp
, solved with exact number type ET
.The models of QuadraticProgram\:
Quadratic_program<NT>
Quadratic_program_from_mps<NT>
Quadratic_program_from_iterators<A_it, B_it, R_it, FL_it, L_it, FU_it, U_it, D_it, C_it>
#include <CGAL/QP_functions.h>