22 #include <boost/timer.hpp>
31 catch(ios_base::failure& except) {
32 cerr <<
"In/out exception in politopixAPI::savePolytope(const string&, boost::shared_ptr<Polytope_Rn>&): " << except.what() << endl;
43 catch(ios_base::failure& except) {
44 cerr <<
"In/out exception in politopixAPI::loadPolytope(const string&, boost::shared_ptr<Polytope_Rn>&): " << except.what() << endl;
55 catch(ios_base::failure& except) {
56 cerr <<
"In/out exception in politopixAPI::savePolyhedralCone(const string&, boost::shared_ptr<PolyhedralCone_Rn>&): " << except.what() << endl;
67 catch(ios_base::failure& except) {
68 cerr <<
"In/out exception in politopixAPI::loadPolyhedralCone(const string&, boost::shared_ptr<PolyhedralCone_Rn>&): " << except.what() << endl;
86 boost::timer this_timer;
89 if (A->numberOfGenerators() == 0) {
91 A->createBoundingBox(bb_size);
95 unsigned int truncationStep = 2*A->dimension();
97 boost::shared_ptr<PolyhedralCone_Rn>,
100 DD(A, lexmin_ite, NRP, truncationStep);
102 else if (A->numberOfHalfSpaces() == 0) {
106 catch(invalid_argument& except) {
107 cerr <<
"TIME=" << this_timer.elapsed() << endl;
108 cerr <<
"Invalid argument exception in politopixAPI::computeDoubleDescriptionWithoutCheck() " << except.what() << endl;
111 catch(out_of_range& except) {
112 cerr <<
"TIME=" << this_timer.elapsed() << endl;
113 cerr <<
"Out of range exception in politopixAPI::computeDoubleDescriptionWithoutCheck() " << except.what() << endl;
116 catch(ios_base::failure& except) {
117 cerr <<
"TIME=" << this_timer.elapsed() << endl;
118 cerr <<
"In/out exception in politopixAPI::computeDoubleDescriptionWithoutCheck() " << except.what() << endl;
121 catch(logic_error& except) {
122 cerr <<
"TIME=" << this_timer.elapsed() << endl;
123 cerr <<
"Logic error exception in politopixAPI::computeDoubleDescriptionWithoutCheck() " << except.what() << endl;
127 cerr <<
"TIME=" << this_timer.elapsed() << endl;
128 cerr <<
"Unexpected exception caught in politopixAPI::computeDoubleDescriptionWithoutCheck() !" << endl;
131 cout <<
"TIME=" << this_timer.elapsed() << endl;
138 int res = computeDoubleDescriptionWithoutCheck(A, bb_size);
145 const boost::shared_ptr<Polytope_Rn>& B)
throw (invalid_argument,out_of_range,ios_base::failure,logic_error)
148 unsigned int truncationStep = A->numberOfHalfSpaces();
149 for (iteHSB.begin(); iteHSB.end()!=
true; iteHSB.next()) {
150 A->addHalfSpace(iteHSB.current());
152 boost::timer this_timer;
157 boost::shared_ptr<PolyhedralCone_Rn>,
160 DD(A, lexmin_ite, NRP, truncationStep);
162 catch(invalid_argument& except) {
163 cerr <<
"TIME=" << this_timer.elapsed() << endl;
164 cerr <<
"Invalid argument exception in politopixAPI::computeIntersection() " << except.what() << endl;
167 catch(out_of_range& except) {
168 cerr <<
"TIME=" << this_timer.elapsed() << endl;
169 cerr <<
"Out of range exception in politopixAPI::computeIntersection() " << except.what() << endl;
172 catch(ios_base::failure& except) {
173 cerr <<
"TIME=" << this_timer.elapsed() << endl;
174 cerr <<
"In/out exception in politopixAPI::computeIntersection() " << except.what() << endl;
177 catch(logic_error& except) {
178 cerr <<
"TIME=" << this_timer.elapsed() << endl;
179 cerr <<
"Logic error exception in politopixAPI::computeIntersection() " << except.what() << endl;
183 cerr <<
"TIME=" << this_timer.elapsed() << endl;
184 cerr <<
"Unexpected exception caught in politopixAPI::computeIntersection() !" << endl;
187 cout <<
"TIME=" << this_timer.elapsed() << endl;
189 if (A->numberOfGenerators()==0 && A->numberOfHalfSpaces()==0) {
190 cout <<
"The result is empty." << std::endl;
198 boost::shared_ptr<Polytope_Rn>& C)
throw (invalid_argument,out_of_range,ios_base::failure,logic_error)
202 for (iteHS0.begin(); iteHS0.end()!=
true; iteHS0.next()) {
203 C->addGenerator(iteHS0.current());
207 for (iteHS1.begin(); iteHS1.end()!=
true; iteHS1.next()) {
208 C->addHalfSpace(iteHS1.current());
212 for (iteHS2.begin(); iteHS2.end()!=
true; iteHS2.next()) {
213 C->addHalfSpace(iteHS2.current());
216 boost::timer this_timer;
219 unsigned int truncationStep = A->numberOfHalfSpaces();
224 boost::shared_ptr<PolyhedralCone_Rn>,
227 DD(C, lexmin_ite, NRP, truncationStep);
229 catch(invalid_argument& except) {
230 cerr <<
"TIME=" << this_timer.elapsed() << endl;
231 cerr <<
"Invalid argument exception in politopixAPI::computeIntersection() " << except.what() << endl;
234 catch(out_of_range& except) {
235 cerr <<
"TIME=" << this_timer.elapsed() << endl;
236 cerr <<
"Out of range exception in politopixAPI::computeIntersection() " << except.what() << endl;
239 catch(ios_base::failure& except) {
240 cerr <<
"TIME=" << this_timer.elapsed() << endl;
241 cerr <<
"In/out exception in politopixAPI::computeIntersection() " << except.what() << endl;
244 catch(logic_error& except) {
245 cerr <<
"TIME=" << this_timer.elapsed() << endl;
246 cerr <<
"Logic error exception in politopixAPI::computeIntersection() " << except.what() << endl;
250 cerr <<
"TIME=" << this_timer.elapsed() << endl;
251 cerr <<
"Unexpected exception caught in politopixAPI::computeIntersection() !" << endl;
254 cout <<
"TIME=" << this_timer.elapsed() << endl;
256 if (C->numberOfGenerators()==0 && C->numberOfHalfSpaces()==0) {
257 cout <<
"The result is empty." << std::endl;
265 const boost::shared_ptr<Polytope_Rn>& B)
throw (invalid_argument,out_of_range,ios_base::failure,logic_error)
267 boost::timer this_timer;
268 bool isInside =
false;
270 isInside = A->isIncluded(B);
272 catch(invalid_argument& except) {
273 cerr <<
"TIME=" << this_timer.elapsed() << endl;
274 cerr <<
"Invalid argument exception in politopixAPI::computeIntersection() " << except.what() << endl;
277 catch(out_of_range& except) {
278 cerr <<
"TIME=" << this_timer.elapsed() << endl;
279 cerr <<
"Out of range exception in politopixAPI::computeIntersection() " << except.what() << endl;
282 catch(ios_base::failure& except) {
283 cerr <<
"TIME=" << this_timer.elapsed() << endl;
284 cerr <<
"In/out exception in politopixAPI::computeIntersection() " << except.what() << endl;
287 catch(logic_error& except) {
288 cerr <<
"TIME=" << this_timer.elapsed() << endl;
289 cerr <<
"Logic error exception in politopixAPI::computeIntersection() " << except.what() << endl;
293 cerr <<
"TIME=" << this_timer.elapsed() << endl;
294 cerr <<
"Unexpected exception caught in politopixAPI::computeIntersection() !" << endl;
297 cout <<
"TIME=" << this_timer.elapsed() << endl;
303 const boost::shared_ptr<Polytope_Rn>& B)
throw (invalid_argument,out_of_range,ios_base::failure,logic_error)
306 unsigned int truncationStep = A->numberOfHalfSpaces();
307 for (iteHSB.begin(); iteHSB.end()!=
true; iteHSB.next()) {
308 A->addHalfSpace(iteHSB.current());
310 boost::timer this_timer;
315 boost::shared_ptr<PolyhedralCone_Rn>,
318 DD(A, lexmin_ite, NRP, truncationStep);
320 catch(invalid_argument& except) {
321 cerr <<
"TIME=" << this_timer.elapsed() << endl;
322 cerr <<
"Invalid argument exception in politopixAPI::computeIntersectionWithoutCheck() " << except.what() << endl;
325 catch(out_of_range& except) {
326 cerr <<
"TIME=" << this_timer.elapsed() << endl;
327 cerr <<
"Out of range exception in politopixAPI::computeIntersectionWithoutCheck() " << except.what() << endl;
330 catch(ios_base::failure& except) {
331 cerr <<
"TIME=" << this_timer.elapsed() << endl;
332 cerr <<
"In/out exception in politopixAPI::computeIntersectionWithoutCheck() " << except.what() << endl;
335 catch(logic_error& except) {
336 cerr <<
"TIME=" << this_timer.elapsed() << endl;
337 cerr <<
"Logic error exception in politopixAPI::computeIntersectionWithoutCheck() " << except.what() << endl;
341 cerr <<
"TIME=" << this_timer.elapsed() << endl;
342 cerr <<
"Unexpected exception caught in politopixAPI::computeIntersectionWithoutCheck() !" << endl;
347 if (A->numberOfGenerators()==0 && A->numberOfHalfSpaces()==0) {
356 const boost::shared_ptr<Polytope_Rn>& B, boost::shared_ptr<Polytope_Rn>& C)
357 throw (invalid_argument,out_of_range,ios_base::failure,logic_error)
359 boost::timer this_timer;
363 catch(invalid_argument& except) {
364 cerr <<
"TIME=" << this_timer.elapsed() << endl;
365 cerr <<
"Invalid argument exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
368 catch(out_of_range& except) {
369 cerr <<
"TIME=" << this_timer.elapsed() << endl;
370 cerr <<
"Out of range exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
373 catch(ios_base::failure& except) {
374 cerr <<
"TIME=" << this_timer.elapsed() << endl;
375 cerr <<
"In/out exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
378 catch(logic_error& except) {
379 cerr <<
"TIME=" << this_timer.elapsed() << endl;
380 cerr <<
"Logic error exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
384 cerr <<
"TIME=" << this_timer.elapsed() << endl;
385 cerr <<
"Unexpected exception caught in politopixAPI::computeMinkowskiSumOfPolytopes() !" << endl;
388 cout <<
"TIME=" << this_timer.elapsed() << endl;
394 const std::vector< boost::shared_ptr<Polytope_Rn> >& arrayOfPol,
395 boost::shared_ptr<Polytope_Rn>& C)
throw (invalid_argument, out_of_range, ios_base::failure, logic_error)
397 boost::timer this_timer;
401 catch(invalid_argument& except) {
402 cerr <<
"TIME=" << this_timer.elapsed() << endl;
403 cerr <<
"Invalid argument exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
406 catch(out_of_range& except) {
407 cerr <<
"TIME=" << this_timer.elapsed() << endl;
408 cerr <<
"Out of range exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
411 catch(ios_base::failure& except) {
412 cerr <<
"TIME=" << this_timer.elapsed() << endl;
413 cerr <<
"In/out exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
416 catch(logic_error& except) {
417 cerr <<
"TIME=" << this_timer.elapsed() << endl;
418 cerr <<
"Logic error exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
422 cerr <<
"TIME=" << this_timer.elapsed() << endl;
423 cerr <<
"Unexpected exception caught in politopixAPI::computeMinkowskiSumOfPolytopes() !" << endl;
426 cout <<
"TIME=" << this_timer.elapsed() << endl;
432 const boost::shared_ptr<Polytope_Rn>& B,
433 boost::shared_ptr<Polytope_Rn>& C,
434 const std::vector< std::vector<int> >& genitorsOfGeneratorsA,
435 const std::vector< std::vector<int> >& genitorsOfGeneratorsB,
436 std::vector< std::vector<int> >& traceGenerators)
437 throw (invalid_argument,out_of_range,ios_base::failure,logic_error)
439 boost::timer this_timer;
441 MinkowskiSum Ope(A,B,C,genitorsOfGeneratorsA,genitorsOfGeneratorsB,traceGenerators);
443 catch(invalid_argument& except) {
444 cerr <<
"TIME=" << this_timer.elapsed() << endl;
445 cerr <<
"Invalid argument exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
448 catch(out_of_range& except) {
449 cerr <<
"TIME=" << this_timer.elapsed() << endl;
450 cerr <<
"Out of range exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
453 catch(ios_base::failure& except) {
454 cerr <<
"TIME=" << this_timer.elapsed() << endl;
455 cerr <<
"In/out exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
458 catch(logic_error& except) {
459 cerr <<
"TIME=" << this_timer.elapsed() << endl;
460 cerr <<
"Logic error exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
464 cerr <<
"TIME=" << this_timer.elapsed() << endl;
465 cerr <<
"Unexpected exception caught in politopixAPI::computeMinkowskiSumOfPolytopes() !" << endl;
468 cout <<
"TIME=" << this_timer.elapsed() << endl;
474 const boost::shared_ptr<Polytope_Rn>& B,
bool getFaceMapping)
throw (invalid_argument,out_of_range,ios_base::failure,logic_error)
477 boost::timer this_timer;
479 ret_val = A->checkEquality(B, getFaceMapping);
481 catch(invalid_argument& except) {
482 cerr <<
"TIME=" << this_timer.elapsed() << endl;
483 cerr <<
"Invalid argument exception in politopixAPI::checkEqualityOfPolytopes() " << except.what() << endl;
486 catch(out_of_range& except) {
487 cerr <<
"TIME=" << this_timer.elapsed() << endl;
488 cerr <<
"Out of range exception in politopixAPI::checkEqualityOfPolytopes() " << except.what() << endl;
491 catch(ios_base::failure& except) {
492 cerr <<
"TIME=" << this_timer.elapsed() << endl;
493 cerr <<
"In/out exception in politopixAPI::checkEqualityOfPolytopes() " << except.what() << endl;
496 catch(logic_error& except) {
497 cerr <<
"TIME=" << this_timer.elapsed() << endl;
498 cerr <<
"Logic error exception in politopixAPI::checkEqualityOfPolytopes() " << except.what() << endl;
502 cerr <<
"TIME=" << this_timer.elapsed() << endl;
503 cerr <<
"Unexpected exception caught in politopixAPI::checkEqualityOfPolytopes() " << endl;
506 cout <<
"TIME=" << this_timer.elapsed() << endl;
512 const boost::shared_ptr<Polytope_Rn>& B)
throw (invalid_argument,out_of_range,ios_base::failure,logic_error)
515 boost::timer this_timer;
517 ret_val = A->checkEqualityOfVertices(B);
519 catch(invalid_argument& except) {
520 cerr <<
"TIME=" << this_timer.elapsed() << endl;
521 cerr <<
"Invalid argument exception in politopixAPI::checkEqualityOfVertices() " << except.what() << endl;
524 catch(out_of_range& except) {
525 cerr <<
"TIME=" << this_timer.elapsed() << endl;
526 cerr <<
"Out of range exception in politopixAPI::checkEqualityOfVertices() " << except.what() << endl;
529 catch(ios_base::failure& except) {
530 cerr <<
"TIME=" << this_timer.elapsed() << endl;
531 cerr <<
"In/out exception in politopixAPI::checkEqualityOfVertices() " << except.what() << endl;
534 catch(logic_error& except) {
535 cerr <<
"TIME=" << this_timer.elapsed() << endl;
536 cerr <<
"Logic error exception in politopixAPI::checkEqualityOfVertices() " << except.what() << endl;
540 cerr <<
"TIME=" << this_timer.elapsed() << endl;
541 cerr <<
"Unexpected exception caught in politopixAPI::checkEqualityOfVertices() " << endl;
544 cout <<
"TIME=" << this_timer.elapsed() << endl;
552 boost::timer this_timer;
556 catch(std::invalid_argument& e) {
557 std::cerr <<
"politopixAPI::computeVolume() : invalid argument exception " << e.what() << std::endl;
560 catch(std::out_of_range& e) {
561 std::cerr <<
"politopixAPI::computeVolume() : out of range exception " << e.what() << std::endl;
564 catch(std::ios_base::failure& e) {
565 std::cerr <<
"politopixAPI::computeVolume() : in/out exception " << e.what() << std::endl;
569 std::cerr <<
"politopixAPI::computeVolume() : unexpected exception caught !" << std::endl;
572 cout <<
"TIME=" << this_timer.elapsed() << endl;
583 A->createBoundingBox(M);
591 {
for (
unsigned int i=0; i<dimension; ++i) {
599 {
for (
unsigned int i=0; i<dimension; ++i) {
600 boost::shared_ptr<HalfSpace_Rn> HS;
603 for (
unsigned int coord_count=0; coord_count<dimension; coord_count++) {
604 if (coord_count == i)
605 HS->setCoefficient(coord_count, 1.);
607 HS->setCoefficient(coord_count, 0.);
612 {
for (
unsigned int i=0; i<dimension; ++i) {
613 boost::shared_ptr<HalfSpace_Rn> HS;
616 for (
unsigned int coord_count=0; coord_count<dimension; coord_count++) {
617 if (coord_count == i)
618 HS->setCoefficient(coord_count, -1.);
620 HS->setCoefficient(coord_count, 0.);
636 const boost::shared_ptr<Polytope_Rn>& A,
637 const boost::shared_ptr<Polytope_Rn>& B,
638 boost::shared_ptr<Polytope_Rn>& C,
639 const std::set< unsigned int >& firstOperandCaps,
640 const std::set< unsigned int >& secondOperandCaps,
641 std::set< unsigned int >& newCaps,
642 double bb_size)
throw (invalid_argument,out_of_range,ios_base::failure,logic_error)
644 boost::timer this_timer;
648 catch(invalid_argument& except) {
649 cerr <<
"TIME=" << this_timer.elapsed() << endl;
650 cerr <<
"Invalid argument exception in politopixAPI::PseudoIntersection() " << except.what() << endl;
653 catch(out_of_range& except) {
654 cerr <<
"TIME=" << this_timer.elapsed() << endl;
655 cerr <<
"Out of range exception in politopixAPI::PseudoIntersection() " << except.what() << endl;
658 catch(ios_base::failure& except) {
659 cerr <<
"TIME=" << this_timer.elapsed() << endl;
660 cerr <<
"In/out exception in politopixAPI::PseudoIntersection() " << except.what() << endl;
663 catch(logic_error& except) {
664 cerr <<
"TIME=" << this_timer.elapsed() << endl;
665 cerr <<
"Logic error exception in politopixAPI::PseudoIntersection() " << except.what() << endl;
669 cerr <<
"TIME=" << this_timer.elapsed() << endl;
670 cerr <<
"Unexpected exception caught in politopixAPI::PseudoIntersection() !" << endl;
673 cout <<
"TIME=" << this_timer.elapsed() << endl;
675 if (C->numberOfGenerators()==0 && C->numberOfHalfSpaces()==0) {
676 cout <<
"The result is empty." << std::endl;
684 const boost::shared_ptr<Polytope_Rn>& A,
685 const boost::shared_ptr<Polytope_Rn>& B,
686 boost::shared_ptr<Polytope_Rn>& C,
687 const std::set< unsigned int >& firstOperandCaps,
688 const std::set< unsigned int >& secondOperandCaps,
689 std::set< unsigned int >& newCaps,
690 double bb_size)
throw (invalid_argument,out_of_range,ios_base::failure,logic_error)
692 boost::timer this_timer;
696 catch(invalid_argument& except) {
697 cerr <<
"TIME=" << this_timer.elapsed() << endl;
698 cerr <<
"Invalid argument exception in politopixAPI::pseudoSum() " << except.what() << endl;
701 catch(out_of_range& except) {
702 cerr <<
"TIME=" << this_timer.elapsed() << endl;
703 cerr <<
"Out of range exception in politopixAPI::pseudoSum() " << except.what() << endl;
706 catch(ios_base::failure& except) {
707 cerr <<
"TIME=" << this_timer.elapsed() << endl;
708 cerr <<
"In/out exception in politopixAPI::pseudoSum() " << except.what() << endl;
711 catch(logic_error& except) {
712 cerr <<
"TIME=" << this_timer.elapsed() << endl;
713 cerr <<
"Logic error exception in politopixAPI::pseudoSum() " << except.what() << endl;
717 cerr <<
"TIME=" << this_timer.elapsed() << endl;
718 cerr <<
"Unexpected exception caught in politopixAPI::pseudoSum() !" << endl;
721 cout <<
"TIME=" << this_timer.elapsed() << endl;
723 if (C->numberOfGenerators()==0 && C->numberOfHalfSpaces()==0) {
724 cout <<
"The result is empty." << std::endl;
732 const boost::shared_ptr<Polytope_Rn>& inputSpace,
733 const std::vector<Point_Rn>& listOfSeeds,
734 std::vector< boost::shared_ptr<Polytope_Rn> >& VoronoiCells,
737 boost::timer this_timer;
738 boost::shared_ptr<Voronoi_Rn> VD;
777 throw std::invalid_argument(
"politopixAPI::computeVoronoiDiagram() vorAlgo out of range");
781 VoronoiCells = VD->getVoronoiCells();
783 catch(length_error& except) {
784 cerr <<
"TIME=" << this_timer.elapsed() << endl;
785 cerr <<
"Length error exception in politopixAPI::computeVoronoiDiagram() " << except.what() << endl;
789 cerr <<
"TIME=" << this_timer.elapsed() << endl;
790 cerr <<
"Unexpected exception caught in politopixAPI::computeVoronoiDiagram() !" << endl;
793 cout <<
"TIME=" << this_timer.elapsed() << endl;
794 return VD->checkTopologyAndGeometry() ==
true ?
TEST_OK :
TEST_KO;
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.
static polito_EXPORT int savePolytope(const string &pathA, boost::shared_ptr< Polytope_Rn > &A)
Save a polytope in the corresponding file name.
Creation of a n-coordinate geometric point designed to be shared by its neighbour faces...
Computes the Voronoi diagram cell by cell.
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.
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)
Compute the intersection between two HV-polytopes with the Double Description algorithm.
Compute the Minkowski sum of two polytopes and then remove all cap half-spaces to truncate again...
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)
Compute the Voronoi Diagram in a n-dimensional space i.e. a partitioning of an input space into regio...
Model a polytope using its two equivalent definitions : the convex hull and the half-space intersecti...
A half-space whose frontier is a linear (n-1) dimension space. _constant + _coefficients[0].x1 + ... + _coefficients[n-1].xn >= 0.
Refers to the class CellByCellDistNgb (100% of the distances between seeds to be sorted) ...
double getCoordinate(unsigned int i) const
static polito_EXPORT double computeVolume(const boost::shared_ptr< Polytope_Rn > P)
Return the volume of the given polytope P with its double description. The implemented algorithm can ...
static polito_EXPORT int checkEqualityOfPolytopes(const boost::shared_ptr< Polytope_Rn > &A, const boost::shared_ptr< Polytope_Rn > &B, bool getFaceMapping=false)
Check whether two HV-polytopes are identical Check whether the vertices of A are inside B half-spaces...
static polito_EXPORT int makeCube(boost::shared_ptr< Polytope_Rn > &A, double M)
Create a cube whose vertices will be (+-M, ..., +-M)
Compute the Minkowski sum of two polytopes.
Do the same as CellByCellVoronoiDist with the neighbourhood computation.
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.
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.)
Compute the Minkowski sum of two polytopes and then remove all cap half-spaces to truncate again...
static polito_EXPORT bool checkEqualityOfVertices(const boost::shared_ptr< Polytope_Rn > &A, const boost::shared_ptr< Polytope_Rn > &B)
Check whether two V-polytopes are identical Check whether the sets of vertices of A and B are equal...
static polito_EXPORT void load(const std::string &filename, boost::shared_ptr< PolyhedralCone_Rn > POLY)
Load the main data format to store polytopes.
static polito_EXPORT void save(const std::string &filename, boost::shared_ptr< PolyhedralCone_Rn > POLY)
Save the polytope to the main data format.
Remove all cap half-spaces and then compute the intersection of two capped polytopes.
Refers to the class CellByCellVoronoi.
static polito_EXPORT int checkTopologyAndGeometry(const boost::shared_ptr< PolyhedralCone_Rn > &A)
Check whether a HV-polytopes is correct.
Refers to the class CellByCellDistNgb ( 1% of the distances between seeds to be sorted) ...
static polito_EXPORT int loadPolytope(const string &pathA, boost::shared_ptr< Polytope_Rn > &A)
Load a polytope from the corresponding file name.
static polito_EXPORT int savePolyhedralCone(const string &pathA, boost::shared_ptr< PolyhedralCone_Rn > &A)
Save a polyhedral cone in the corresponding file name.
static polito_EXPORT int PolarPolytope(const boost::shared_ptr< Polytope_Rn > &original_pol, boost::shared_ptr< Polytope_Rn > &polar_pol)
Compute the polar polytope.
static polito_EXPORT int computeDoubleDescription(boost::shared_ptr< Polytope_Rn > &A, double bb_size=1000.)
Compute the HV-description for a given H-polytope or V-polytope with the Double Description algorithm...
Refers to the class CellByCellDistNgb ( 30% of the distances between seeds to be sorted) ...
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.
static polito_EXPORT unsigned int getDimension()
Return the dimension of the cartesian space we work in.
static double compute(const boost::shared_ptr< Polytope_Rn > P)
Return the volume of the given polytope P.
This class can be more time-consuming than WeakRedundancyProcessing or NoRedundancyProcessing because...
static polito_EXPORT int loadPolyhedralCone(const string &pathA, boost::shared_ptr< PolyhedralCone_Rn > &A)
Load a polyhedral cone from the corresponding file name.
Refers to the class CellByCellVoronoiDist.
Refers to the class IncrementalVoronoi.
This class is designed to run the list of all geometric objects representing a polytope.
static polito_EXPORT int computeDoubleDescriptionWithoutCheck(boost::shared_ptr< Polytope_Rn > &A, double bb_size=1000.)
Compute the HV-description for a given H-polytope or V-polytope with the Double Description algorithm...
static polito_EXPORT int computeIntersectionWithoutCheck(boost::shared_ptr< Polytope_Rn > &A, const boost::shared_ptr< Polytope_Rn > &B)
Compute the intersection between two HV-polytopes with the Double Description algorithm.
TypeOfAlgorithm
Choose the kind of algorithm to be run.
static polito_EXPORT bool isIncluded(const boost::shared_ptr< Polytope_Rn > &A, const boost::shared_ptr< Polytope_Rn > &B)
Test whether the polytope A V-description is inside the polytope B H-description. ...
At the end of the ith iteration, the Voronoi diagram of i seeds is computed.
Refers to the class CellByCellDistNgb ( 10% of the distances between seeds to be sorted) ...
void setAverageNumberOfNeighbours(double pc)
When we want to do a partial sort so try to assess the average number of facet neighbours for each ce...
static int Compute(boost::shared_ptr< Polytope_Rn > &pol, double bb_size=1000.)
Use the polarity to get the facets from the generators.
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.)
Remove all cap half-spaces and then compute the intersection of two capped polytopes.
Refers to the class CellByCellDistNgb ( 20% of the distances between seeds to be sorted) ...
The algorithm implemented here is an incremental algorithm as mentioned in How Good are Convex Hull ...
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)
Compute the Minkowski sum between two HV-polytopes.