\( \newcommand{\E}{\mathrm{E}} \) \( \newcommand{\A}{\mathrm{A}} \) \( \newcommand{\R}{\mathrm{R}} \) \( \newcommand{\N}{\mathrm{N}} \) \( \newcommand{\Q}{\mathrm{Q}} \) \( \newcommand{\Z}{\mathrm{Z}} \) \( \def\ccSum #1#2#3{ \sum_{#1}^{#2}{#3} } \def\ccProd #1#2#3{ \sum_{#1}^{#2}{#3} }\)
CGAL 5.0 - Linear and Quadratic Programming Solver
QP_solver/solve_convex_hull_containment_lp.h
// example: function to check whether a point is in the convex
// hull of other points
#include <CGAL/QP_models.h>
#include <CGAL/QP_functions.h>
// unary function to get homogeneous begin-iterator of point
template <class Point_d>
struct Homogeneous_begin {
typedef typename Point_d::Homogeneous_const_iterator result_type;
result_type operator() (const Point_d& p) const {
return p.homogeneous_begin();
}
};
// function to solve the LP that tests whether a point is in the
// convex hull of other points; the type ET is an exact type used
// for the internal computations
template <class Point_d, class RandomAccessIterator, class ET>
solve_convex_hull_containment_lp (const Point_d& p,
RandomAccessIterator end, const ET& dummy)
{
// Constraint matrix type: A[j][i] is the i-th homogeneous coordinate of p_j
typedef boost::transform_iterator
<Homogeneous_begin<Point_d>, RandomAccessIterator> A_it;
// Right-hand side type: b[i] is the i-th homogeneous coordinate of p
typedef typename Point_d::Homogeneous_const_iterator B_it;
// Relation type ("=")
// input number type
// Linear objective function type (c=0: we only test feasibility)
// the nonnegative linear program type
typedef
Program;
// ok, we are prepared now: construct program and solve it
Program lp (static_cast<int>(end-begin), // number of variables
p.dimension()+1, // number of constraints
A_it (begin), B_it (p.homogeneous_begin()),
R_it (CGAL::EQUAL), C_it (0));
}