politopix  5.0.0
politopixAPI.cpp
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) 2015-2021 : 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 #include <stdio.h>
22 #include <boost/timer.hpp>
23 #include "politopixAPI.h"
24 #ifdef BOOST_PYTHON
25 #include "politopy.h"
26 #endif
27 
28 
29 int politopixAPI::savePolytope(const string& pathA, boost::shared_ptr<Polytope_Rn>& A) {
30  try {
31  IO_Polytope::save(pathA, A);
32  }
33  catch(ios_base::failure& except) {
34  cerr << "In/out exception in politopixAPI::savePolytope(const string&, boost::shared_ptr<Polytope_Rn>&): " << except.what() << endl;
35  throw except;
36  }
37  return TEST_OK;
38 }
39 
40 int politopixAPI::loadPolytope(const string& pathA, boost::shared_ptr<Polytope_Rn>& A)
41 {
42  try {
43  IO_Polytope::load(pathA, A);
44  }
45  catch(ios_base::failure& except) {
46  cerr << "In/out exception in politopixAPI::loadPolytope(const string&, boost::shared_ptr<Polytope_Rn>&): " << except.what() << endl;
47  throw except;
48  }
49  return TEST_OK;
50 }
51 
52 int politopixAPI::savePolyhedralCone(const string& pathA, boost::shared_ptr<PolyhedralCone_Rn>& A)
53 {
54  try {
55  IO_Polytope::save(pathA, A);
56  }
57  catch(ios_base::failure& except) {
58  cerr << "In/out exception in politopixAPI::savePolyhedralCone(const string&, boost::shared_ptr<PolyhedralCone_Rn>&): " << except.what() << endl;
59  throw except;
60  }
61  return TEST_OK;
62 }
63 
64 int politopixAPI::loadPolyhedralCone(const string& pathA, boost::shared_ptr<PolyhedralCone_Rn>& A) {
65  try {
66  IO_Polytope::load(pathA, A);
67  }
68  catch(ios_base::failure& except) {
69  cerr << "In/out exception in politopixAPI::loadPolyhedralCone(const string&, boost::shared_ptr<PolyhedralCone_Rn>&): " << except.what() << endl;
70  throw except;
71  }
72  return TEST_OK;
73 }
74 
75 int politopixAPI::addHalfspace(boost::shared_ptr<Polytope_Rn>& A, const boost::shared_ptr<HalfSpace_Rn>& HS) {
76  A->addHalfSpace(HS);
77  return TEST_OK;
78 }
79 
80 int politopixAPI::addGenerator(boost::shared_ptr<Polytope_Rn>& A, const boost::shared_ptr<Generator_Rn>& GN) {
81  A->addGenerator(GN);
82  return TEST_OK;
83 }
84 
85 int politopixAPI::computeDoubleDescriptionWithoutCheck(boost::shared_ptr<Polytope_Rn>& A, double bb_size) {
86  boost::timer this_timer;
87  try {
88  // Do we compute the generators or the half-spaces?
89  if (A->numberOfGenerators() == 0) {
90  // Create the bounding box centered on the origin including the H-polytope A.
91  A->createBoundingBox(bb_size);
93  //StrongRedundancyProcessing< boost::shared_ptr<PolyhedralCone_Rn> > NRP(A->neigbourhoodCondition());
94  // The number of facets for a cube.
95  unsigned int truncationStep = 2*A->dimension();
97  boost::shared_ptr<PolyhedralCone_Rn>,
99  //StrongRedundancyProcessing< boost::shared_ptr<PolyhedralCone_Rn> > >
100  DD(A, lexmin_ite, truncationStep);
101  }
102  else if (A->numberOfHalfSpaces() == 0) {
104  }
105  }
106  catch(invalid_argument& except) {
107  cerr << "TIME=" << this_timer.elapsed() << endl;
108  cerr << "Invalid argument exception in politopixAPI::computeDoubleDescriptionWithoutCheck() " << except.what() << endl;
109  return TEST_KO;
110  }
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;
114  return TEST_KO;
115  }
116  catch(ios_base::failure& except) {
117  cerr << "TIME=" << this_timer.elapsed() << endl;
118  cerr << "In/out exception in politopixAPI::computeDoubleDescriptionWithoutCheck() " << except.what() << endl;
119  return TEST_KO;
120  }
121  catch(logic_error& except) {
122  cerr << "TIME=" << this_timer.elapsed() << endl;
123  cerr << "Logic error exception in politopixAPI::computeDoubleDescriptionWithoutCheck() " << except.what() << endl;
124  return TEST_KO;
125  }
126  catch(...) {
127  cerr << "TIME=" << this_timer.elapsed() << endl;
128  cerr << "Unexpected exception caught in politopixAPI::computeDoubleDescriptionWithoutCheck() !" << endl;
129  return TEST_KO;
130  }
131  cout << "TIME=" << this_timer.elapsed() << endl;
132 
133  return TEST_OK;
134 }
135 
136 int politopixAPI::computeDoubleDescription(boost::shared_ptr<Polytope_Rn>& A, double bb_size) {
137  int res = computeDoubleDescriptionWithoutCheck(A, bb_size);
138  if (res == TEST_KO)
139  return TEST_KO;
140  return A->checkTopologyAndGeometry() == true ? TEST_OK : TEST_KO;
141 }
142 
143 int politopixAPI::computeIntersection(boost::shared_ptr<Polytope_Rn>& A,
144  const boost::shared_ptr<Polytope_Rn>& B) {
146  unsigned int truncationStep = A->numberOfHalfSpaces();
147  for (iteHSB.begin(); iteHSB.end()!=true; iteHSB.next()) {
148  A->addHalfSpace(iteHSB.current());
149  }
150  boost::timer this_timer;
151  try {
152  constIteratorOfListOfGeometricObjects< boost::shared_ptr<HalfSpace_Rn> > lexmin_ite(A->getListOfHalfSpaces());
153  //StrongRedundancyProcessing< boost::shared_ptr<PolyhedralCone_Rn> > NRP(A->neigbourhoodCondition());
155  boost::shared_ptr<PolyhedralCone_Rn>,
157  //StrongRedundancyProcessing< boost::shared_ptr<PolyhedralCone_Rn> > >
158  DD(A, lexmin_ite, truncationStep);
159  }
160  catch(invalid_argument& except) {
161  cerr << "TIME=" << this_timer.elapsed() << endl;
162  cerr << "Invalid argument exception in politopixAPI::computeIntersection() " << except.what() << endl;
163  return TEST_KO;
164  }
165  catch(out_of_range& except) {
166  cerr << "TIME=" << this_timer.elapsed() << endl;
167  cerr << "Out of range exception in politopixAPI::computeIntersection() " << except.what() << endl;
168  return TEST_KO;
169  }
170  catch(ios_base::failure& except) {
171  cerr << "TIME=" << this_timer.elapsed() << endl;
172  cerr << "In/out exception in politopixAPI::computeIntersection() " << except.what() << endl;
173  return TEST_KO;
174  }
175  catch(logic_error& except) {
176  cerr << "TIME=" << this_timer.elapsed() << endl;
177  cerr << "Logic error exception in politopixAPI::computeIntersection() " << except.what() << endl;
178  return TEST_KO;
179  }
180  catch(...) {
181  cerr << "TIME=" << this_timer.elapsed() << endl;
182  cerr << "Unexpected exception caught in politopixAPI::computeIntersection() !" << endl;
183  return TEST_KO;
184  }
185  cout << "TIME=" << this_timer.elapsed() << endl;
186 
187  if (A->numberOfGenerators()==0 && A->numberOfHalfSpaces()==0) {
188  cout << "The result is empty." << std::endl;
189  return TEST_OK;
190  }
191 
192  return A->checkTopologyAndGeometry() == true ? TEST_OK : TEST_KO;
193 }
194 
195 int politopixAPI::computeIntersection(const boost::shared_ptr<Polytope_Rn>& A, const boost::shared_ptr<Polytope_Rn>& B,
196  boost::shared_ptr<Polytope_Rn>& C) {
197  // Copy the generators of A.
199  for (iteHS0.begin(); iteHS0.end()!=true; iteHS0.next()) {
200  C->addGenerator(iteHS0.current());
201  }
202  // Copy the half-spaces of A.
204  for (iteHS1.begin(); iteHS1.end()!=true; iteHS1.next()) {
205  C->addHalfSpace(iteHS1.current());
206  }
207  // Copy the half-spaces of B.
209  for (iteHS2.begin(); iteHS2.end()!=true; iteHS2.next()) {
210  C->addHalfSpace(iteHS2.current());
211  }
212 
213  boost::timer this_timer;
214  try {
215  // Now the truncation.
216  unsigned int truncationStep = A->numberOfHalfSpaces();
217  constIteratorOfListOfGeometricObjects< boost::shared_ptr<HalfSpace_Rn> > lexmin_ite(C->getListOfHalfSpaces());
218  //StrongRedundancyProcessing< boost::shared_ptr<PolyhedralCone_Rn> > NRP(C->neigbourhoodCondition());
219  // The number of facets for a cube.
221  boost::shared_ptr<PolyhedralCone_Rn>,
223  //StrongRedundancyProcessing< boost::shared_ptr<PolyhedralCone_Rn> > >
224  DD(C, lexmin_ite, truncationStep);
225  }
226  catch(invalid_argument& except) {
227  cerr << "TIME=" << this_timer.elapsed() << endl;
228  cerr << "Invalid argument exception in politopixAPI::computeIntersection() " << except.what() << endl;
229  return TEST_KO;
230  }
231  catch(out_of_range& except) {
232  cerr << "TIME=" << this_timer.elapsed() << endl;
233  cerr << "Out of range exception in politopixAPI::computeIntersection() " << except.what() << endl;
234  return TEST_KO;
235  }
236  catch(ios_base::failure& except) {
237  cerr << "TIME=" << this_timer.elapsed() << endl;
238  cerr << "In/out exception in politopixAPI::computeIntersection() " << except.what() << endl;
239  return TEST_KO;
240  }
241  catch(logic_error& except) {
242  cerr << "TIME=" << this_timer.elapsed() << endl;
243  cerr << "Logic error exception in politopixAPI::computeIntersection() " << except.what() << endl;
244  return TEST_KO;
245  }
246  catch(...) {
247  cerr << "TIME=" << this_timer.elapsed() << endl;
248  cerr << "Unexpected exception caught in politopixAPI::computeIntersection() !" << endl;
249  return TEST_KO;
250  }
251  cout << "TIME=" << this_timer.elapsed() << endl;
252 
253  if (C->numberOfGenerators()==0 && C->numberOfHalfSpaces()==0) {
254  cout << "The result is empty." << std::endl;
255  return TEST_OK;
256  }
257 
258  return C->checkTopologyAndGeometry() == true ? TEST_OK : TEST_KO;
259 }
260 
261 bool politopixAPI::isIncluded(const boost::shared_ptr<Polytope_Rn>& A,
262  const boost::shared_ptr<Polytope_Rn>& B) {
263  boost::timer this_timer;
264  bool isInside = false;
265  try {
266  isInside = A->isIncluded(B);
267  }
268  catch(invalid_argument& except) {
269  cerr << "TIME=" << this_timer.elapsed() << endl;
270  cerr << "Invalid argument exception in politopixAPI::isIncluded() " << except.what() << endl;
271  return TEST_KO;
272  }
273  catch(out_of_range& except) {
274  cerr << "TIME=" << this_timer.elapsed() << endl;
275  cerr << "Out of range exception in politopixAPI::isIncluded() " << except.what() << endl;
276  return TEST_KO;
277  }
278  catch(ios_base::failure& except) {
279  cerr << "TIME=" << this_timer.elapsed() << endl;
280  cerr << "In/out exception in politopixAPI::isIncluded() " << except.what() << endl;
281  return TEST_KO;
282  }
283  catch(logic_error& except) {
284  cerr << "TIME=" << this_timer.elapsed() << endl;
285  cerr << "Logic error exception in politopixAPI::isIncluded() " << except.what() << endl;
286  return TEST_KO;
287  }
288  catch(...) {
289  cerr << "TIME=" << this_timer.elapsed() << endl;
290  cerr << "Unexpected exception caught in politopixAPI::isIncluded() !" << endl;
291  return TEST_KO;
292  }
293  cout << "TIME=" << this_timer.elapsed() << endl;
294 
295  return isInside;
296 }
297 
298 int politopixAPI::computeIntersectionWithoutCheck(boost::shared_ptr<Polytope_Rn>& A,
299  const boost::shared_ptr<Polytope_Rn>& B) {
301  unsigned int truncationStep = A->numberOfHalfSpaces();
302  for (iteHSB.begin(); iteHSB.end()!=true; iteHSB.next()) {
303  A->addHalfSpace(iteHSB.current());
304  }
305  boost::timer this_timer;
306  try {
307  constIteratorOfListOfGeometricObjects< boost::shared_ptr<HalfSpace_Rn> > lexmin_ite(A->getListOfHalfSpaces());
308  //StrongRedundancyProcessing< boost::shared_ptr<PolyhedralCone_Rn> > NRP(A->neigbourhoodCondition());
310  boost::shared_ptr<PolyhedralCone_Rn>,
312  //StrongRedundancyProcessing< boost::shared_ptr<PolyhedralCone_Rn> > >
313  DD(A, lexmin_ite, truncationStep);
314  }
315  catch(invalid_argument& except) {
316  cerr << "TIME=" << this_timer.elapsed() << endl;
317  cerr << "Invalid argument exception in politopixAPI::computeIntersectionWithoutCheck() " << except.what() << endl;
318  return TEST_KO;
319  }
320  catch(out_of_range& except) {
321  cerr << "TIME=" << this_timer.elapsed() << endl;
322  cerr << "Out of range exception in politopixAPI::computeIntersectionWithoutCheck() " << except.what() << endl;
323  return TEST_KO;
324  }
325  catch(ios_base::failure& except) {
326  cerr << "TIME=" << this_timer.elapsed() << endl;
327  cerr << "In/out exception in politopixAPI::computeIntersectionWithoutCheck() " << except.what() << endl;
328  return TEST_KO;
329  }
330  catch(logic_error& except) {
331  cerr << "TIME=" << this_timer.elapsed() << endl;
332  cerr << "Logic error exception in politopixAPI::computeIntersectionWithoutCheck() " << except.what() << endl;
333  return TEST_KO;
334  }
335  catch(...) {
336  cerr << "TIME=" << this_timer.elapsed() << endl;
337  cerr << "Unexpected exception caught in politopixAPI::computeIntersectionWithoutCheck() !" << endl;
338  return TEST_KO;
339  }
340  //cout << "TIME=" << this_timer.elapsed() << endl;
341 
342  if (A->numberOfGenerators()==0 && A->numberOfHalfSpaces()==0) {
343  //cout << "The result is empty." << std::endl;
344  return TEST_OK;
345  }
346 
347  return TEST_OK;
348 }
349 
350 int politopixAPI::computeMinkowskiSumOfPolytopes(const boost::shared_ptr<Polytope_Rn>& A,
351  const boost::shared_ptr<Polytope_Rn>& B, boost::shared_ptr<Polytope_Rn>& C) {
352  boost::timer this_timer;
353  try {
354  MinkowskiSum Ope(A,B,C);
355  }
356  catch(invalid_argument& except) {
357  cerr << "TIME=" << this_timer.elapsed() << endl;
358  cerr << "Invalid argument exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
359  return TEST_KO;
360  }
361  catch(out_of_range& except) {
362  cerr << "TIME=" << this_timer.elapsed() << endl;
363  cerr << "Out of range exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
364  return TEST_KO;
365  }
366  catch(ios_base::failure& except) {
367  cerr << "TIME=" << this_timer.elapsed() << endl;
368  cerr << "In/out exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
369  return TEST_KO;
370  }
371  catch(logic_error& except) {
372  cerr << "TIME=" << this_timer.elapsed() << endl;
373  cerr << "Logic error exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
374  return TEST_KO;
375  }
376  catch(...) {
377  cerr << "TIME=" << this_timer.elapsed() << endl;
378  cerr << "Unexpected exception caught in politopixAPI::computeMinkowskiSumOfPolytopes() !" << endl;
379  return TEST_KO;
380  }
381  cout << "TIME=" << this_timer.elapsed() << endl;
382 
383  return C->checkTopologyAndGeometry() == true ? TEST_OK : TEST_KO;
384 }
385 
387  const std::vector< boost::shared_ptr<Polytope_Rn> >& arrayOfPol,
388  boost::shared_ptr<Polytope_Rn>& C) {
389  boost::timer this_timer;
390  try {
391  MinkowskiSum Ope(arrayOfPol,C);
392  }
393  catch(invalid_argument& except) {
394  cerr << "TIME=" << this_timer.elapsed() << endl;
395  cerr << "Invalid argument exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
396  return TEST_KO;
397  }
398  catch(out_of_range& except) {
399  cerr << "TIME=" << this_timer.elapsed() << endl;
400  cerr << "Out of range exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
401  return TEST_KO;
402  }
403  catch(ios_base::failure& except) {
404  cerr << "TIME=" << this_timer.elapsed() << endl;
405  cerr << "In/out exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
406  return TEST_KO;
407  }
408  catch(logic_error& except) {
409  cerr << "TIME=" << this_timer.elapsed() << endl;
410  cerr << "Logic error exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
411  return TEST_KO;
412  }
413  catch(...) {
414  cerr << "TIME=" << this_timer.elapsed() << endl;
415  cerr << "Unexpected exception caught in politopixAPI::computeMinkowskiSumOfPolytopes() !" << endl;
416  return TEST_KO;
417  }
418  cout << "TIME=" << this_timer.elapsed() << endl;
419 
420  return C->checkTopologyAndGeometry() == true ? TEST_OK : TEST_KO;
421 }
422 
423 int politopixAPI::computeMinkowskiSumOfPolytopes(const boost::shared_ptr<Polytope_Rn>& A,
424  const boost::shared_ptr<Polytope_Rn>& B,
425  boost::shared_ptr<Polytope_Rn>& C,
426  const std::vector< std::vector<int> >& genitorsOfGeneratorsA,
427  const std::vector< std::vector<int> >& genitorsOfGeneratorsB,
428  std::vector< std::vector<int> >& traceGenerators) {
429  boost::timer this_timer;
430  try {
431  MinkowskiSum Ope(A,B,C,genitorsOfGeneratorsA,genitorsOfGeneratorsB,traceGenerators);
432  }
433  catch(invalid_argument& except) {
434  cerr << "TIME=" << this_timer.elapsed() << endl;
435  cerr << "Invalid argument exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
436  return TEST_KO;
437  }
438  catch(out_of_range& except) {
439  cerr << "TIME=" << this_timer.elapsed() << endl;
440  cerr << "Out of range exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
441  return TEST_KO;
442  }
443  catch(ios_base::failure& except) {
444  cerr << "TIME=" << this_timer.elapsed() << endl;
445  cerr << "In/out exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
446  return TEST_KO;
447  }
448  catch(logic_error& except) {
449  cerr << "TIME=" << this_timer.elapsed() << endl;
450  cerr << "Logic error exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
451  return TEST_KO;
452  }
453  catch(...) {
454  cerr << "TIME=" << this_timer.elapsed() << endl;
455  cerr << "Unexpected exception caught in politopixAPI::computeMinkowskiSumOfPolytopes() !" << endl;
456  return TEST_KO;
457  }
458  cout << "TIME=" << this_timer.elapsed() << endl;
459 
460  return C->checkTopologyAndGeometry() == true ? TEST_OK : TEST_KO;
461 }
462 
463 int politopixAPI::checkEqualityOfPolytopes(const boost::shared_ptr<Polytope_Rn>& A,
464  const boost::shared_ptr<Polytope_Rn>& B, bool getFaceMapping) {
465  bool ret_val;
466  boost::timer this_timer;
467  try {
468  ret_val = A->checkEquality(B, getFaceMapping);
469  }
470  catch(invalid_argument& except) {
471  cerr << "TIME=" << this_timer.elapsed() << endl;
472  cerr << "Invalid argument exception in politopixAPI::checkEqualityOfPolytopes() " << except.what() << endl;
473  return TEST_KO;
474  }
475  catch(out_of_range& except) {
476  cerr << "TIME=" << this_timer.elapsed() << endl;
477  cerr << "Out of range exception in politopixAPI::checkEqualityOfPolytopes() " << except.what() << endl;
478  return TEST_KO;
479  }
480  catch(ios_base::failure& except) {
481  cerr << "TIME=" << this_timer.elapsed() << endl;
482  cerr << "In/out exception in politopixAPI::checkEqualityOfPolytopes() " << except.what() << endl;
483  return TEST_KO;
484  }
485  catch(logic_error& except) {
486  cerr << "TIME=" << this_timer.elapsed() << endl;
487  cerr << "Logic error exception in politopixAPI::checkEqualityOfPolytopes() " << except.what() << endl;
488  return TEST_KO;
489  }
490  catch(...) {
491  cerr << "TIME=" << this_timer.elapsed() << endl;
492  cerr << "Unexpected exception caught in politopixAPI::checkEqualityOfPolytopes() " << endl;
493  return TEST_KO;
494  }
495  cout << "TIME=" << this_timer.elapsed() << endl;
496 
497  return ret_val == true ? TEST_OK : TEST_KO;
498 }
499 
500 bool politopixAPI::checkEqualityOfVertices(const boost::shared_ptr<Polytope_Rn>& A,
501  const boost::shared_ptr<Polytope_Rn>& B) {
502  bool ret_val;
503  boost::timer this_timer;
504  try {
505  ret_val = A->checkEqualityOfVertices(B);
506  }
507  catch(invalid_argument& except) {
508  cerr << "TIME=" << this_timer.elapsed() << endl;
509  cerr << "Invalid argument exception in politopixAPI::checkEqualityOfVertices() " << except.what() << endl;
510  return TEST_KO;
511  }
512  catch(out_of_range& except) {
513  cerr << "TIME=" << this_timer.elapsed() << endl;
514  cerr << "Out of range exception in politopixAPI::checkEqualityOfVertices() " << except.what() << endl;
515  return TEST_KO;
516  }
517  catch(ios_base::failure& except) {
518  cerr << "TIME=" << this_timer.elapsed() << endl;
519  cerr << "In/out exception in politopixAPI::checkEqualityOfVertices() " << except.what() << endl;
520  return TEST_KO;
521  }
522  catch(logic_error& except) {
523  cerr << "TIME=" << this_timer.elapsed() << endl;
524  cerr << "Logic error exception in politopixAPI::checkEqualityOfVertices() " << except.what() << endl;
525  return TEST_KO;
526  }
527  catch(...) {
528  cerr << "TIME=" << this_timer.elapsed() << endl;
529  cerr << "Unexpected exception caught in politopixAPI::checkEqualityOfVertices() " << endl;
530  return TEST_KO;
531  }
532  cout << "TIME=" << this_timer.elapsed() << endl;
533 
534  return ret_val == true ? TEST_OK : TEST_KO;
535 }
536 
537 double politopixAPI::computeVolume(const boost::shared_ptr<Polytope_Rn> Pol) {
538  double vol = -1.;
539  boost::timer this_timer;
540  try {
542  }
543  catch(std::invalid_argument& e) {
544  std::cerr << "politopixAPI::computeVolume() : invalid argument exception " << e.what() << std::endl;
545  return -1;
546  }
547  catch(std::out_of_range& e) {
548  std::cerr << "politopixAPI::computeVolume() : out of range exception " << e.what() << std::endl;
549  return -1;
550  }
551  catch(std::ios_base::failure& e) {
552  std::cerr << "politopixAPI::computeVolume() : in/out exception " << e.what() << std::endl;
553  return -1;
554  }
555  catch(...) {
556  std::cerr << "politopixAPI::computeVolume() : unexpected exception caught !" << std::endl;
557  return -1;
558  }
559  cout << "TIME=" << this_timer.elapsed() << endl;
560 
561  return vol;
562 }
563 
564 int politopixAPI::checkTopologyAndGeometry(const boost::shared_ptr<PolyhedralCone_Rn>& A, bool check_all) {
565  return A->checkTopologyAndGeometry(check_all) == true ? TEST_OK : TEST_KO;
566 }
567 
568 int politopixAPI::makeCube(boost::shared_ptr<Polytope_Rn>& A, double M) {
569  A.reset(new Polytope_Rn());
570  A->createBoundingBox(M);
571  return TEST_OK;
572 }
573 
574 int politopixAPI::makeBox(boost::shared_ptr<Polytope_Rn>& A, Point_Rn Pmin, Point_Rn Pmax) {
575  A.reset(new Polytope_Rn());
576  unsigned int dimension = Rn::getDimension();
577  double M = Pmin.getCoordinate(0);
578  {for (unsigned int i=0; i<dimension; ++i) {
579  if (M < fabs(Pmin.getCoordinate(i)))
580  M = fabs(Pmin.getCoordinate(i));
581  if (M < fabs(Pmax.getCoordinate(i)))
582  M = fabs(Pmax.getCoordinate(i));
583  }}
584  //A->createBoundingBox(M+1); // This will be done in computeDoubleDescriptionWithoutCheck()
585  // Now create the half-space of the shape: x_i >= Pmin.setCoordinate(i) <=> -Pmin.setCoordinate(i) + x_i >= 0.
586  {for (unsigned int i=0; i<dimension; ++i) {
587  boost::shared_ptr<HalfSpace_Rn> HS;
588  HS.reset(new HalfSpace_Rn(dimension));
589  HS->setConstant(-Pmin.getCoordinate(i));
590  for (unsigned int coord_count=0; coord_count<dimension; coord_count++) {
591  if (coord_count == i)
592  HS->setCoefficient(coord_count, 1.);
593  else
594  HS->setCoefficient(coord_count, 0.);
595  }
596  A->addHalfSpace(HS);
597  }}
598  // Now create the half-space of the shape: x_i <= Pmax.setCoordinate(i) <=> Pmax.setCoordinate(i) - x_i >= 0.
599  {for (unsigned int i=0; i<dimension; ++i) {
600  boost::shared_ptr<HalfSpace_Rn> HS;
601  HS.reset(new HalfSpace_Rn(dimension));
602  HS->setConstant(Pmax.getCoordinate(i));
603  for (unsigned int coord_count=0; coord_count<dimension; coord_count++) {
604  if (coord_count == i)
605  HS->setCoefficient(coord_count, -1.);
606  else
607  HS->setCoefficient(coord_count, 0.);
608  }
609  A->addHalfSpace(HS);
610  }}
612 }
613 
614 int politopixAPI::PolarPolytope(const boost::shared_ptr<Polytope_Rn>& original_pol, boost::shared_ptr<Polytope_Rn>& polar_pol) {
615  return TopGeomTools::PolarPolytope(original_pol, polar_pol);
616 }
617 
618 int politopixAPI::Translate(boost::shared_ptr<Polytope_Rn>& pol, const boost::numeric::ublas::vector<double>& v2t) {
619  return TopGeomTools::Translate(pol, v2t);
620 }
621 
623  const boost::shared_ptr<Polytope_Rn>& A,
624  const boost::shared_ptr<Polytope_Rn>& B,
625  boost::shared_ptr<Polytope_Rn>& C,
626  const std::set< unsigned int >& firstOperandCaps,
627  const std::set< unsigned int >& secondOperandCaps,
628  std::set< unsigned int >& newCaps,
629  double bb_size) {
630  boost::timer this_timer;
631  try {
632  PseudoIntersectionWithoutCaps PIWC(A, B, C, firstOperandCaps, secondOperandCaps, newCaps, bb_size);
633  }
634  catch(invalid_argument& except) {
635  cerr << "TIME=" << this_timer.elapsed() << endl;
636  cerr << "Invalid argument exception in politopixAPI::PseudoIntersection() " << except.what() << endl;
637  return TEST_KO;
638  }
639  catch(out_of_range& except) {
640  cerr << "TIME=" << this_timer.elapsed() << endl;
641  cerr << "Out of range exception in politopixAPI::PseudoIntersection() " << except.what() << endl;
642  return TEST_KO;
643  }
644  catch(ios_base::failure& except) {
645  cerr << "TIME=" << this_timer.elapsed() << endl;
646  cerr << "In/out exception in politopixAPI::PseudoIntersection() " << except.what() << endl;
647  return TEST_KO;
648  }
649  catch(logic_error& except) {
650  cerr << "TIME=" << this_timer.elapsed() << endl;
651  cerr << "Logic error exception in politopixAPI::PseudoIntersection() " << except.what() << endl;
652  return TEST_KO;
653  }
654  catch(...) {
655  cerr << "TIME=" << this_timer.elapsed() << endl;
656  cerr << "Unexpected exception caught in politopixAPI::PseudoIntersection() !" << endl;
657  return TEST_KO;
658  }
659  cout << "TIME=" << this_timer.elapsed() << endl;
660 
661  if (C->numberOfGenerators()==0 && C->numberOfHalfSpaces()==0) {
662  cout << "The result is empty." << std::endl;
663  return TEST_OK;
664  }
665 
666  return C->checkTopologyAndGeometry() == true ? TEST_OK : TEST_KO;
667 }
668 
670  const boost::shared_ptr<Polytope_Rn>& A,
671  const boost::shared_ptr<Polytope_Rn>& B,
672  boost::shared_ptr<Polytope_Rn>& C,
673  const std::set< unsigned int >& firstOperandCaps,
674  const std::set< unsigned int >& secondOperandCaps,
675  std::set< unsigned int >& newCaps,
676  double bb_size) {
677  boost::timer this_timer;
678  try {
679  PseudoSumWithoutCaps PSWC(A, B, C, firstOperandCaps, secondOperandCaps, newCaps, bb_size);
680  }
681  catch(invalid_argument& except) {
682  cerr << "TIME=" << this_timer.elapsed() << endl;
683  cerr << "Invalid argument exception in politopixAPI::pseudoSum() " << except.what() << endl;
684  return TEST_KO;
685  }
686  catch(out_of_range& except) {
687  cerr << "TIME=" << this_timer.elapsed() << endl;
688  cerr << "Out of range exception in politopixAPI::pseudoSum() " << except.what() << endl;
689  return TEST_KO;
690  }
691  catch(ios_base::failure& except) {
692  cerr << "TIME=" << this_timer.elapsed() << endl;
693  cerr << "In/out exception in politopixAPI::pseudoSum() " << except.what() << endl;
694  return TEST_KO;
695  }
696  catch(logic_error& except) {
697  cerr << "TIME=" << this_timer.elapsed() << endl;
698  cerr << "Logic error exception in politopixAPI::pseudoSum() " << except.what() << endl;
699  return TEST_KO;
700  }
701  catch(...) {
702  cerr << "TIME=" << this_timer.elapsed() << endl;
703  cerr << "Unexpected exception caught in politopixAPI::pseudoSum() !" << endl;
704  return TEST_KO;
705  }
706  cout << "TIME=" << this_timer.elapsed() << endl;
707 
708  if (C->numberOfGenerators()==0 && C->numberOfHalfSpaces()==0) {
709  cout << "The result is empty." << std::endl;
710  return TEST_OK;
711  }
712 
713  return C->checkTopologyAndGeometry() == true ? TEST_OK : TEST_KO;
714 }
715 
717  const boost::shared_ptr<Polytope_Rn>& inputSpace,
718  const std::vector<Point_Rn>& listOfSeeds,
719  std::vector< boost::shared_ptr<Polytope_Rn> >& VoronoiCells,
720  Voronoi_Rn::TypeOfAlgorithm vorAlgo) {
721  boost::timer this_timer;
722  boost::shared_ptr<Voronoi_Rn> VD;
723  if (vorAlgo == Voronoi_Rn::Incremental) {
724  VD.reset(new IncrementalVoronoi(inputSpace, listOfSeeds));
725  }
726  else if (vorAlgo == Voronoi_Rn::CellByCell) {
727  VD.reset(new CellByCellVoronoi(inputSpace, listOfSeeds));
728  }
729  else if (vorAlgo == Voronoi_Rn::CellByCellDist) {
730  VD.reset(new CellByCellVoronoiDist(inputSpace, listOfSeeds));
731  }
732  else if (vorAlgo == Voronoi_Rn::CellByCellDistNgb) {
733  VD.reset(new CellByCellVoronoiDistNgb(inputSpace, listOfSeeds));
734  }
735  else if (vorAlgo == Voronoi_Rn::CellByCellDistNgb_1pc) {
736  VD.reset(new CellByCellVoronoiDistNgb(inputSpace, listOfSeeds));
737  boost::shared_ptr<CellByCellVoronoiDistNgb> VDd = boost::static_pointer_cast<CellByCellVoronoiDistNgb>(VD);
738  VDd->setAverageNumberOfNeighbours(1.);
739  }
740  else if (vorAlgo == Voronoi_Rn::CellByCellDistNgb_10pc) {
741  VD.reset(new CellByCellVoronoiDistNgb(inputSpace, listOfSeeds));
742  boost::shared_ptr<CellByCellVoronoiDistNgb> VDd = boost::static_pointer_cast<CellByCellVoronoiDistNgb>(VD);
743  VDd->setAverageNumberOfNeighbours(10.);
744  }
745  else if (vorAlgo == Voronoi_Rn::CellByCellDistNgb_20pc) {
746  VD.reset(new CellByCellVoronoiDistNgb(inputSpace, listOfSeeds));
747  boost::shared_ptr<CellByCellVoronoiDistNgb> VDd = boost::static_pointer_cast<CellByCellVoronoiDistNgb>(VD);
748  VDd->setAverageNumberOfNeighbours(20.);
749  }
750  else if (vorAlgo == Voronoi_Rn::CellByCellDistNgb_30pc) {
751  VD.reset(new CellByCellVoronoiDistNgb(inputSpace, listOfSeeds));
752  boost::shared_ptr<CellByCellVoronoiDistNgb> VDd = boost::static_pointer_cast<CellByCellVoronoiDistNgb>(VD);
753  VDd->setAverageNumberOfNeighbours(30.);
754  }
755  else if (vorAlgo == Voronoi_Rn::CellByCellDistNgb_40pc) {
756  VD.reset(new CellByCellVoronoiDistNgb(inputSpace, listOfSeeds));
757  boost::shared_ptr<CellByCellVoronoiDistNgb> VDd = boost::static_pointer_cast<CellByCellVoronoiDistNgb>(VD);
758  VDd->setAverageNumberOfNeighbours(40.);
759  }
760  else
761  throw std::invalid_argument("politopixAPI::computeVoronoiDiagram() vorAlgo out of range");
762 
763  try {
764  VD->compute();
765  VoronoiCells = VD->getVoronoiCells();
766  }
767  catch(length_error& except) {
768  cerr << "TIME=" << this_timer.elapsed() << endl;
769  cerr << "Length error exception in politopixAPI::computeVoronoiDiagram() " << except.what() << endl;
770  return TEST_KO;
771  }
772  catch(...) {
773  cerr << "TIME=" << this_timer.elapsed() << endl;
774  cerr << "Unexpected exception caught in politopixAPI::computeVoronoiDiagram() !" << endl;
775  return TEST_KO;
776  }
777  cout << "TIME=" << this_timer.elapsed() << endl;
778  return VD->checkTopologyAndGeometry() == true ? TEST_OK : TEST_KO;
779 }
CellByCellVoronoiDist
Constructor.
Definition: Voronoi_Rn.h:185
politopixAPI::savePolytope
static polito_EXPORT int savePolytope(const string &pathA, boost::shared_ptr< Polytope_Rn > &A)
Save a polytope in the corresponding file name.
Definition: politopixAPI.cpp:29
PseudoSumWithoutCaps
Compute the Minkowski sum of two polytopes and then remove all cap half-spaces to truncate again.
Definition: PolyhedralAlgorithms_Rn.h:267
politopixAPI::makeBox
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.
Definition: politopixAPI.cpp:574
Point_Rn::getCoordinate
double getCoordinate(unsigned int i) const
Definition: Point_Rn.cpp:70
TEST_OK
#define TEST_OK
Definition: politopixAPI.h:38
politopixAPI.h
politopixAPI::checkEqualityOfVertices
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.
Definition: politopixAPI.cpp:500
Voronoi_Rn::CellByCellDist
@ CellByCellDist
Refers to the class CellByCellVoronoi.
Definition: Voronoi_Rn.h:47
politopixAPI::isIncluded
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.
Definition: politopixAPI.cpp:261
Voronoi_Rn::CellByCellDistNgb_1pc
@ CellByCellDistNgb_1pc
Refers to the class CellByCellDistNgb (100% of the distances between seeds to be sorted)
Definition: Voronoi_Rn.h:49
CellByCellVoronoi
Computes the Voronoi diagram cell by cell.
Definition: Voronoi_Rn.h:127
politopixAPI::computeMinkowskiSumOfPolytopes
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.
Definition: politopixAPI.cpp:350
IncrementalVoronoi
At the end of the ith iteration, the Voronoi diagram of i seeds is computed.
Definition: Voronoi_Rn.h:102
politopixAPI::pseudoIntersection
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.
Definition: politopixAPI.cpp:622
Voronoi_Rn::CellByCellDistNgb_40pc
@ CellByCellDistNgb_40pc
Refers to the class CellByCellDistNgb ( 30% of the distances between seeds to be sorted)
Definition: Voronoi_Rn.h:53
politopixAPI::computeIntersectionWithoutCheck
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.
Definition: politopixAPI.cpp:298
politopixAPI::computeVoronoiDiagram
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...
Definition: politopixAPI.cpp:716
Rn::getDimension
static polito_EXPORT unsigned int getDimension()
Return the dimension of the cartesian space we work in.
Definition: Rn.cpp:29
politopixAPI::computeIntersection
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.
Definition: politopixAPI.cpp:195
constIteratorOfListOfGeometricObjects::current
const GEOMETRIC_OBJECT current()
Return the current geometric element.
Definition: GeometricObjectIterator_Rn.h:177
Voronoi_Rn::CellByCellDistNgb_10pc
@ CellByCellDistNgb_10pc
Refers to the class CellByCellDistNgb ( 1% of the distances between seeds to be sorted)
Definition: Voronoi_Rn.h:50
constIteratorOfListOfGeometricObjects::end
bool end() const
Tell whether we have reached the end of the list.
Definition: GeometricObjectIterator_Rn.h:174
HalfSpace_Rn
A half-space whose frontier is a linear (n-1) dimension space. _constant + _coefficients[0]....
Definition: HalfSpace_Rn.h:40
politopixAPI::PolarPolytope
static polito_EXPORT int PolarPolytope(const boost::shared_ptr< Polytope_Rn > &original_pol, boost::shared_ptr< Polytope_Rn > &polar_pol)
Compute the polar polytope.
Definition: politopixAPI.cpp:614
Voronoi_Rn::CellByCellDistNgb
@ CellByCellDistNgb
Refers to the class CellByCellVoronoiDist.
Definition: Voronoi_Rn.h:48
politopixAPI::addHalfspace
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.
Definition: politopixAPI.cpp:75
TEST_KO
#define TEST_KO
Definition: politopixAPI.h:39
politopixAPI::addGenerator
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.
Definition: politopixAPI.cpp:80
Point_Rn
Creation of a n-coordinate geometric point designed to be shared by its neighbour faces.
Definition: Point_Rn.h:34
DoubleDescriptionFromGenerators::Compute
static int Compute(boost::shared_ptr< Polytope_Rn > &pol, double bb_size=1000.)
Use the polarity to get the facets from the generators.
Definition: PolyhedralAlgorithms_Rn.cpp:1634
constIteratorOfListOfGeometricObjects
This class is designed to run the list of all geometric objects representing a polytope.
Definition: GeometricObjectIterator_Rn.h:142
Voronoi_Rn::Incremental
@ Incremental
Definition: Voronoi_Rn.h:45
DoubleDescription
The algorithm implemented here is an incremental algorithm as mentioned in How Good are Convex Hull ...
Definition: DoubleDescription_Rn.h:86
politopixAPI::loadPolyhedralCone
static polito_EXPORT int loadPolyhedralCone(const string &pathA, boost::shared_ptr< PolyhedralCone_Rn > &A)
Load a polyhedral cone from the corresponding file name.
Definition: politopixAPI.cpp:64
politopixAPI::computeVolume
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 ...
Definition: politopixAPI.cpp:537
TopGeomTools::PolarPolytope
static int PolarPolytope(const boost::shared_ptr< Polytope_Rn > &original_pol, boost::shared_ptr< Polytope_Rn > &polar_pol, bool forceComputation=true, double bb_size=1000.)
Compute the polar polytope.
Definition: PolyhedralAlgorithms_Rn.cpp:1428
VolumeOfPolytopes_Rn::compute
static double compute(const boost::shared_ptr< Polytope_Rn > P)
Return the volume of the given polytope P.
Definition: VolumeOfPolytopes_Rn.h:153
constIteratorOfListOfGeometricObjects::begin
void begin()
Move the iterator at the beginning of the list.
Definition: GeometricObjectIterator_Rn.h:150
TopGeomTools::Translate
static int Translate(boost::shared_ptr< Polytope_Rn > &pol, const boost::numeric::ublas::vector< double > &v2t)
Translate a polytope by the given vector.
Definition: PolyhedralAlgorithms_Rn.cpp:1248
politopixAPI::computeDoubleDescription
static polito_EXPORT int computeDoubleDescription(boost::shared_ptr< Polytope_Rn > &A, double bb_size)
Compute the HV-description for a given H-polytope or V-polytope with the Double Description algorithm...
Definition: politopixAPI.cpp:136
PseudoIntersectionWithoutCaps
Remove all cap half-spaces and then compute the intersection of two capped polytopes.
Definition: PolyhedralAlgorithms_Rn.h:304
politopixAPI::checkTopologyAndGeometry
static polito_EXPORT int checkTopologyAndGeometry(const boost::shared_ptr< PolyhedralCone_Rn > &A, bool check_all=false)
Check whether a HV-polytopes is correct.
Definition: politopixAPI.cpp:564
politopy.h
politopixAPI::pseudoSum
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.
Definition: politopixAPI.cpp:669
politopixAPI::makeCube
static polito_EXPORT int makeCube(boost::shared_ptr< Polytope_Rn > &A, double M)
Create a cube whose vertices will be (+-M, ..., +-M)
Definition: politopixAPI.cpp:568
politopixAPI::loadPolytope
static polito_EXPORT int loadPolytope(const string &pathA, boost::shared_ptr< Polytope_Rn > &A)
Load a polytope from the corresponding file name.
Definition: politopixAPI.cpp:40
IO_Polytope::load
static polito_EXPORT void load(const std::string &filename, boost::shared_ptr< PolyhedralCone_Rn > POLY)
Load the main data format to store polytopes.
Definition: IO_Polytope.cpp:26
Voronoi_Rn::TypeOfAlgorithm
TypeOfAlgorithm
Choose the kind of algorithm to be run.
Definition: Voronoi_Rn.h:44
Voronoi_Rn::CellByCellDistNgb_30pc
@ CellByCellDistNgb_30pc
Refers to the class CellByCellDistNgb ( 20% of the distances between seeds to be sorted)
Definition: Voronoi_Rn.h:52
constIteratorOfListOfGeometricObjects::next
void next()
Move the iterator one step forward.
Definition: GeometricObjectIterator_Rn.h:153
CellByCellVoronoiDistNgb
Do the same as CellByCellVoronoiDist with the neighbourhood computation.
Definition: Voronoi_Rn.h:231
politopixAPI::Translate
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.
Definition: politopixAPI.cpp:618
IO_Polytope::save
static polito_EXPORT void save(const std::string &filename, boost::shared_ptr< PolyhedralCone_Rn > POLY)
Save the polytope to the main data format.
Definition: IO_Polytope.cpp:132
Voronoi_Rn::CellByCell
@ CellByCell
Refers to the class IncrementalVoronoi.
Definition: Voronoi_Rn.h:46
politopixAPI::computeDoubleDescriptionWithoutCheck
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...
Definition: politopixAPI.cpp:85
Voronoi_Rn::CellByCellDistNgb_20pc
@ CellByCellDistNgb_20pc
Refers to the class CellByCellDistNgb ( 10% of the distances between seeds to be sorted)
Definition: Voronoi_Rn.h:51
Polytope_Rn
Model a polytope using its two equivalent definitions : the convex hull and the half-space intersecti...
Definition: Polytope_Rn.h:36
politopixAPI::checkEqualityOfPolytopes
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...
Definition: politopixAPI.cpp:463
politopixAPI::savePolyhedralCone
static polito_EXPORT int savePolyhedralCone(const string &pathA, boost::shared_ptr< PolyhedralCone_Rn > &A)
Save a polyhedral cone in the corresponding file name.
Definition: politopixAPI.cpp:52
MinkowskiSum
Compute the Minkowski sum of two polytopes.
Definition: PolyhedralAlgorithms_Rn.h:138