politopix  4.1.0
 All Classes Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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-2017 : Delos Vincent
3 //
4 // This program is free software: you can redistribute it and/or modify
5 // it under the terms of the GNU Lesser 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 Lesser General Public License for more details.
13 //
14 // You should have received a copy of the GNU Lesser 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 
25 
26 int politopixAPI::savePolytope(const string& pathA, boost::shared_ptr<Polytope_Rn>& A) throw (ios_base::failure)
27 {
28  try {
29  IO_Polytope::save(pathA, A);
30  }
31  catch(ios_base::failure& except) {
32  cerr << "In/out exception in politopixAPI::savePolytope(const string&, boost::shared_ptr<Polytope_Rn>&): " << except.what() << endl;
33  throw except;
34  }
35  return TEST_OK;
36 }
37 
38 int politopixAPI::loadPolytope(const string& pathA, boost::shared_ptr<Polytope_Rn>& A) throw (ios_base::failure)
39 {
40  try {
41  IO_Polytope::load(pathA, A);
42  }
43  catch(ios_base::failure& except) {
44  cerr << "In/out exception in politopixAPI::loadPolytope(const string&, boost::shared_ptr<Polytope_Rn>&): " << except.what() << endl;
45  throw except;
46  }
47  return TEST_OK;
48 }
49 
50 int politopixAPI::savePolyhedralCone(const string& pathA, boost::shared_ptr<PolyhedralCone_Rn>& A) throw (ios_base::failure)
51 {
52  try {
53  IO_Polytope::save(pathA, A);
54  }
55  catch(ios_base::failure& except) {
56  cerr << "In/out exception in politopixAPI::savePolyhedralCone(const string&, boost::shared_ptr<PolyhedralCone_Rn>&): " << except.what() << endl;
57  throw except;
58  }
59  return TEST_OK;
60 }
61 
62 int politopixAPI::loadPolyhedralCone(const string& pathA, boost::shared_ptr<PolyhedralCone_Rn>& A) throw (ios_base::failure)
63 {
64  try {
65  IO_Polytope::load(pathA, A);
66  }
67  catch(ios_base::failure& except) {
68  cerr << "In/out exception in politopixAPI::loadPolyhedralCone(const string&, boost::shared_ptr<PolyhedralCone_Rn>&): " << except.what() << endl;
69  throw except;
70  }
71  return TEST_OK;
72 }
73 
74 int politopixAPI::addHalfspace(boost::shared_ptr<Polytope_Rn>& A, const boost::shared_ptr<HalfSpace_Rn>& HS) {
75  A->addHalfSpace(HS);
76  return TEST_OK;
77 }
78 
79 int politopixAPI::addGenerator(boost::shared_ptr<Polytope_Rn>& A, const boost::shared_ptr<Generator_Rn>& GN) {
80  A->addGenerator(GN);
81  return TEST_OK;
82 }
83 
84 int politopixAPI::computeDoubleDescriptionWithoutCheck(boost::shared_ptr<Polytope_Rn>& A, double bb_size) throw (invalid_argument,out_of_range,ios_base::failure,logic_error)
85 {
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);
94  // The number of facets for a cube.
95  unsigned int truncationStep = 2*A->dimension();
97  boost::shared_ptr<PolyhedralCone_Rn>,
100  DD(A, lexmin_ite, NRP, 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) throw (invalid_argument,out_of_range,ios_base::failure,logic_error)
137 {
138  int res = computeDoubleDescriptionWithoutCheck(A, bb_size);
139  if (res == TEST_KO)
140  return TEST_KO;
141  return A->checkTopologyAndGeometry() == true ? TEST_OK : TEST_KO;
142 }
143 
144 int politopixAPI::computeIntersection(boost::shared_ptr<Polytope_Rn>& A,
145  const boost::shared_ptr<Polytope_Rn>& B) throw (invalid_argument,out_of_range,ios_base::failure,logic_error)
146 {
148  unsigned int truncationStep = A->numberOfHalfSpaces();
149  for (iteHSB.begin(); iteHSB.end()!=true; iteHSB.next()) {
150  A->addHalfSpace(iteHSB.current());
151  }
152  boost::timer this_timer;
153  try {
154  constIteratorOfListOfGeometricObjects< boost::shared_ptr<HalfSpace_Rn> > lexmin_ite(A->getListOfHalfSpaces());
157  boost::shared_ptr<PolyhedralCone_Rn>,
160  DD(A, lexmin_ite, NRP, truncationStep);
161  }
162  catch(invalid_argument& except) {
163  cerr << "TIME=" << this_timer.elapsed() << endl;
164  cerr << "Invalid argument exception in politopixAPI::computeIntersection() " << except.what() << endl;
165  return TEST_KO;
166  }
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;
170  return TEST_KO;
171  }
172  catch(ios_base::failure& except) {
173  cerr << "TIME=" << this_timer.elapsed() << endl;
174  cerr << "In/out exception in politopixAPI::computeIntersection() " << except.what() << endl;
175  return TEST_KO;
176  }
177  catch(logic_error& except) {
178  cerr << "TIME=" << this_timer.elapsed() << endl;
179  cerr << "Logic error exception in politopixAPI::computeIntersection() " << except.what() << endl;
180  return TEST_KO;
181  }
182  catch(...) {
183  cerr << "TIME=" << this_timer.elapsed() << endl;
184  cerr << "Unexpected exception caught in politopixAPI::computeIntersection() !" << endl;
185  return TEST_KO;
186  }
187  cout << "TIME=" << this_timer.elapsed() << endl;
188 
189  if (A->numberOfGenerators()==0 && A->numberOfHalfSpaces()==0) {
190  cout << "The result is empty." << std::endl;
191  return TEST_OK;
192  }
193 
194  return A->checkTopologyAndGeometry() == true ? TEST_OK : TEST_KO;
195 }
196 
197 int politopixAPI::computeIntersection(const boost::shared_ptr<Polytope_Rn>& A, const boost::shared_ptr<Polytope_Rn>& B,
198  boost::shared_ptr<Polytope_Rn>& C) throw (invalid_argument,out_of_range,ios_base::failure,logic_error)
199 {
200  // Copy the generators of A.
202  for (iteHS0.begin(); iteHS0.end()!=true; iteHS0.next()) {
203  C->addGenerator(iteHS0.current());
204  }
205  // Copy the half-spaces of A.
207  for (iteHS1.begin(); iteHS1.end()!=true; iteHS1.next()) {
208  C->addHalfSpace(iteHS1.current());
209  }
210  // Copy the half-spaces of B.
212  for (iteHS2.begin(); iteHS2.end()!=true; iteHS2.next()) {
213  C->addHalfSpace(iteHS2.current());
214  }
215 
216  boost::timer this_timer;
217  try {
218  // Now the truncation.
219  unsigned int truncationStep = A->numberOfHalfSpaces();
220  constIteratorOfListOfGeometricObjects< boost::shared_ptr<HalfSpace_Rn> > lexmin_ite(C->getListOfHalfSpaces());
222  // The number of facets for a cube.
224  boost::shared_ptr<PolyhedralCone_Rn>,
227  DD(C, lexmin_ite, NRP, truncationStep);
228  }
229  catch(invalid_argument& except) {
230  cerr << "TIME=" << this_timer.elapsed() << endl;
231  cerr << "Invalid argument exception in politopixAPI::computeIntersection() " << except.what() << endl;
232  return TEST_KO;
233  }
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;
237  return TEST_KO;
238  }
239  catch(ios_base::failure& except) {
240  cerr << "TIME=" << this_timer.elapsed() << endl;
241  cerr << "In/out exception in politopixAPI::computeIntersection() " << except.what() << endl;
242  return TEST_KO;
243  }
244  catch(logic_error& except) {
245  cerr << "TIME=" << this_timer.elapsed() << endl;
246  cerr << "Logic error exception in politopixAPI::computeIntersection() " << except.what() << endl;
247  return TEST_KO;
248  }
249  catch(...) {
250  cerr << "TIME=" << this_timer.elapsed() << endl;
251  cerr << "Unexpected exception caught in politopixAPI::computeIntersection() !" << endl;
252  return TEST_KO;
253  }
254  cout << "TIME=" << this_timer.elapsed() << endl;
255 
256  if (C->numberOfGenerators()==0 && C->numberOfHalfSpaces()==0) {
257  cout << "The result is empty." << std::endl;
258  return TEST_OK;
259  }
260 
261  return C->checkTopologyAndGeometry() == true ? TEST_OK : TEST_KO;
262 }
263 
264 bool politopixAPI::isIncluded(const boost::shared_ptr<Polytope_Rn>& A,
265  const boost::shared_ptr<Polytope_Rn>& B) throw (invalid_argument,out_of_range,ios_base::failure,logic_error)
266 {
267  boost::timer this_timer;
268  bool isInside = false;
269  try {
270  isInside = A->isIncluded(B);
271  }
272  catch(invalid_argument& except) {
273  cerr << "TIME=" << this_timer.elapsed() << endl;
274  cerr << "Invalid argument exception in politopixAPI::computeIntersection() " << except.what() << endl;
275  return TEST_KO;
276  }
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;
280  return TEST_KO;
281  }
282  catch(ios_base::failure& except) {
283  cerr << "TIME=" << this_timer.elapsed() << endl;
284  cerr << "In/out exception in politopixAPI::computeIntersection() " << except.what() << endl;
285  return TEST_KO;
286  }
287  catch(logic_error& except) {
288  cerr << "TIME=" << this_timer.elapsed() << endl;
289  cerr << "Logic error exception in politopixAPI::computeIntersection() " << except.what() << endl;
290  return TEST_KO;
291  }
292  catch(...) {
293  cerr << "TIME=" << this_timer.elapsed() << endl;
294  cerr << "Unexpected exception caught in politopixAPI::computeIntersection() !" << endl;
295  return TEST_KO;
296  }
297  cout << "TIME=" << this_timer.elapsed() << endl;
298 
299  return isInside;
300 }
301 
302 int politopixAPI::computeIntersectionWithoutCheck(boost::shared_ptr<Polytope_Rn>& A,
303  const boost::shared_ptr<Polytope_Rn>& B) throw (invalid_argument,out_of_range,ios_base::failure,logic_error)
304 {
306  unsigned int truncationStep = A->numberOfHalfSpaces();
307  for (iteHSB.begin(); iteHSB.end()!=true; iteHSB.next()) {
308  A->addHalfSpace(iteHSB.current());
309  }
310  boost::timer this_timer;
311  try {
312  constIteratorOfListOfGeometricObjects< boost::shared_ptr<HalfSpace_Rn> > lexmin_ite(A->getListOfHalfSpaces());
315  boost::shared_ptr<PolyhedralCone_Rn>,
318  DD(A, lexmin_ite, NRP, truncationStep);
319  }
320  catch(invalid_argument& except) {
321  cerr << "TIME=" << this_timer.elapsed() << endl;
322  cerr << "Invalid argument exception in politopixAPI::computeIntersectionWithoutCheck() " << except.what() << endl;
323  return TEST_KO;
324  }
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;
328  return TEST_KO;
329  }
330  catch(ios_base::failure& except) {
331  cerr << "TIME=" << this_timer.elapsed() << endl;
332  cerr << "In/out exception in politopixAPI::computeIntersectionWithoutCheck() " << except.what() << endl;
333  return TEST_KO;
334  }
335  catch(logic_error& except) {
336  cerr << "TIME=" << this_timer.elapsed() << endl;
337  cerr << "Logic error exception in politopixAPI::computeIntersectionWithoutCheck() " << except.what() << endl;
338  return TEST_KO;
339  }
340  catch(...) {
341  cerr << "TIME=" << this_timer.elapsed() << endl;
342  cerr << "Unexpected exception caught in politopixAPI::computeIntersectionWithoutCheck() !" << endl;
343  return TEST_KO;
344  }
345  //cout << "TIME=" << this_timer.elapsed() << endl;
346 
347  if (A->numberOfGenerators()==0 && A->numberOfHalfSpaces()==0) {
348  //cout << "The result is empty." << std::endl;
349  return TEST_OK;
350  }
351 
352  return TEST_OK;
353 }
354 
355 int politopixAPI::computeMinkowskiSumOfPolytopes(const boost::shared_ptr<Polytope_Rn>& A,
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)
358 {
359  boost::timer this_timer;
360  try {
361  MinkowskiSum Ope(A,B,C);
362  }
363  catch(invalid_argument& except) {
364  cerr << "TIME=" << this_timer.elapsed() << endl;
365  cerr << "Invalid argument exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
366  return TEST_KO;
367  }
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;
371  return TEST_KO;
372  }
373  catch(ios_base::failure& except) {
374  cerr << "TIME=" << this_timer.elapsed() << endl;
375  cerr << "In/out exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
376  return TEST_KO;
377  }
378  catch(logic_error& except) {
379  cerr << "TIME=" << this_timer.elapsed() << endl;
380  cerr << "Logic error exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
381  return TEST_KO;
382  }
383  catch(...) {
384  cerr << "TIME=" << this_timer.elapsed() << endl;
385  cerr << "Unexpected exception caught in politopixAPI::computeMinkowskiSumOfPolytopes() !" << endl;
386  return TEST_KO;
387  }
388  cout << "TIME=" << this_timer.elapsed() << endl;
389 
390  return C->checkTopologyAndGeometry() == true ? TEST_OK : TEST_KO;
391 }
392 
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)
396 {
397  boost::timer this_timer;
398  try {
399  MinkowskiSum Ope(arrayOfPol,C);
400  }
401  catch(invalid_argument& except) {
402  cerr << "TIME=" << this_timer.elapsed() << endl;
403  cerr << "Invalid argument exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
404  return TEST_KO;
405  }
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;
409  return TEST_KO;
410  }
411  catch(ios_base::failure& except) {
412  cerr << "TIME=" << this_timer.elapsed() << endl;
413  cerr << "In/out exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
414  return TEST_KO;
415  }
416  catch(logic_error& except) {
417  cerr << "TIME=" << this_timer.elapsed() << endl;
418  cerr << "Logic error exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
419  return TEST_KO;
420  }
421  catch(...) {
422  cerr << "TIME=" << this_timer.elapsed() << endl;
423  cerr << "Unexpected exception caught in politopixAPI::computeMinkowskiSumOfPolytopes() !" << endl;
424  return TEST_KO;
425  }
426  cout << "TIME=" << this_timer.elapsed() << endl;
427 
428  return C->checkTopologyAndGeometry() == true ? TEST_OK : TEST_KO;
429 }
430 
431 int politopixAPI::computeMinkowskiSumOfPolytopes(const boost::shared_ptr<Polytope_Rn>& A,
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)
438 {
439  boost::timer this_timer;
440  try {
441  MinkowskiSum Ope(A,B,C,genitorsOfGeneratorsA,genitorsOfGeneratorsB,traceGenerators);
442  }
443  catch(invalid_argument& except) {
444  cerr << "TIME=" << this_timer.elapsed() << endl;
445  cerr << "Invalid argument exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
446  return TEST_KO;
447  }
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;
451  return TEST_KO;
452  }
453  catch(ios_base::failure& except) {
454  cerr << "TIME=" << this_timer.elapsed() << endl;
455  cerr << "In/out exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
456  return TEST_KO;
457  }
458  catch(logic_error& except) {
459  cerr << "TIME=" << this_timer.elapsed() << endl;
460  cerr << "Logic error exception in politopixAPI::computeMinkowskiSumOfPolytopes() " << except.what() << endl;
461  return TEST_KO;
462  }
463  catch(...) {
464  cerr << "TIME=" << this_timer.elapsed() << endl;
465  cerr << "Unexpected exception caught in politopixAPI::computeMinkowskiSumOfPolytopes() !" << endl;
466  return TEST_KO;
467  }
468  cout << "TIME=" << this_timer.elapsed() << endl;
469 
470  return C->checkTopologyAndGeometry() == true ? TEST_OK : TEST_KO;
471 }
472 
473 int politopixAPI::checkEqualityOfPolytopes(const boost::shared_ptr<Polytope_Rn>& A,
474  const boost::shared_ptr<Polytope_Rn>& B, bool getFaceMapping) throw (invalid_argument,out_of_range,ios_base::failure,logic_error)
475 {
476  bool ret_val;
477  boost::timer this_timer;
478  try {
479  ret_val = A->checkEquality(B, getFaceMapping);
480  }
481  catch(invalid_argument& except) {
482  cerr << "TIME=" << this_timer.elapsed() << endl;
483  cerr << "Invalid argument exception in politopixAPI::checkEqualityOfPolytopes() " << except.what() << endl;
484  return TEST_KO;
485  }
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;
489  return TEST_KO;
490  }
491  catch(ios_base::failure& except) {
492  cerr << "TIME=" << this_timer.elapsed() << endl;
493  cerr << "In/out exception in politopixAPI::checkEqualityOfPolytopes() " << except.what() << endl;
494  return TEST_KO;
495  }
496  catch(logic_error& except) {
497  cerr << "TIME=" << this_timer.elapsed() << endl;
498  cerr << "Logic error exception in politopixAPI::checkEqualityOfPolytopes() " << except.what() << endl;
499  return TEST_KO;
500  }
501  catch(...) {
502  cerr << "TIME=" << this_timer.elapsed() << endl;
503  cerr << "Unexpected exception caught in politopixAPI::checkEqualityOfPolytopes() " << endl;
504  return TEST_KO;
505  }
506  cout << "TIME=" << this_timer.elapsed() << endl;
507 
508  return ret_val == true ? TEST_OK : TEST_KO;
509 }
510 
511 bool politopixAPI::checkEqualityOfVertices(const boost::shared_ptr<Polytope_Rn>& A,
512  const boost::shared_ptr<Polytope_Rn>& B) throw (invalid_argument,out_of_range,ios_base::failure,logic_error)
513 {
514  bool ret_val;
515  boost::timer this_timer;
516  try {
517  ret_val = A->checkEqualityOfVertices(B);
518  }
519  catch(invalid_argument& except) {
520  cerr << "TIME=" << this_timer.elapsed() << endl;
521  cerr << "Invalid argument exception in politopixAPI::checkEqualityOfVertices() " << except.what() << endl;
522  return TEST_KO;
523  }
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;
527  return TEST_KO;
528  }
529  catch(ios_base::failure& except) {
530  cerr << "TIME=" << this_timer.elapsed() << endl;
531  cerr << "In/out exception in politopixAPI::checkEqualityOfVertices() " << except.what() << endl;
532  return TEST_KO;
533  }
534  catch(logic_error& except) {
535  cerr << "TIME=" << this_timer.elapsed() << endl;
536  cerr << "Logic error exception in politopixAPI::checkEqualityOfVertices() " << except.what() << endl;
537  return TEST_KO;
538  }
539  catch(...) {
540  cerr << "TIME=" << this_timer.elapsed() << endl;
541  cerr << "Unexpected exception caught in politopixAPI::checkEqualityOfVertices() " << endl;
542  return TEST_KO;
543  }
544  cout << "TIME=" << this_timer.elapsed() << endl;
545 
546  return ret_val == true ? TEST_OK : TEST_KO;
547 }
548 
549 double politopixAPI::computeVolume(const boost::shared_ptr<Polytope_Rn> Pol) throw (invalid_argument, out_of_range, ios_base::failure)
550 {
551  double vol = -1.;
552  boost::timer this_timer;
553  try {
555  }
556  catch(std::invalid_argument& e) {
557  std::cerr << "politopixAPI::computeVolume() : invalid argument exception " << e.what() << std::endl;
558  return -1;
559  }
560  catch(std::out_of_range& e) {
561  std::cerr << "politopixAPI::computeVolume() : out of range exception " << e.what() << std::endl;
562  return -1;
563  }
564  catch(std::ios_base::failure& e) {
565  std::cerr << "politopixAPI::computeVolume() : in/out exception " << e.what() << std::endl;
566  return -1;
567  }
568  catch(...) {
569  std::cerr << "politopixAPI::computeVolume() : unexpected exception caught !" << std::endl;
570  return -1;
571  }
572  cout << "TIME=" << this_timer.elapsed() << endl;
573 
574  return vol;
575 }
576 
577 int politopixAPI::checkTopologyAndGeometry(const boost::shared_ptr<PolyhedralCone_Rn>& A) {
578  return A->checkTopologyAndGeometry() == true ? TEST_OK : TEST_KO;
579 }
580 
581 int politopixAPI::makeCube(boost::shared_ptr<Polytope_Rn>& A, double M) {
582  A.reset(new Polytope_Rn());
583  A->createBoundingBox(M);
584  return TEST_OK;
585 }
586 
587 int politopixAPI::makeBox(boost::shared_ptr<Polytope_Rn>& A, Point_Rn Pmin, Point_Rn Pmax) {
588  A.reset(new Polytope_Rn());
589  unsigned int dimension = Rn::getDimension();
590  double M = Pmin.getCoordinate(0);
591  {for (unsigned int i=0; i<dimension; ++i) {
592  if (M < fabs(Pmin.getCoordinate(i)))
593  M = fabs(Pmin.getCoordinate(i));
594  if (M < fabs(Pmax.getCoordinate(i)))
595  M = fabs(Pmax.getCoordinate(i));
596  }}
597  //A->createBoundingBox(M+1); // This will be done in computeDoubleDescriptionWithoutCheck()
598  // Now create the half-space of the shape: x_i >= Pmin.setCoordinate(i) <=> -Pmin.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(-Pmin.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  }}
611  // Now create the half-space of the shape: x_i <= Pmax.setCoordinate(i) <=> Pmax.setCoordinate(i) - x_i >= 0.
612  {for (unsigned int i=0; i<dimension; ++i) {
613  boost::shared_ptr<HalfSpace_Rn> HS;
614  HS.reset(new HalfSpace_Rn(dimension));
615  HS->setConstant(Pmax.getCoordinate(i));
616  for (unsigned int coord_count=0; coord_count<dimension; coord_count++) {
617  if (coord_count == i)
618  HS->setCoefficient(coord_count, -1.);
619  else
620  HS->setCoefficient(coord_count, 0.);
621  }
622  A->addHalfSpace(HS);
623  }}
625 }
626 
627 int politopixAPI::PolarPolytope(const boost::shared_ptr<Polytope_Rn>& original_pol, boost::shared_ptr<Polytope_Rn>& polar_pol) {
628  return TopGeomTools::PolarPolytope(original_pol, polar_pol);
629 }
630 
631 int politopixAPI::Translate(boost::shared_ptr<Polytope_Rn>& pol, const boost::numeric::ublas::vector<double>& v2t) {
632  return TopGeomTools::Translate(pol, v2t);
633 }
634 
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)
643 {
644  boost::timer this_timer;
645  try {
646  PseudoIntersectionWithoutCaps PIWC(A, B, C, firstOperandCaps, secondOperandCaps, newCaps, bb_size);
647  }
648  catch(invalid_argument& except) {
649  cerr << "TIME=" << this_timer.elapsed() << endl;
650  cerr << "Invalid argument exception in politopixAPI::PseudoIntersection() " << except.what() << endl;
651  return TEST_KO;
652  }
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;
656  return TEST_KO;
657  }
658  catch(ios_base::failure& except) {
659  cerr << "TIME=" << this_timer.elapsed() << endl;
660  cerr << "In/out exception in politopixAPI::PseudoIntersection() " << except.what() << endl;
661  return TEST_KO;
662  }
663  catch(logic_error& except) {
664  cerr << "TIME=" << this_timer.elapsed() << endl;
665  cerr << "Logic error exception in politopixAPI::PseudoIntersection() " << except.what() << endl;
666  return TEST_KO;
667  }
668  catch(...) {
669  cerr << "TIME=" << this_timer.elapsed() << endl;
670  cerr << "Unexpected exception caught in politopixAPI::PseudoIntersection() !" << endl;
671  return TEST_KO;
672  }
673  cout << "TIME=" << this_timer.elapsed() << endl;
674 
675  if (C->numberOfGenerators()==0 && C->numberOfHalfSpaces()==0) {
676  cout << "The result is empty." << std::endl;
677  return TEST_OK;
678  }
679 
680  return C->checkTopologyAndGeometry() == true ? TEST_OK : TEST_KO;
681 }
682 
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)
691 {
692  boost::timer this_timer;
693  try {
694  PseudoSumWithoutCaps PSWC(A, B, C, firstOperandCaps, secondOperandCaps, newCaps, bb_size);
695  }
696  catch(invalid_argument& except) {
697  cerr << "TIME=" << this_timer.elapsed() << endl;
698  cerr << "Invalid argument exception in politopixAPI::pseudoSum() " << except.what() << endl;
699  return TEST_KO;
700  }
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;
704  return TEST_KO;
705  }
706  catch(ios_base::failure& except) {
707  cerr << "TIME=" << this_timer.elapsed() << endl;
708  cerr << "In/out exception in politopixAPI::pseudoSum() " << except.what() << endl;
709  return TEST_KO;
710  }
711  catch(logic_error& except) {
712  cerr << "TIME=" << this_timer.elapsed() << endl;
713  cerr << "Logic error exception in politopixAPI::pseudoSum() " << except.what() << endl;
714  return TEST_KO;
715  }
716  catch(...) {
717  cerr << "TIME=" << this_timer.elapsed() << endl;
718  cerr << "Unexpected exception caught in politopixAPI::pseudoSum() !" << endl;
719  return TEST_KO;
720  }
721  cout << "TIME=" << this_timer.elapsed() << endl;
722 
723  if (C->numberOfGenerators()==0 && C->numberOfHalfSpaces()==0) {
724  cout << "The result is empty." << std::endl;
725  return TEST_OK;
726  }
727 
728  return C->checkTopologyAndGeometry() == true ? TEST_OK : TEST_KO;
729 }
730 
732  const boost::shared_ptr<Polytope_Rn>& inputSpace,
733  const std::vector<Point_Rn>& listOfSeeds,
734  std::vector< boost::shared_ptr<Polytope_Rn> >& VoronoiCells,
735  Voronoi_Rn::TypeOfAlgorithm vorAlgo) throw (std::invalid_argument)
736 {
737  boost::timer this_timer;
738  boost::shared_ptr<Voronoi_Rn> VD;
739  if (vorAlgo == Voronoi_Rn::Incremental) {
740  VD.reset(new IncrementalVoronoi(inputSpace, listOfSeeds));
741  }
742  else if (vorAlgo == Voronoi_Rn::CellByCell) {
743  VD.reset(new CellByCellVoronoi(inputSpace, listOfSeeds));
744  }
745  else if (vorAlgo == Voronoi_Rn::CellByCellDist) {
746  VD.reset(new CellByCellVoronoiDist(inputSpace, listOfSeeds));
747  }
748  else if (vorAlgo == Voronoi_Rn::CellByCellDistNgb) {
749  VD.reset(new CellByCellVoronoiDistNgb(inputSpace, listOfSeeds));
750  }
751  else if (vorAlgo == Voronoi_Rn::CellByCellDistNgb_1pc) {
752  VD.reset(new CellByCellVoronoiDistNgb(inputSpace, listOfSeeds));
753  boost::shared_ptr<CellByCellVoronoiDistNgb> VDd = boost::static_pointer_cast<CellByCellVoronoiDistNgb>(VD);
755  }
756  else if (vorAlgo == Voronoi_Rn::CellByCellDistNgb_10pc) {
757  VD.reset(new CellByCellVoronoiDistNgb(inputSpace, listOfSeeds));
758  boost::shared_ptr<CellByCellVoronoiDistNgb> VDd = boost::static_pointer_cast<CellByCellVoronoiDistNgb>(VD);
760  }
761  else if (vorAlgo == Voronoi_Rn::CellByCellDistNgb_20pc) {
762  VD.reset(new CellByCellVoronoiDistNgb(inputSpace, listOfSeeds));
763  boost::shared_ptr<CellByCellVoronoiDistNgb> VDd = boost::static_pointer_cast<CellByCellVoronoiDistNgb>(VD);
765  }
766  else if (vorAlgo == Voronoi_Rn::CellByCellDistNgb_30pc) {
767  VD.reset(new CellByCellVoronoiDistNgb(inputSpace, listOfSeeds));
768  boost::shared_ptr<CellByCellVoronoiDistNgb> VDd = boost::static_pointer_cast<CellByCellVoronoiDistNgb>(VD);
770  }
771  else if (vorAlgo == Voronoi_Rn::CellByCellDistNgb_40pc) {
772  VD.reset(new CellByCellVoronoiDistNgb(inputSpace, listOfSeeds));
773  boost::shared_ptr<CellByCellVoronoiDistNgb> VDd = boost::static_pointer_cast<CellByCellVoronoiDistNgb>(VD);
775  }
776  else
777  throw std::invalid_argument("politopixAPI::computeVoronoiDiagram() vorAlgo out of range");
778 
779  try {
780  VD->compute();
781  VoronoiCells = VD->getVoronoiCells();
782  }
783  catch(length_error& except) {
784  cerr << "TIME=" << this_timer.elapsed() << endl;
785  cerr << "Length error exception in politopixAPI::computeVoronoiDiagram() " << except.what() << endl;
786  return TEST_KO;
787  }
788  catch(...) {
789  cerr << "TIME=" << this_timer.elapsed() << endl;
790  cerr << "Unexpected exception caught in politopixAPI::computeVoronoiDiagram() !" << endl;
791  return TEST_KO;
792  }
793  cout << "TIME=" << this_timer.elapsed() << endl;
794  return VD->checkTopologyAndGeometry() == true ? TEST_OK : TEST_KO;
795 }
796 
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...
Definition: Point_Rn.h:34
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.
Computes the Voronoi diagram cell by cell.
Definition: Voronoi_Rn.h:127
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...
Constructor.
Definition: Voronoi_Rn.h:185
Model a polytope using its two equivalent definitions : the convex hull and the half-space intersecti...
Definition: Polytope_Rn.h:34
A half-space whose frontier is a linear (n-1) dimension space. _constant + _coefficients[0].x1 + ... + _coefficients[n-1].xn >= 0.
Definition: HalfSpace_Rn.h:37
Refers to the class CellByCellDistNgb (100% of the distances between seeds to be sorted) ...
Definition: Voronoi_Rn.h:49
double getCoordinate(unsigned int i) const
Definition: Point_Rn.cpp:70
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.
Definition: Voronoi_Rn.h:230
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.
Definition: IO_Polytope.cpp:26
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.
Definition: Voronoi_Rn.h:47
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) ...
Definition: Voronoi_Rn.h:50
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) ...
Definition: Voronoi_Rn.h:53
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.
Definition: Rn.cpp:29
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.
Definition: Voronoi_Rn.h:48
Refers to the class IncrementalVoronoi.
Definition: Voronoi_Rn.h:46
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.
Definition: Voronoi_Rn.h:44
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.
Definition: Voronoi_Rn.h:102
Refers to the class CellByCellDistNgb ( 10% of the distances between seeds to be sorted) ...
Definition: Voronoi_Rn.h:51
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...
Definition: Voronoi_Rn.h:219
static int Translate(boost::shared_ptr< Polytope_Rn > &pol, const boost::numeric::ublas::vector< double > &v2t)
Translate a polytope by the given vector.
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) ...
Definition: Voronoi_Rn.h:52
#define TEST_KO
Definition: politopixAPI.h:39
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.
#define TEST_OK
Definition: politopixAPI.h:38