politopix  4.1.0
 All Classes Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
politopixAPI Class Reference

The main class to run the basic API. More...

#include <politopixAPI.h>

Collaboration diagram for politopixAPI:
Collaboration graph

Static Public Member Functions

static polito_EXPORT int savePolytope (const string &pathA, boost::shared_ptr< Polytope_Rn > &A) throw (ios_base::failure)
 Save a polytope in the corresponding file name. More...
 
static polito_EXPORT int savePolyhedralCone (const string &pathA, boost::shared_ptr< PolyhedralCone_Rn > &A) throw (ios_base::failure)
 Save a polyhedral cone in the corresponding file name. More...
 
static polito_EXPORT int loadPolytope (const string &pathA, boost::shared_ptr< Polytope_Rn > &A) throw (ios_base::failure)
 Load a polytope from the corresponding file name. More...
 
static polito_EXPORT int loadPolyhedralCone (const string &pathA, boost::shared_ptr< PolyhedralCone_Rn > &A) throw (ios_base::failure)
 Load a polyhedral cone from the corresponding file name. More...
 
static polito_EXPORT int addHalfspace (boost::shared_ptr< Polytope_Rn > &A, const boost::shared_ptr< HalfSpace_Rn > &HS)
 Add a new half-space into the polytope data structure. More...
 
static polito_EXPORT int addGenerator (boost::shared_ptr< Polytope_Rn > &A, const boost::shared_ptr< Generator_Rn > &GN)
 Add a new generator into the polytope data structure. More...
 
static polito_EXPORT int computeDoubleDescription (boost::shared_ptr< Polytope_Rn > &A, double bb_size=1000.) throw (invalid_argument, out_of_range, ios_base::failure, logic_error)
 Compute the HV-description for a given H-polytope or V-polytope with the Double Description algorithm. More...
 
static polito_EXPORT int computeDoubleDescriptionWithoutCheck (boost::shared_ptr< Polytope_Rn > &A, double bb_size=1000.) throw (invalid_argument, out_of_range, ios_base::failure, logic_error)
 Compute the HV-description for a given H-polytope or V-polytope with the Double Description algorithm. More...
 
static polito_EXPORT int computeIntersection (const boost::shared_ptr< Polytope_Rn > &A, const boost::shared_ptr< Polytope_Rn > &B, boost::shared_ptr< Polytope_Rn > &C) throw (invalid_argument, out_of_range, ios_base::failure, logic_error)
 Compute the intersection between two HV-polytopes with the Double Description algorithm. More...
 
static polito_EXPORT int computeIntersection (boost::shared_ptr< Polytope_Rn > &A, const boost::shared_ptr< Polytope_Rn > &B) throw (invalid_argument, out_of_range, ios_base::failure, logic_error)
 Compute the intersection between two HV-polytopes with the Double Description algorithm. More...
 
static polito_EXPORT int computeIntersectionWithoutCheck (boost::shared_ptr< Polytope_Rn > &A, const boost::shared_ptr< Polytope_Rn > &B) throw (invalid_argument, out_of_range, ios_base::failure, logic_error)
 Compute the intersection between two HV-polytopes with the Double Description algorithm. More...
 
static polito_EXPORT bool isIncluded (const boost::shared_ptr< Polytope_Rn > &A, const boost::shared_ptr< Polytope_Rn > &B) throw (invalid_argument, out_of_range, ios_base::failure, logic_error)
 Test whether the polytope A V-description is inside the polytope B H-description. More...
 
static polito_EXPORT int computeMinkowskiSumOfPolytopes (const boost::shared_ptr< Polytope_Rn > &A, const boost::shared_ptr< Polytope_Rn > &B, boost::shared_ptr< Polytope_Rn > &C) throw (invalid_argument, out_of_range, ios_base::failure, logic_error)
 Compute the Minkowski sum between two HV-polytopes. More...
 
static polito_EXPORT int computeMinkowskiSumOfPolytopes (const std::vector< boost::shared_ptr< Polytope_Rn > > &arrayOfPol, boost::shared_ptr< Polytope_Rn > &C) throw (invalid_argument, out_of_range, ios_base::failure, logic_error)
 Compute the Minkowski sum of several HV-polytopes stored in an array. More...
 
static polito_EXPORT int computeMinkowskiSumOfPolytopes (const boost::shared_ptr< Polytope_Rn > &A, const boost::shared_ptr< Polytope_Rn > &B, boost::shared_ptr< Polytope_Rn > &C, const std::vector< std::vector< int > > &genitorsOfGeneratorsA, const std::vector< std::vector< int > > &genitorsOfGeneratorsB, std::vector< std::vector< int > > &traceGenerators) throw (invalid_argument, out_of_range, ios_base::failure, logic_error)
 Compute the Minkowski sum between two HV-polytopes tracing the generators. More...
 
static polito_EXPORT int checkEqualityOfPolytopes (const boost::shared_ptr< Polytope_Rn > &A, const boost::shared_ptr< Polytope_Rn > &B, bool getFaceMapping=false) throw (invalid_argument, out_of_range, ios_base::failure, logic_error)
 Check whether two HV-polytopes are identical Check whether the vertices of A are inside B half-spaces and vice-versa. Perform also some topological verifications. More...
 
static polito_EXPORT bool checkEqualityOfVertices (const boost::shared_ptr< Polytope_Rn > &A, const boost::shared_ptr< Polytope_Rn > &B) throw (invalid_argument, out_of_range, ios_base::failure, logic_error)
 Check whether two V-polytopes are identical Check whether the sets of vertices of A and B are equal. More...
 
static polito_EXPORT int checkTopologyAndGeometry (const boost::shared_ptr< PolyhedralCone_Rn > &A)
 Check whether a HV-polytopes is correct. More...
 
static polito_EXPORT int makeCube (boost::shared_ptr< Polytope_Rn > &A, double M)
 Create a cube whose vertices will be (+-M, ..., +-M) More...
 
static polito_EXPORT int makeBox (boost::shared_ptr< Polytope_Rn > &A, Point_Rn Pmin, Point_Rn Pmax)
 Create a box or rectangle parallelepiped whose corners are Pmin and Pmax. More...
 
static polito_EXPORT int PolarPolytope (const boost::shared_ptr< Polytope_Rn > &original_pol, boost::shared_ptr< Polytope_Rn > &polar_pol)
 Compute the polar polytope. More...
 
static polito_EXPORT int Translate (boost::shared_ptr< Polytope_Rn > &pol, const boost::numeric::ublas::vector< double > &v2t)
 Translate a polytope or polyhedral cone by the given vector. More...
 
static polito_EXPORT double computeVolume (const boost::shared_ptr< Polytope_Rn > P) throw (invalid_argument, out_of_range, ios_base::failure)
 Return the volume of the given polytope P with its double description. The implemented algorithm can be found in Volume Computation for Polytopes: Strategies and Performances by Andreas Enge in Encyclopedia of Optimization 2nd edition, p 4032-4073. More...
 
static polito_EXPORT int pseudoIntersection (const boost::shared_ptr< Polytope_Rn > &A, const boost::shared_ptr< Polytope_Rn > &B, boost::shared_ptr< Polytope_Rn > &C, const std::set< unsigned int > &firstOperandCaps, const std::set< unsigned int > &secondOperandCaps, std::set< unsigned int > &newCaps, double bb_size=1000.) throw (invalid_argument,out_of_range,ios_base::failure,logic_error)
 Remove all cap half-spaces and then compute the intersection of two capped polytopes. More...
 
static polito_EXPORT int pseudoSum (const boost::shared_ptr< Polytope_Rn > &A, const boost::shared_ptr< Polytope_Rn > &B, boost::shared_ptr< Polytope_Rn > &C, const std::set< unsigned int > &firstOperandCaps, const std::set< unsigned int > &secondOperandCaps, std::set< unsigned int > &newCaps, double bb_size=1000.) throw (invalid_argument,out_of_range,ios_base::failure,logic_error)
 Compute the Minkowski sum of two polytopes and then remove all cap half-spaces to truncate again. More...
 
static polito_EXPORT int computeVoronoiDiagram (const boost::shared_ptr< Polytope_Rn > &inputSpace, const std::vector< Point_Rn > &listOfSeeds, std::vector< boost::shared_ptr< Polytope_Rn > > &VoronoiCells, Voronoi_Rn::TypeOfAlgorithm vorAlgo=Voronoi_Rn::Incremental) throw (std::invalid_argument)
 Compute the Voronoi Diagram in a n-dimensional space i.e. a partitioning of an input space into regions based on distance to points called seeds. More...
 

Detailed Description

The main class to run the basic API.

Definition at line 43 of file politopixAPI.h.

Member Function Documentation

int politopixAPI::addGenerator ( boost::shared_ptr< Polytope_Rn > &  A,
const boost::shared_ptr< Generator_Rn > &  GN 
)
static

Add a new generator into the polytope data structure.

Parameters
AThe current polytope
GNThe corresponding generator
Returns
TEST_OK or 0 if the process was successful, TEST_KO or -1 if something went wrong.

Definition at line 79 of file politopixAPI.cpp.

int politopixAPI::addHalfspace ( boost::shared_ptr< Polytope_Rn > &  A,
const boost::shared_ptr< HalfSpace_Rn > &  HS 
)
static

Add a new half-space into the polytope data structure.

Parameters
AThe current polytope
HSThe corresponding half-space
Returns
TEST_OK or 0 if the process was successful, TEST_KO or -1 if something went wrong.

Definition at line 74 of file politopixAPI.cpp.

int politopixAPI::checkEqualityOfPolytopes ( const boost::shared_ptr< Polytope_Rn > &  A,
const boost::shared_ptr< Polytope_Rn > &  B,
bool  getFaceMapping = false 
)
throw (invalid_argument,
out_of_range,
ios_base::failure,
logic_error
)
static

Check whether two HV-polytopes are identical Check whether the vertices of A are inside B half-spaces and vice-versa. Perform also some topological verifications.

Parameters
AThe 1st HV-polytope
BThe 2nd HV-polytope
getFaceMappingIf true, print the mapping between the generators and faces of both polytopes in case of equality.
Returns
TEST_OK or 0 if the process was successful, TEST_KO or -1 if something went wrong.

Definition at line 473 of file politopixAPI.cpp.

Here is the caller graph for this function:

bool politopixAPI::checkEqualityOfVertices ( const boost::shared_ptr< Polytope_Rn > &  A,
const boost::shared_ptr< Polytope_Rn > &  B 
)
throw (invalid_argument,
out_of_range,
ios_base::failure,
logic_error
)
static

Check whether two V-polytopes are identical Check whether the sets of vertices of A and B are equal.

Parameters
AThe 1st V-polytope
BThe 2nd V-polytope
Returns
true if

\[ \mathcal{V}_A = \mathcal{V}_B \]

, false otherwise.

Definition at line 511 of file politopixAPI.cpp.

Here is the caller graph for this function:

int politopixAPI::checkTopologyAndGeometry ( const boost::shared_ptr< PolyhedralCone_Rn > &  A)
static

Check whether a HV-polytopes is correct.

Parameters
AA HV-polytope
Returns
TEST_OK or 0 if the process was successful, TEST_KO or -1 if something went wrong.

Definition at line 577 of file politopixAPI.cpp.

Here is the caller graph for this function:

int politopixAPI::computeDoubleDescription ( boost::shared_ptr< Polytope_Rn > &  A,
double  bb_size = 1000. 
)
throw (invalid_argument,
out_of_range,
ios_base::failure,
logic_error
)
static

Compute the HV-description for a given H-polytope or V-polytope with the Double Description algorithm.

Parameters
AThe H-polytope or V-polytope
bb_sizeThe origin centered bounding box size providing the V-description
Returns
TEST_OK or 0 if the process was successful, TEST_KO or -1 if something went wrong.

Definition at line 136 of file politopixAPI.cpp.

Here is the caller graph for this function:

int politopixAPI::computeDoubleDescriptionWithoutCheck ( boost::shared_ptr< Polytope_Rn > &  A,
double  bb_size = 1000. 
)
throw (invalid_argument,
out_of_range,
ios_base::failure,
logic_error
)
static

Compute the HV-description for a given H-polytope or V-polytope with the Double Description algorithm.

Parameters
AThe H-polytope or V-polytope
Returns
TEST_OK or 0 if the process was successful, TEST_KO or -1 if something went wrong.

Definition at line 84 of file politopixAPI.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

int politopixAPI::computeIntersection ( const boost::shared_ptr< Polytope_Rn > &  A,
const boost::shared_ptr< Polytope_Rn > &  B,
boost::shared_ptr< Polytope_Rn > &  C 
)
throw (invalid_argument,
out_of_range,
ios_base::failure,
logic_error
)
static

Compute the intersection between two HV-polytopes with the Double Description algorithm.

Parameters
AThe 1st HV-polytope
BThe 2nd HV-polytope
CThe previously allocated result
Returns
TEST_OK or 0 if the process was successful, TEST_KO or -1 if something went wrong.

Definition at line 197 of file politopixAPI.cpp.

Here is the caller graph for this function:

int politopixAPI::computeIntersection ( boost::shared_ptr< Polytope_Rn > &  A,
const boost::shared_ptr< Polytope_Rn > &  B 
)
throw (invalid_argument,
out_of_range,
ios_base::failure,
logic_error
)
static

Compute the intersection between two HV-polytopes with the Double Description algorithm.

Parameters
AThe 1st HV-polytope, then the result
BThe 2nd HV-polytope
Returns
TEST_OK or 0 if the process was successful, TEST_KO or -1 if something went wrong.

Definition at line 144 of file politopixAPI.cpp.

int politopixAPI::computeIntersectionWithoutCheck ( boost::shared_ptr< Polytope_Rn > &  A,
const boost::shared_ptr< Polytope_Rn > &  B 
)
throw (invalid_argument,
out_of_range,
ios_base::failure,
logic_error
)
static

Compute the intersection between two HV-polytopes with the Double Description algorithm.

Parameters
AThe 1st HV-polytope, then the result
BThe 2nd HV-polytope
Returns
TEST_OK or 0 if the process was successful, TEST_KO or -1 if something went wrong.

Definition at line 302 of file politopixAPI.cpp.

int politopixAPI::computeMinkowskiSumOfPolytopes ( const boost::shared_ptr< Polytope_Rn > &  A,
const boost::shared_ptr< Polytope_Rn > &  B,
boost::shared_ptr< Polytope_Rn > &  C 
)
throw (invalid_argument,
out_of_range,
ios_base::failure,
logic_error
)
static

Compute the Minkowski sum between two HV-polytopes.

Parameters
AThe 1st HV-polytope
BThe 2nd HV-polytope
CThe previously allocated sum
Returns
TEST_OK or 0 if the process was successful, TEST_KO or -1 if something went wrong.

Definition at line 355 of file politopixAPI.cpp.

Here is the caller graph for this function:

int politopixAPI::computeMinkowskiSumOfPolytopes ( const std::vector< boost::shared_ptr< Polytope_Rn > > &  arrayOfPol,
boost::shared_ptr< Polytope_Rn > &  C 
)
throw (invalid_argument,
out_of_range,
ios_base::failure,
logic_error
)
static

Compute the Minkowski sum of several HV-polytopes stored in an array.

Parameters
arrayOfPolThe array of HV-polytopes
CThe previously allocated sum
Returns
TEST_OK or 0 if the process was successful, TEST_KO or -1 if something went wrong.

Definition at line 393 of file politopixAPI.cpp.

int politopixAPI::computeMinkowskiSumOfPolytopes ( const boost::shared_ptr< Polytope_Rn > &  A,
const boost::shared_ptr< Polytope_Rn > &  B,
boost::shared_ptr< Polytope_Rn > &  C,
const std::vector< std::vector< int > > &  genitorsOfGeneratorsA,
const std::vector< std::vector< int > > &  genitorsOfGeneratorsB,
std::vector< std::vector< int > > &  traceGenerators 
)
throw (invalid_argument,
out_of_range,
ios_base::failure,
logic_error
)
static

Compute the Minkowski sum between two HV-polytopes tracing the generators.

Parameters
AThe 1st HV-polytope
BThe 2nd HV-polytope
CThe previously allocated sum
traceGeneratorsGive for each generator of C, the list of numbers identifying its genitors in A and B
Returns
TEST_OK or 0 if the process was successful, TEST_KO or -1 if something went wrong.

Definition at line 431 of file politopixAPI.cpp.

double politopixAPI::computeVolume ( const boost::shared_ptr< Polytope_Rn P)
throw (invalid_argument,
out_of_range,
ios_base::failure
)
static

Return the volume of the given polytope P with its double description. The implemented algorithm can be found in Volume Computation for Polytopes: Strategies and Performances by Andreas Enge in Encyclopedia of Optimization 2nd edition, p 4032-4073.

Definition at line 549 of file politopixAPI.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

int politopixAPI::computeVoronoiDiagram ( const boost::shared_ptr< Polytope_Rn > &  inputSpace,
const std::vector< Point_Rn > &  listOfSeeds,
std::vector< boost::shared_ptr< Polytope_Rn > > &  VoronoiCells,
Voronoi_Rn::TypeOfAlgorithm  vorAlgo = Voronoi_Rn::Incremental 
)
throw (std::invalid_argument
)
static

Compute the Voronoi Diagram in a n-dimensional space i.e. a partitioning of an input space into regions based on distance to points called seeds.

Parameters
inputSpaceThe input HV-polytope, most of the times a parallelepiped
listOfSeedsThe list of points to be considered as seeds
VoronoiCellsThe list of the returned HV-polytopes partitioning the input space
vorAlgoThe type of algorithm to be runned: 0 = incremental, 1 = cell by cell
Returns
TEST_OK or 0 if the process was successful, TEST_KO or -1 if something went wrong.

Definition at line 731 of file politopixAPI.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

bool politopixAPI::isIncluded ( const boost::shared_ptr< Polytope_Rn > &  A,
const boost::shared_ptr< Polytope_Rn > &  B 
)
throw (invalid_argument,
out_of_range,
ios_base::failure,
logic_error
)
static

Test whether the polytope A V-description is inside the polytope B H-description.

Parameters
AThe 1st V-polytope
BThe 2nd H-polytope
Returns
true if

\[ A \subset B \]

, false otherwise.

Definition at line 264 of file politopixAPI.cpp.

int politopixAPI::loadPolyhedralCone ( const string &  pathA,
boost::shared_ptr< PolyhedralCone_Rn > &  A 
)
throw (ios_base::failure
)
static

Load a polyhedral cone from the corresponding file name.

Parameters
pathAThe name of the current file
AThe previously allocated polyhedral cone
Returns
TEST_OK or 0 if the process was successful, TEST_KO or -1 if something went wrong.

Definition at line 62 of file politopixAPI.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

int politopixAPI::loadPolytope ( const string &  pathA,
boost::shared_ptr< Polytope_Rn > &  A 
)
throw (ios_base::failure
)
static

Load a polytope from the corresponding file name.

Parameters
pathAThe name of the current file
AThe previously allocated polytope
Returns
TEST_OK or 0 if the process was successful, TEST_KO or -1 if something went wrong.

Definition at line 38 of file politopixAPI.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

int politopixAPI::makeBox ( boost::shared_ptr< Polytope_Rn > &  A,
Point_Rn  Pmin,
Point_Rn  Pmax 
)
static

Create a box or rectangle parallelepiped whose corners are Pmin and Pmax.

Parameters
AThe box variable
Pminlower corner
Pminupper corner
Returns
TEST_OK or 0 if the process was successful, TEST_KO or -1 if something went wrong.

Definition at line 587 of file politopixAPI.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

int politopixAPI::makeCube ( boost::shared_ptr< Polytope_Rn > &  A,
double  M 
)
static

Create a cube whose vertices will be (+-M, ..., +-M)

Parameters
AThe cube variable
MThe cube half side length.
Returns
TEST_OK or 0 if the process was successful, TEST_KO or -1 if something went wrong.

Definition at line 581 of file politopixAPI.cpp.

Here is the caller graph for this function:

int politopixAPI::PolarPolytope ( const boost::shared_ptr< Polytope_Rn > &  original_pol,
boost::shared_ptr< Polytope_Rn > &  polar_pol 
)
static

Compute the polar polytope.

Parameters
original_polThe input polytope
polar_polThe polar polytope
Returns
TEST_OK or 0 if the process was successful, TEST_KO or -1 if something went wrong.

Definition at line 627 of file politopixAPI.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

int politopixAPI::pseudoIntersection ( const boost::shared_ptr< Polytope_Rn > &  A,
const boost::shared_ptr< Polytope_Rn > &  B,
boost::shared_ptr< Polytope_Rn > &  C,
const std::set< unsigned int > &  firstOperandCaps,
const std::set< unsigned int > &  secondOperandCaps,
std::set< unsigned int > &  newCaps,
double  bb_size = 1000. 
)
throw (invalid_argument,
out_of_range,
ios_base::failure,
logic_error
)
static

Remove all cap half-spaces and then compute the intersection of two capped polytopes.

Parameters
AThe 1st HV-polytope
BThe 2nd HV-polytope
CThe previously allocated intersection
Returns
TEST_OK or 0 if the process was successful, TEST_KO or -1 if something went wrong.

Definition at line 635 of file politopixAPI.cpp.

Here is the caller graph for this function:

int politopixAPI::pseudoSum ( const boost::shared_ptr< Polytope_Rn > &  A,
const boost::shared_ptr< Polytope_Rn > &  B,
boost::shared_ptr< Polytope_Rn > &  C,
const std::set< unsigned int > &  firstOperandCaps,
const std::set< unsigned int > &  secondOperandCaps,
std::set< unsigned int > &  newCaps,
double  bb_size = 1000. 
)
throw (invalid_argument,
out_of_range,
ios_base::failure,
logic_error
)
static

Compute the Minkowski sum of two polytopes and then remove all cap half-spaces to truncate again.

Parameters
AThe 1st HV-polytope
BThe 2nd HV-polytope
CThe previously allocated sum
Returns
TEST_OK or 0 if the process was successful, TEST_KO or -1 if something went wrong.

Definition at line 683 of file politopixAPI.cpp.

Here is the caller graph for this function:

int politopixAPI::savePolyhedralCone ( const string &  pathA,
boost::shared_ptr< PolyhedralCone_Rn > &  A 
)
throw (ios_base::failure
)
static

Save a polyhedral cone in the corresponding file name.

Parameters
pathAThe name of the current file
AThe corresponding polyhedral cone
Returns
TEST_OK or 0 if the process was successful, TEST_KO or -1 if something went wrong.

Definition at line 50 of file politopixAPI.cpp.

Here is the call graph for this function:

int politopixAPI::savePolytope ( const string &  pathA,
boost::shared_ptr< Polytope_Rn > &  A 
)
throw (ios_base::failure
)
static

Save a polytope in the corresponding file name.

Parameters
pathAThe name of the current file
AThe corresponding polytope
Returns
TEST_OK or 0 if the process was successful, TEST_KO or -1 if something went wrong.

Definition at line 26 of file politopixAPI.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

int politopixAPI::Translate ( boost::shared_ptr< Polytope_Rn > &  pol,
const boost::numeric::ublas::vector< double > &  v2t 
)
static

Translate a polytope or polyhedral cone by the given vector.

Parameters
polThe corresponding polytope
v2tThe translation vector
Returns
TEST_OK or 0 if the process was successful, TEST_KO or -1 if something went wrong.

Definition at line 631 of file politopixAPI.cpp.

Here is the call graph for this function:


The documentation for this class was generated from the following files: