politopix  5.0.0
politopixAPI.h
Go to the documentation of this file.
1 // politopix allows to make computations on polytopes such as finding vertices, intersecting, Minkowski sums, ...
2 // Copyright (C) 2014-2020 : Delos Vincent
3 //
4 // This program is free software: you can redistribute it and/or modify
5 // it under the terms of the GNU General Public License as published by
6 // the Free Software Foundation, either version 3 of the License, or
7 // (at your option) any later version.
8 //
9 // This program is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 // GNU General Public License for more details.
13 //
14 // You should have received a copy of the GNU General Public License
15 // along with this program. If not, see <http://www.gnu.org/licenses/>.
16 //
19 // I2M (UMR CNRS 5295 / University of Bordeaux)
20 
21 #ifndef INC_polAPI
22 #define INC_polAPI
23 
24 #include "polito_Export.h"
25 #include "Rn.h"
26 #include "Voronoi_Rn.h"
27 #include "Polytope_Rn.h"
28 #include "IO_Polytope.h"
29 #include "Generator_Rn.h"
30 #include "HalfSpace_Rn.h"
31 #include "VolumeOfPolytopes_Rn.h"
32 #include "DoubleDescription_Rn.h"
35 
36 using namespace std;
37 
38 #define TEST_OK 0
39 #define TEST_KO -1
40 
43 class politopixAPI {
44 
45  public:
46 
51  static polito_EXPORT int savePolytope(const string& pathA, boost::shared_ptr<Polytope_Rn>& A);
52 
57  static polito_EXPORT int savePolyhedralCone(const string& pathA, boost::shared_ptr<PolyhedralCone_Rn>& A);
58 
63  static polito_EXPORT int loadPolytope(const string& pathA, boost::shared_ptr<Polytope_Rn>& A);
64 
69  static polito_EXPORT int loadPolyhedralCone(const string& pathA, boost::shared_ptr<PolyhedralCone_Rn>& A);
70 
75  static polito_EXPORT int addHalfspace(boost::shared_ptr<Polytope_Rn>& A, const boost::shared_ptr<HalfSpace_Rn>& HS);
76 
81  static polito_EXPORT int addGenerator(boost::shared_ptr<Polytope_Rn>& A, const boost::shared_ptr<Generator_Rn>& GN);
82 
87  static polito_EXPORT int computeDoubleDescription(boost::shared_ptr<Polytope_Rn>& A, double bb_size);
88 
90  static polito_EXPORT int computeDoubleDescription(boost::shared_ptr<Polytope_Rn>& A) { return computeDoubleDescription(A,1000.); }
91 
95  static polito_EXPORT int computeDoubleDescriptionWithoutCheck(boost::shared_ptr<Polytope_Rn>& A, double bb_size=1000.);
96 
102  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);
103 
108  static polito_EXPORT int computeIntersection(boost::shared_ptr<Polytope_Rn>& A, const boost::shared_ptr<Polytope_Rn>& B);
109 
114  static polito_EXPORT int computeIntersectionWithoutCheck(boost::shared_ptr<Polytope_Rn>& A, const boost::shared_ptr<Polytope_Rn>& B);
115 
120  static polito_EXPORT bool isIncluded(const boost::shared_ptr<Polytope_Rn>& A, const boost::shared_ptr<Polytope_Rn>& B);
121 
127  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);
128 
133  static polito_EXPORT int computeMinkowskiSumOfPolytopes(const std::vector< boost::shared_ptr<Polytope_Rn> >& arrayOfPol, boost::shared_ptr<Polytope_Rn>& C);
134 
141  static polito_EXPORT int computeMinkowskiSumOfPolytopes(const boost::shared_ptr<Polytope_Rn>& A,
142  const boost::shared_ptr<Polytope_Rn>& B,
143  boost::shared_ptr<Polytope_Rn>& C,
144  const std::vector< std::vector<int> >& genitorsOfGeneratorsA,
145  const std::vector< std::vector<int> >& genitorsOfGeneratorsB,
146  std::vector< std::vector<int> >& traceGenerators);
147 
154  static polito_EXPORT int checkEqualityOfPolytopes(const boost::shared_ptr<Polytope_Rn>& A,const boost::shared_ptr<Polytope_Rn>& B, bool getFaceMapping=false);
155 
161  static polito_EXPORT bool checkEqualityOfVertices(const boost::shared_ptr<Polytope_Rn>& A, const boost::shared_ptr<Polytope_Rn>& B);
162 
166  static polito_EXPORT int checkTopologyAndGeometry(const boost::shared_ptr<PolyhedralCone_Rn>& A, bool check_all=false);
167 
172  static polito_EXPORT int makeCube(boost::shared_ptr<Polytope_Rn>& A, double M);
173 
179  static polito_EXPORT int makeBox(boost::shared_ptr<Polytope_Rn>& A, Point_Rn Pmin, Point_Rn Pmax);
180 
185  static polito_EXPORT int PolarPolytope(const boost::shared_ptr<Polytope_Rn>& original_pol, boost::shared_ptr<Polytope_Rn>& polar_pol);
186 
191  static polito_EXPORT int Translate(boost::shared_ptr<Polytope_Rn>& pol, const boost::numeric::ublas::vector<double>& v2t);
192 
196  static polito_EXPORT double computeVolume(const boost::shared_ptr<Polytope_Rn> P);
197 
203  static polito_EXPORT int pseudoIntersection(
204  const boost::shared_ptr<Polytope_Rn>& A,
205  const boost::shared_ptr<Polytope_Rn>& B,
206  boost::shared_ptr<Polytope_Rn>& C,
207  const std::set< unsigned int >& firstOperandCaps,
208  const std::set< unsigned int >& secondOperandCaps,
209  std::set< unsigned int >& newCaps,
210  double bb_size=1000.);
211 
217  static polito_EXPORT int pseudoSum(
218  const boost::shared_ptr<Polytope_Rn>& A,
219  const boost::shared_ptr<Polytope_Rn>& B,
220  boost::shared_ptr<Polytope_Rn>& C,
221  const std::set< unsigned int >& firstOperandCaps,
222  const std::set< unsigned int >& secondOperandCaps,
223  std::set< unsigned int >& newCaps,
224  double bb_size=1000.);
225 
232  static polito_EXPORT int computeVoronoiDiagram(
233  const boost::shared_ptr<Polytope_Rn>& inputSpace,
234  const std::vector<Point_Rn>& listOfSeeds,
235  std::vector< boost::shared_ptr<Polytope_Rn> >& VoronoiCells,
237 
238 };
239 
242 template<class T> class PolyhedralBody {
243 
244  public:
245 
246  static int sum(const boost::shared_ptr<T>& firstBody, const boost::shared_ptr<T>& secondBody, boost::shared_ptr<T>& result) {
247  result = firstBody->sum(secondBody); return TEST_OK;
248  }
249 
250  static int sum(const std::vector< boost::shared_ptr<T> >& allBodies, boost::shared_ptr<T>& result) {
251  result->sum(allBodies); return TEST_OK;
252  }
253 
254  static int intersect(const boost::shared_ptr<T>& firstBody, const boost::shared_ptr<T>& secondBody, boost::shared_ptr<T>& result) {
255  result = firstBody->intersect(secondBody); return TEST_OK;
256  }
257 
258  static int intersect(const std::vector< boost::shared_ptr<T> >& allBodies, boost::shared_ptr<T>& result) {
259  result->intersect(allBodies); return TEST_OK;
260  }
261 
266  static bool isIncluded(const boost::shared_ptr<T>& firstBody, const boost::shared_ptr<T>& secondBody) {
267  return firstBody->isIncluded(secondBody);
268  }
269 
275  static bool isIncluded(const boost::shared_ptr<T>& firstBody, const boost::shared_ptr<T>& secondBody, double& minmaxDistance) {
276  return firstBody->isIncluded(secondBody, minmaxDistance);
277  }
278 
279  static void resetConstantsOfAllHalfSpaces(boost::shared_ptr<T>& body, double cst) {
280  body->resetConstantsOfAllHalfSpaces(cst);
281  }
282 
283  static int scalingFactor(boost::shared_ptr<T>& body, double factor) {
284  return body->scalingFactor(factor);
285  }
286 
287  static bool checkTopologyAndGeometry(const boost::shared_ptr<T>& body) {
288  return body->checkTopologyAndGeometry();
289  }
290 
291  static std::string getName() {
292  return T::getName();
293  }
294 
295  static std::string getFileExtension() {
296  return T::getFileExtension();
297  }
298 
299  static void save(const std::string& pathA, boost::shared_ptr<T>& body) {
300  body->save(pathA);
301  }
302 
303  static void load(const std::string& pathA, boost::shared_ptr<T>& body, double bb_size=10000.) {
304  body->load(pathA, bb_size);
305  }
306 
307  static void dump(std::ostream &current_ostream, boost::shared_ptr<T>& body) {
308  body->dump(current_ostream);
309  }
310 
311 };
312 
313 
314 #endif // INC_polAPI
GeometricObjectIterator_Rn.h
TEST_OK
#define TEST_OK
Definition: politopixAPI.h:38
PolyhedralBody::load
static void load(const std::string &pathA, boost::shared_ptr< T > &body, double bb_size=10000.)
Definition: politopixAPI.h:303
politopixAPI::computeDoubleDescription
static polito_EXPORT int computeDoubleDescription(boost::shared_ptr< Polytope_Rn > &A)
The same as computeDoubleDescription(A,bb_size) with bb_size=1000. for the boost python interface.
Definition: politopixAPI.h:90
PolyhedralBody::sum
static int sum(const boost::shared_ptr< T > &firstBody, const boost::shared_ptr< T > &secondBody, boost::shared_ptr< T > &result)
Definition: politopixAPI.h:246
PolyhedralBody::checkTopologyAndGeometry
static bool checkTopologyAndGeometry(const boost::shared_ptr< T > &body)
Definition: politopixAPI.h:287
polito_EXPORT
#define polito_EXPORT
Definition: polito_Export.h:15
Voronoi_Rn.h
PolyhedralBody::save
static void save(const std::string &pathA, boost::shared_ptr< T > &body)
Definition: politopixAPI.h:299
politopixAPI
The main class to run the basic API.
Definition: politopixAPI.h:43
VolumeOfPolytopes_Rn.h
PolyhedralAlgorithms_Rn.h
PolyhedralBody::resetConstantsOfAllHalfSpaces
static void resetConstantsOfAllHalfSpaces(boost::shared_ptr< T > &body, double cst)
Definition: politopixAPI.h:279
PolyhedralBody::isIncluded
static bool isIncluded(const boost::shared_ptr< T > &firstBody, const boost::shared_ptr< T > &secondBody)
Test whether firstBody is inside secondBody.
Definition: politopixAPI.h:266
PolyhedralBody::sum
static int sum(const std::vector< boost::shared_ptr< T > > &allBodies, boost::shared_ptr< T > &result)
Definition: politopixAPI.h:250
Rn.h
Point_Rn
Creation of a n-coordinate geometric point designed to be shared by its neighbour faces.
Definition: Point_Rn.h:34
Generator_Rn.h
PolyhedralBody::getFileExtension
static std::string getFileExtension()
Definition: politopixAPI.h:295
PolyhedralBody::dump
static void dump(std::ostream &current_ostream, boost::shared_ptr< T > &body)
Definition: politopixAPI.h:307
Voronoi_Rn::Incremental
@ Incremental
Definition: Voronoi_Rn.h:45
PolyhedralBody::isIncluded
static bool isIncluded(const boost::shared_ptr< T > &firstBody, const boost::shared_ptr< T > &secondBody, double &minmaxDistance)
Test whether firstBody is inside secondBody.
Definition: politopixAPI.h:275
DoubleDescription_Rn.h
PolyhedralBody::getName
static std::string getName()
Definition: politopixAPI.h:291
Polytope_Rn.h
Voronoi_Rn::TypeOfAlgorithm
TypeOfAlgorithm
Choose the kind of algorithm to be run.
Definition: Voronoi_Rn.h:44
IO_Polytope.h
computeIntersection
int(* computeIntersection)(const boost::shared_ptr< Polytope_Rn > &, const boost::shared_ptr< Polytope_Rn > &, boost::shared_ptr< Polytope_Rn > &)
Definition: politopy.h:38
PolyhedralBody::intersect
static int intersect(const std::vector< boost::shared_ptr< T > > &allBodies, boost::shared_ptr< T > &result)
Definition: politopixAPI.h:258
PolyhedralBody::scalingFactor
static int scalingFactor(boost::shared_ptr< T > &body, double factor)
Definition: politopixAPI.h:283
PolyhedralBody
Generic class to run computations whether we have, polytopes, capped polytopes, prismatic polyhedra,...
Definition: politopixAPI.h:242
polito_Export.h
HalfSpace_Rn.h
PolyhedralBody::intersect
static int intersect(const boost::shared_ptr< T > &firstBody, const boost::shared_ptr< T > &secondBody, boost::shared_ptr< T > &result)
Definition: politopixAPI.h:254