politopix  4.1.0
 All Classes Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
Voronoi_Rn.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) 2011-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 
23 #include "DoubleDescription_Rn.h"
24 #include "Voronoi_Rn.h"
25 #include "IO_Polytope.h"
26 #include <math.h>
27 #include <numeric>
28 #include <sstream>
29 #include <iostream>
30 #include <fstream>
31 #include <sstream>
32 
33 
34 // Constructor
35 Voronoi_Rn::Voronoi_Rn(const boost::shared_ptr<Polytope_Rn>& inputSpace,
36  const std::vector<Point_Rn>& listOfSeeds):_inputSpace(inputSpace),_listOfSeeds(listOfSeeds) {
37  _listOfVoronoiCells.resize(_listOfSeeds.size());
38 }
39 
40 boost::shared_ptr<HalfSpace_Rn> Voronoi_Rn::computeMidPlane(std::vector<Point_Rn>::const_iterator iteSeed1, std::vector<Point_Rn>::const_iterator iteSeed2) {
41  unsigned int dimension = iteSeed1->dimension();
42  boost::shared_ptr<HalfSpace_Rn> HS;
43  HS.reset(new HalfSpace_Rn(dimension));
44  double sum_square = 0.;
45  for (unsigned int coord_count=0; coord_count<dimension; coord_count++) {
46  HS->setCoefficient(coord_count, iteSeed1->getCoordinate(coord_count)-iteSeed2->getCoordinate(coord_count));
47  double sq = iteSeed1->getCoordinate(coord_count)*iteSeed1->getCoordinate(coord_count);
48  sq -= iteSeed2->getCoordinate(coord_count)*iteSeed2->getCoordinate(coord_count);
49  sq /= 2.;
50  sum_square += sq;
51  }
52  HS->setConstant(-sum_square);
53  return HS;
54 }
55 
56 // CHECK POLYHEDRON
57 bool Voronoi_Rn::checkTopologyAndGeometry() const throw (std::domain_error) {
58  return true;
59 }
60 
61 // Dump the cell structure on the given output.
62 void Voronoi_Rn::dump(std::ostream& out) const {
63  out << std::endl << "#VORONOI DIAGRAM" << std::endl;
64  unsigned int nb=0;
65  std::vector< boost::shared_ptr<Polytope_Rn> >::const_iterator itePC=_listOfVoronoiCells.begin();
66  {for (; itePC!=_listOfVoronoiCells.end(); ++itePC) {
67  out << "# Cell "<< nb << std::endl;
68  (*itePC)->dump(out);
69  out << std::endl;
70 #ifdef DEBUG
71  //std::ostringstream stream_C;
72  //stream_C << "cell";
73  //stream_C << nb;
74  //stream_C << ".ptop";
75  //std::string cellfile = stream_C.str();
76  //IO_Polytope::save(cellfile, *itePC);
77 #endif
78  ++nb;
79  }}
80 }
81 
82 void Voronoi_Rn::gnuplot(std::ostream& out) const throw (std::domain_error) {
83  if (Rn::getDimension() != 2)
84  throw std::domain_error("Voronoi_Rn::gnuplot(std::ostream& out) dimension is not 2 ");
85  unsigned int nb = 1;
86  // To code an unique color for each polygon.
87  std::vector< boost::shared_ptr<Polytope_Rn> >::const_iterator iteVC=_listOfVoronoiCells.begin();
88  {for (; iteVC!=_listOfVoronoiCells.end(); ++iteVC) {
89  double frac = (double)nb / (double)(*iteVC)->numberOfGenerators();
90  std::ostringstream stream_;
91  stream_ << nb;
92  std::string valString = stream_.str();
93  Visualization::gnuplot2D(*iteVC, valString, frac, out);
94  ++nb;
95  }}
96 }
97 
98 bool IncrementalVoronoi::compute() throw (std::length_error) {
99 
100  if (_listOfSeeds.empty() == true)
101  throw std::length_error("IncrementalVoronoi::compute() needs seeds as input!");
102 
103  if (_listOfSeeds.size() == 1) {
104  // Make a deep copy of the current input space
105  boost::shared_ptr<Polytope_Rn> is(new Polytope_Rn( *_inputSpace ));
106  _listOfVoronoiCells[0] = is;
107  return true;
108  }
109 
110  _listOfVoronoiCells.clear();
111  // The first cell is the whole input space.
112  boost::shared_ptr<Polytope_Rn> is(new Polytope_Rn( *_inputSpace ));
113  _listOfVoronoiCells.push_back(is);
114  //unsigned int alreadyProcessedSeeds = 1;
115  unsigned int counter = 0;
116  {for (std::vector<Point_Rn>::const_iterator newVoronoiSeed = (_listOfSeeds.begin()+1); newVoronoiSeed!=_listOfSeeds.end(); ++newVoronoiSeed) {
117 
118 #ifdef DEBUG
119  std::cout << std::endl;
120  std::cout << "Seed " << newVoronoiSeed-_listOfSeeds.begin() << ": ";
121  (newVoronoiSeed)->save(std::cout);
122  std::cout << std::endl;
123 #endif
124  // NEW CELL //
127  // Construction of the current cell as a first copy of the input space, then to be chopped.
128  boost::shared_ptr<Polytope_Rn> newVoronoiCell(new Polytope_Rn(*_inputSpace));
129  unsigned int newTruncationStep = newVoronoiCell->numberOfHalfSpaces();
130  {for (std::vector<Point_Rn>::const_iterator processedSeed = _listOfSeeds.begin(); processedSeed!=newVoronoiSeed; ++processedSeed) {
131  // Compute all the half-spaces of the current cell associated to newVoronoiSeed.
132  boost::shared_ptr<HalfSpace_Rn> HS = computeMidPlane(newVoronoiSeed, processedSeed);
133  newVoronoiCell->addHalfSpace(HS);
134  }}
135  // Now the truncation of the new cell
136  constIteratorOfListOfGeometricObjects< boost::shared_ptr<HalfSpace_Rn> > lexmin_ite(newVoronoiCell->getListOfHalfSpaces());
137  StrongRedundancyProcessing< boost::shared_ptr<PolyhedralCone_Rn> > NRP(newVoronoiCell->neigbourhoodCondition());
139  boost::shared_ptr<PolyhedralCone_Rn>,
142  DD(newVoronoiCell, lexmin_ite, NRP, newTruncationStep);
143 
145  // OLD CELLS //
147  // Chopping of the already created cells.
148  std::vector<Point_Rn>::const_iterator processedSeed;
149  std::vector< boost::shared_ptr<Polytope_Rn> >::iterator oldVoronoiCell;
150  {for (oldVoronoiCell = _listOfVoronoiCells.begin(), processedSeed = _listOfSeeds.begin();
151  oldVoronoiCell!=_listOfVoronoiCells.end(); ++oldVoronoiCell, ++processedSeed) {
152 #ifdef DEBUG
153  //gnuplot(std::cout);
154  //(newVoronoiSeed)->save(std::cout);
155  //std::cout << "### Voronoi cell ###";
156  //(*oldVoronoiCell)->dump(std::cout);
157  //std::cout << std::endl;
158 #endif
159  // Build the half-spaces the other way
160  boost::shared_ptr<HalfSpace_Rn> HS = computeMidPlane(processedSeed, newVoronoiSeed);
161  unsigned int oldTruncationStep = (*oldVoronoiCell)->numberOfHalfSpaces();
162  (*oldVoronoiCell)->addHalfSpace(HS);
163  // Now the truncation of the old cells
164  constIteratorOfListOfGeometricObjects< boost::shared_ptr<HalfSpace_Rn> > lexmin_ite2((*oldVoronoiCell)->getListOfHalfSpaces());
165  StrongRedundancyProcessing< boost::shared_ptr<PolyhedralCone_Rn> > NRP2((*oldVoronoiCell)->neigbourhoodCondition());
167  boost::shared_ptr<PolyhedralCone_Rn>,
168  constIteratorOfListOfGeometricObjects< boost::shared_ptr<HalfSpace_Rn> >,
170  DD2(*oldVoronoiCell, lexmin_ite2, NRP2, oldTruncationStep);
171  }}
172  _listOfVoronoiCells.push_back(newVoronoiCell);
173  //std::cout << counter << std::endl;
174  ++counter;
175  }}
176 #ifdef DEBUG
177  //gnuplot(std::cout);
178  //(newVoronoiSeed)->save(std::cout);
179  //std::cout << "### Voronoi cell ###";
180  //(*oldVoronoiCell)->dump(std::cout);
181  //std::cout << std::endl;
182  //dump(std::cout);
183 #endif
184 
185  return true;
186 }
187 
188 bool CellByCellVoronoi::compute() throw (std::length_error) {
189 
190  if (_listOfSeeds.empty() == true)
191  throw std::length_error("CellByCellVoronoi::compute() needs seeds as input!");
192 
193  if (_listOfSeeds.size() == 1) {
194  // Make a deep copy of the current input space
195  boost::shared_ptr<Polytope_Rn> is(new Polytope_Rn( *_inputSpace ));
196  _listOfVoronoiCells[0] = is;
197  return true;
198  }
199 
200  //unsigned int alreadyProcessedSeeds = 1;
201  unsigned int counter = 0;
202  {for (std::vector<Point_Rn>::const_iterator newVoronoiSeed = _listOfSeeds.begin(); newVoronoiSeed!=_listOfSeeds.end(); ++newVoronoiSeed) {
203 
204 #ifdef DEBUG
205  std::cout << std::endl;
206  std::cout << "Seed " << newVoronoiSeed-_listOfSeeds.begin() << ": ";
207  (newVoronoiSeed)->save(std::cout);
208  std::cout << std::endl;
209 #endif
210  // NEW CELL //
213  // Construction of the current cell as a first copy of the input space, then to be chopped.
214  boost::shared_ptr<Polytope_Rn> newVoronoiCell(new Polytope_Rn(*_inputSpace));
215  unsigned int newTruncationStep = newVoronoiCell->numberOfHalfSpaces();
216  {for (std::vector<Point_Rn>::const_iterator currentSeed = _listOfSeeds.begin(); currentSeed!=_listOfSeeds.end(); ++currentSeed) {
217  if (currentSeed != newVoronoiSeed) {
218  // Compute all the half-spaces of the current cell associated to newVoronoiSeed.
219  boost::shared_ptr<HalfSpace_Rn> HS = computeMidPlane(newVoronoiSeed, currentSeed);
220  newVoronoiCell->addHalfSpace(HS);
221  }
222  }}
223  // Now the truncation of the new cell
224  constIteratorOfListOfGeometricObjects< boost::shared_ptr<HalfSpace_Rn> > lexmin_ite(newVoronoiCell->getListOfHalfSpaces());
225  StrongRedundancyProcessing< boost::shared_ptr<PolyhedralCone_Rn> > NRP(newVoronoiCell->neigbourhoodCondition());
227  boost::shared_ptr<PolyhedralCone_Rn>,
230  DD(newVoronoiCell, lexmin_ite, NRP, newTruncationStep);
231 
232  _listOfVoronoiCells[counter] = newVoronoiCell;
233  //std::cout << counter << std::endl;
234  ++counter;
235  }}
236 #ifdef DEBUG
237  //gnuplot(std::cout);
238  //(newVoronoiSeed)->save(std::cout);
239  //std::cout << "### Voronoi cell ###";
240  //(*oldVoronoiCell)->dump(std::cout);
241  //std::cout << std::endl;
242  //dump(std::cout);
243 #endif
244 
245  return true;
246 }
247 
248 //bool CellByCellVoronoiNgb::compute() throw (std::length_error) {
249 
250  //if (_listOfSeeds.empty() == true)
251  //throw std::length_error("ParallelVoronoi::compute() needs seeds as input!");
252 
253  //if (_listOfSeeds.size() == 1) {
254  // Make a deep copy of the current input space
255  //boost::shared_ptr<Polytope_Rn> is(new Polytope_Rn( *_inputSpace ));
256  //_listOfVoronoiCells.push_back(is);
257  //return true;
258  //}
259  //unsigned int RnDIM = Rn::getDimension();
260  //std::vector< std::pair<unsigned int,unsigned int> > v(RnDIM);
261  //_facetNeighbours.assign(_listOfSeeds.size(), v);
262  //_facetNeighbours.reserve(_listOfSeeds.size());
263 
264  //unsigned int alreadyProcessedSeeds = 1;
265  //unsigned int counter = 0;
266  //{for (std::vector<Point_Rn>::const_iterator newVoronoiSeed = _listOfSeeds.begin(); newVoronoiSeed!=_listOfSeeds.end(); ++newVoronoiSeed) {
267 
268 //#ifdef DEBUG
269  //std::cout << std::endl;
270  //std::cout << "Seed " << newVoronoiSeed-_listOfSeeds.begin() << ": ";
271  //(newVoronoiSeed)->save(std::cout);
272  //std::cout << std::endl;
273 //#endif
275  // NEW CELL //
277  // Build the array of all half-spaces.
278  //std::vector< boost::shared_ptr<HalfSpace_Rn> > allHalfSpaces;
279  // Construction of the current cell as a first copy of the input space, then to be chopped.
280  //boost::shared_ptr<Polytope_Rn> newVoronoiCell(new Polytope_Rn(*_inputSpace));
281  //unsigned int newTruncationStep = newVoronoiCell->numberOfHalfSpaces();
282  //{for (std::vector<Point_Rn>::const_iterator currentSeed=_listOfSeeds.begin(); currentSeed!=_listOfSeeds.end(); ++currentSeed) {
283  //if (currentSeed != newVoronoiSeed) {
284  // Compute all the half-spaces of the current cell associated to newVoronoiSeed.
285  //boost::shared_ptr<HalfSpace_Rn> HS = computeMidPlane(newVoronoiSeed, currentSeed);
286  //newVoronoiCell->addHalfSpace(HS);
287  //allHalfSpaces.push_back(HS);
288  //}
289  //else {
290  // Build an artificial half-space that will always be respected: 1 + 0.x1 + ... + 0.xn >= 0.
291  // This is not to change the numbering and should have no influence on the process.
292  //boost::shared_ptr<HalfSpace_Rn> HS;
293  //HS.reset(new HalfSpace_Rn(RnDIM));
294  //HS->setConstant(1.);
295  // The normal vector
296  //for (unsigned int coord_count=0; coord_count<RnDIM; coord_count++)
297  //HS->setCoefficient(coord_count, 0.);
298  //allHalfSpaces.push_back(HS);
299  //}
300  //}}
301  // Before we truncate the cell we want to swap half-spaces: put first the half-spaces that are known to give genuine facets.
302  //if (!_facetNeighbours[counter].empty()) {
303  // Half-spaces numbered from 0 to (newTruncationStep-1) belong to the input space only
304  //std::vector< boost::shared_ptr<HalfSpace_Rn> >::iterator hs2translate = allHalfSpaces.begin()+newTruncationStep;
305  //for (unsigned int i=0; i<_facetNeighbours[counter].size(); ++i) {
306  // Reminder _facetNeighbours[counter][i].first refers to a half-space number and _facetNeighbours[counter][i].second to a neighbour polytope number
307  //std::vector< boost::shared_ptr<HalfSpace_Rn> >::iterator hs_ngb = allHalfSpaces.begin() + _facetNeighbours[counter][i].first;
308  //iter_swap(hs2translate, hs_ngb);
309  //++hs2translate;
310  //}
311  //}
312  // Now the truncation of the new cell
313  //constIteratorOfListOfGeometricObjects< boost::shared_ptr<HalfSpace_Rn> > lexmin_ite(newVoronoiCell->getListOfHalfSpaces());
314  //StrongRedundancyProcessing< boost::shared_ptr<PolyhedralCone_Rn> > NRP;
315  //DoubleDescription<
316  //boost::shared_ptr<PolyhedralCone_Rn>,
317  //constIteratorOfListOfGeometricObjects< boost::shared_ptr<HalfSpace_Rn> >,
318  //StrongRedundancyProcessing< boost::shared_ptr<PolyhedralCone_Rn> > >
319  //DD(newVoronoiCell, lexmin_ite, NRP, newTruncationStep);
320 
321  //_listOfVoronoiCells.push_back(newVoronoiCell);
322 
324  //std::cout << counter << std::endl;
325  //std::ostringstream stream_;
326  //stream_ << "./ptop/cellVor";
327  //stream_ << counter;
328  //stream_ << ".ptop";
329  //std::string valString = stream_.str();
330  //std::cout << valString << std::endl;
331  //IO_Polytope::save(valString, newVoronoiCell);
333 
334  // Compute all the facet neighbors of newVoronoiCell.
335  // Don't forget that the first newTruncationStep half-spaces belong to the input space.
336  // Half-spaces numbered from 0 to (newTruncationStep-1) belong to the input space only
337  // and are not relevant from the neighborhood point of view.
338  //updateFacetNeighbours(allHalfSpaces, newVoronoiCell, counter, newTruncationStep);
339 
340  //++counter;
341 
342  //}}
343 //#ifdef DEBUG
344  //gnuplot(std::cout);
345  //(newVoronoiSeed)->save(std::cout);
346  //std::cout << "### Voronoi cell ###";
347  //(*oldVoronoiCell)->dump(std::cout);
348  //std::cout << std::endl;
349  //dump(std::cout);
350 //#endif
351 
352  //return true;
353 //}
354 
355 //void CellByCellVoronoiNgb::updateFacetNeighbours(
356  //const std::vector< boost::shared_ptr<HalfSpace_Rn> >& allHalfSpaces,
357  //boost::shared_ptr<Polytope_Rn> newVoronoiCell,
358  //unsigned int counter,
359  //unsigned int newTruncationStep) {
360  //constIteratorOfListOfGeometricObjects< boost::shared_ptr<HalfSpace_Rn> > iteHS_Vorcell(newVoronoiCell->getListOfHalfSpaces());
361  //unsigned int hs_vorcell_counter = 0;
362  //unsigned int hs_all_counter = 0;
363  //{for (iteHS_Vorcell.begin(); iteHS_Vorcell.end()!=true; iteHS_Vorcell.next()) {
364  //const boost::shared_ptr<HalfSpace_Rn>& currentHalfSpace = iteHS_Vorcell.current();
365  //std::vector< boost::shared_ptr<HalfSpace_Rn> >::const_iterator iteAll;
366  //{for (iteAll=allHalfSpaces.begin(); iteAll!=allHalfSpaces.end(); ++iteAll) {
367  //if (iteHS_Vorcell.current() == *iteAll)
368  //if (hs_all_counter >= newTruncationStep) {
369  //std::pair< unsigned int, unsigned int > curHalfspace_ngbPolytope(hs_vorcell_counter, hs_all_counter);
370  // counter has already been incremented here so decrease it by one.
371  //_facetNeighbours[counter].push_back(curHalfspace_ngbPolytope);//-newTruncationStep);
372  //}
373  //++hs_all_counter;
374  //}}
375  //++hs_vorcell_counter;
376  //}}
377 //}
378 
380  bool operator() (std::pair<double, unsigned int> pt1, std::pair<double, unsigned int> pt2) { return (pt1.first < pt2.first);}
382 
384  std::vector< std::pair< double, unsigned int > >& allDistances,
385  std::vector<Point_Rn>::const_iterator newVoronoiSeed,
386  boost::shared_ptr<Polytope_Rn> newVoronoiCell)
387 {
388  // Determine all distances from currentSeed to the others
389  allDistances.clear();
390  std::pair< double, unsigned int > v;
391  allDistances.assign(_listOfSeeds.size(), v);
392 
393  unsigned int distCounter = 0;
394  {for (std::vector<Point_Rn>::const_iterator currentSeed = _listOfSeeds.begin(); currentSeed!=_listOfSeeds.end(); ++currentSeed) {
395  if (currentSeed != newVoronoiSeed) {
396  boost::numeric::ublas::vector<double> v1 = currentSeed->vect();
397  v1 -= newVoronoiSeed->vect();
398  allDistances[distCounter].first = norm_2(v1); // The distance from newVoronoiSeed to currentSeed
399  }
400  else {
401  allDistances[distCounter].first = 0.; // The distance from newVoronoiSeed to itself
402  }
403  // We number the seeds because the array is going to be shuffled.
404  allDistances[distCounter].second = distCounter;
405  ++distCounter;
406  }}
407 
408  // Sort allDistances using thisSortOnFirstCoord as comparator
409  if (_percentageOfSortedHalfSpaces == 100.)
410  std::sort(allDistances.begin(), allDistances.end(), thisSortOnFirstCoord);
411  else {
412  unsigned int nbOfSortedSeeds = _listOfSeeds.size() * (_percentageOfSortedHalfSpaces/100.);
413  std::partial_sort(allDistances.begin(), allDistances.begin()+nbOfSortedSeeds, allDistances.end(), thisSortOnFirstCoord);
414  }
415 
417  //std::cout << std::endl;
418  //{for (unsigned int i=1; i<allDistances.size(); ++i) {
419  //std::cout << allDistances[i][0] << " " << allDistances[i][1] << std::endl;
420  //}}
421  //std::cout << std::endl;
423 
424  // Put inside _allHalfSpaces all the half-space that could play a role in shaping newVoronoiSeed
425  _allHalfSpaces.clear();
426  constIteratorOfListOfGeometricObjects< boost::shared_ptr<HalfSpace_Rn> > iteHSnewVoronoiSeed(newVoronoiCell->getListOfHalfSpaces());
427  for (iteHSnewVoronoiSeed.begin(); iteHSnewVoronoiSeed.end()!=true; iteHSnewVoronoiSeed.next())
428  _allHalfSpaces.push_back(iteHSnewVoronoiSeed.current());
429  {for (unsigned int i=1; i<allDistances.size(); ++i) {
430  // i=0 <=> newVoronoiSeed
431  std::vector<Point_Rn>::const_iterator currentSeed = _listOfSeeds.begin()+allDistances[i].second;
432  // Compute all the half-spaces of the current cell associated to newVoronoiSeed.
433  boost::shared_ptr<HalfSpace_Rn> HS = computeMidPlane(newVoronoiSeed, currentSeed);
434  newVoronoiCell->addHalfSpace(HS);
435  _allHalfSpaces.push_back(HS);
436  }}
437 }
438 
439 bool CellByCellVoronoiDist::compute() throw (std::length_error) {
440 
441  if (_listOfSeeds.empty() == true)
442  throw std::length_error("CellByCellVoronoiDist::compute() needs seeds as input!");
443 
444  if (_listOfSeeds.size() == 1) {
445  // Make a deep copy of the current input space
446  boost::shared_ptr<Polytope_Rn> is(new Polytope_Rn( *_inputSpace ));
447  _listOfVoronoiCells[0] = is;
448  return true;
449  }
450 
451  //unsigned int alreadyProcessedSeeds = 1;
452  unsigned int counter = 0;
453  {for (std::vector<Point_Rn>::const_iterator newVoronoiSeed = _listOfSeeds.begin(); newVoronoiSeed!=_listOfSeeds.end(); ++newVoronoiSeed) {
454 
455 #ifdef DEBUG
456  std::cout << std::endl;
457  std::cout << "Seed " << newVoronoiSeed-_listOfSeeds.begin() << ": ";
458  (newVoronoiSeed)->save(std::cout);
459  std::cout << std::endl;
460 #endif
461  // NEW CELL //
464  // Construction of the current cell as a first copy of the input space, then to be chopped.
465  boost::shared_ptr<Polytope_Rn> newVoronoiCell(new Polytope_Rn(*_inputSpace));
466  unsigned int newTruncationStep = newVoronoiCell->numberOfHalfSpaces();
467 
468  // Order the half-spaces inside newVoronoiCell according to sorting criteria based on distances.
469  std::vector< std::pair< double, unsigned int > > allDistances;
470  processDistances(allDistances, newVoronoiSeed, newVoronoiCell);
471 
472  // Now the truncation of the new cell
473  constIteratorOfListOfGeometricObjects< boost::shared_ptr<HalfSpace_Rn> > lexmin_ite(newVoronoiCell->getListOfHalfSpaces());
474  StrongRedundancyProcessing< boost::shared_ptr<PolyhedralCone_Rn> > NRP(newVoronoiCell->neigbourhoodCondition());
476  boost::shared_ptr<PolyhedralCone_Rn>,
479  DD(newVoronoiCell, lexmin_ite, NRP, newTruncationStep);
480 
481  _listOfVoronoiCells[counter] = newVoronoiCell;
482  //std::cout << counter << std::endl;
483 
484  ++counter;
485  }}
486 #ifdef DEBUG
487  //gnuplot(std::cout);
488  //(newVoronoiSeed)->save(std::cout);
489  //std::cout << "### Voronoi cell ###";
490  //(*oldVoronoiCell)->dump(std::cout);
491  //std::cout << std::endl;
492  //dump(std::cout);
493 #endif
494 
495  return true;
496 }
497 
498 bool CellByCellVoronoiDistNgb::compute() throw (std::length_error) {
499  std::vector< Point_Rn >::const_iterator iteBeg = _listOfSeeds.begin();
500  std::vector< Point_Rn >::const_iterator iteEnd = _listOfSeeds.end();
501  return CellByCellVoronoiDistNgb::compute(iteBeg, iteEnd);
502 }
503 
505  std::vector< Point_Rn >::const_iterator iteBeg,
506  std::vector< Point_Rn >::const_iterator iteEnd) throw (std::length_error) {
507 
508  if (_listOfSeeds.empty() == true)
509  throw std::length_error("CellByCellVoronoiDistNgb::compute() needs seeds as input!");
510 
511  if (_listOfSeeds.size() == 1) {
512  // Make a deep copy of the current input space
513  boost::shared_ptr<Polytope_Rn> is(new Polytope_Rn( *_inputSpace ));
514  _listOfVoronoiCells[0] = is;
515  return true;
516  }
517 
518  allocateFacetNeighbours();
519 
520  //unsigned int alreadyProcessedSeeds = 1;
521  unsigned int counter = iteBeg - _listOfSeeds.begin();
522  {for (std::vector<Point_Rn>::const_iterator newVoronoiSeed = iteBeg; newVoronoiSeed!=iteEnd; ++newVoronoiSeed) {
523 
524 #ifdef DEBUG
525  std::cout << std::endl;
526  std::cout << "Seed " << newVoronoiSeed-_listOfSeeds.begin() << ": ";
527  (newVoronoiSeed)->save(std::cout);
528  std::cout << std::endl;
529 #endif
530  // NEW CELL //
533  // Construction of the current cell as a first copy of the input space, then to be chopped.
534  boost::shared_ptr<Polytope_Rn> newVoronoiCell(new Polytope_Rn(*_inputSpace));
535  unsigned int newTruncationStep = newVoronoiCell->numberOfHalfSpaces();
536 
537  // Order the half-spaces inside newVoronoiCell according to sorting criteria based on distances.
538  std::vector< std::pair< double, unsigned int > > allDistances;
539  processDistances(allDistances, newVoronoiSeed, newVoronoiCell);
540 
541  // Now the truncation of the new cell
542  constIteratorOfListOfGeometricObjects< boost::shared_ptr<HalfSpace_Rn> > lexmin_ite(newVoronoiCell->getListOfHalfSpaces());
543  StrongRedundancyProcessing< boost::shared_ptr<PolyhedralCone_Rn> > NRP(newVoronoiCell->neigbourhoodCondition());
545  boost::shared_ptr<PolyhedralCone_Rn>,
548  DD(newVoronoiCell, lexmin_ite, NRP, newTruncationStep);
549 
550  _listOfVoronoiCells[counter] = newVoronoiCell;
551  //if (counter%100 == 0)
552  //std::cout << counter << std::endl;
553  // To stop the computations before the end.
554  //if (counter == 500)
555  //return true;
556 
557  updateFacetNeighbours(allDistances, newVoronoiCell, counter, newTruncationStep);
558 
559  ++counter;
560  }}
561 #ifdef DEBUG
562  //gnuplot(std::cout);
563  //(newVoronoiSeed)->save(std::cout);
564  //std::cout << "### Voronoi cell ###";
565  //(*oldVoronoiCell)->dump(std::cout);
566  //std::cout << std::endl;
567  //dump(std::cout);
568 #endif
569 
570  return true;
571 }
572 
573 void CellByCellVoronoiDistNgb::dumpFacetNeighbours(std::ostream &this_ostream) const {
574  std::vector< std::vector< std::pair<boost::shared_ptr<HalfSpace_Rn>, unsigned int> > >::const_iterator iteSeed;
575  // The list of seeds and their associated cells.
576  unsigned int cellCounter = 0;
577  for (iteSeed=_facetNeighbours.begin(); iteSeed!=_facetNeighbours.end(); ++iteSeed) {
578  this_ostream << iteSeed->size() << " ";
579  std::vector< std::pair<boost::shared_ptr<HalfSpace_Rn>, unsigned int> >::const_iterator iteCurHS_NgbNb;
580  for (iteCurHS_NgbNb=iteSeed->begin(); iteCurHS_NgbNb!=iteSeed->end(); ++iteCurHS_NgbNb) {
581  this_ostream << _listOfVoronoiCells[cellCounter]->getHalfSpaceNumber(iteCurHS_NgbNb->first) << " " << iteCurHS_NgbNb->second << " ";
582  }
583  this_ostream << std::endl;
584  ++cellCounter;
585  }
586 }
587 
589  const std::vector< std::pair< double, unsigned int > >& allDistances,
590  boost::shared_ptr<Polytope_Rn> newVoronoiCell,
591  unsigned int counter, unsigned int newTruncationStep) {
592  constIteratorOfListOfGeometricObjects< boost::shared_ptr<HalfSpace_Rn> > iteHS_Vorcell(newVoronoiCell->getListOfHalfSpaces());
593  {for (iteHS_Vorcell.begin(); iteHS_Vorcell.end()!=true; iteHS_Vorcell.next()) {
594  //const boost::shared_ptr<HalfSpace_Rn>& currentHalfSpace = iteHS_Vorcell.current();
595  std::vector< boost::shared_ptr<HalfSpace_Rn> >::const_iterator iteAll;
596  unsigned int hs_all_counter = newTruncationStep;
597  {for (iteAll=_allHalfSpaces.begin()+newTruncationStep; iteAll!=_allHalfSpaces.end(); ++iteAll) {
598  if (iteHS_Vorcell.current() == *iteAll) {
599  // allHalfSpaces has been shuffled by the sorting process so we lost the information about which seeds are connected
600  // to newVoronoiCell. However we can retrieve it in allDistances[hs_all_counter].second
601  // Add -newTruncationStep+1 to hs_all_counter because there a shift between allDistances and _allHalfSpaces
602  // _allHalfSpaces contains first all the input space half-spaces.
603  // Example with a cube as the input space:
604  // _allHalfSpaces: H0 ... H7 H8 ... Hn
605  // allDistances: (0.,counter) (d1,n1)
606  unsigned int ngbSeedNumber = allDistances[hs_all_counter-newTruncationStep+1].second;
607  std::pair< boost::shared_ptr<HalfSpace_Rn>, unsigned int > curHalfspace_ngbPolytope(iteHS_Vorcell.current(), ngbSeedNumber);
608  _facetNeighbours[counter].push_back(curHalfspace_ngbPolytope);//-newTruncationStep);
609  break;
610  }
611  ++hs_all_counter;
612  }}
613  }}
614 }
615 
616 //========================================================================================================================
617 
618 void VoronoiWithProperties::loadProp(const std::string& filename, std::vector< int >& propList) throw (std::ios_base::failure) {
619  std::ifstream file(filename.c_str(), std::ifstream::in);
620  if (!file) {
621  std::string s("VoronoiWithProperties::loadProp() Unable to open ");
622  s += filename;
623  s += "\n";
624  throw std::ios_base::failure(s);
625  }
626 
627  propList.clear();
628  int prop;
629  unsigned int numberOfSeeds, counter=0;
630  std::string line;
631  std::getline(file, line);
632  std::istringstream iline(line);
633  iline >> numberOfSeeds;
634  while (counter != numberOfSeeds) {
635  std::getline(file, line);
636  std::istringstream iline2(line);
637  iline2 >> prop;
638  propList.push_back(prop);
639  ++counter;
640  }
641  file.close();
642 }
643 
644 void VoronoiWithProperties::loadProperties(const std::string& filename) throw (std::ios_base::failure) {
645  std::ifstream file(filename.c_str(), std::ifstream::in);
646  if (!file) {
647  std::string s("VoronoiWithProperties::loadProperties() Unable to open ");
648  s += filename;
649  s += "\n";
650  throw std::ios_base::failure(s);
651  }
652 
653  _listOfSeedProperties.clear();
654  int prop;
655  unsigned int numberOfSeeds, counter=0;
656  std::string line;
657  std::getline(file, line);
658  std::istringstream iline(line);
659  iline >> numberOfSeeds;
660  while (counter != numberOfSeeds) {
661  std::getline(file, line);
662  std::istringstream iline2(line);
663  iline2 >> prop;
664  _listOfSeedProperties.push_back(prop);
665  ++counter;
666  }
667  file.close();
668 }
669 
670 void VoronoiWithProperties::loadSeeds(const std::string& filename) throw (std::ios_base::failure) {
671  boost::shared_ptr<Polytope_Rn> listOfPointsAsVDesc(new Polytope_Rn());
672  IO_Polytope::load(filename, listOfPointsAsVDesc);
673  constIteratorOfListOfGeometricObjects< boost::shared_ptr<Generator_Rn> > iteGN(listOfPointsAsVDesc->getListOfGenerators());
674  {for (iteGN.begin(); iteGN.end()!=true; iteGN.next()) {
675  _listOfSeeds.push_back(Point_Rn(
676  iteGN.current()->getCoordinate(0),
677  iteGN.current()->getCoordinate(1),
678  iteGN.current()->getCoordinate(2)) );
679  }}
680 }
681 
682 void VoronoiWithProperties::loadCell(const std::string& filename) throw (std::ios_base::failure) {
683  boost::shared_ptr<Polytope_Rn> thisCell(new Polytope_Rn());
684  IO_Polytope::load(filename, thisCell);
685  _listOfVoronoiCells.push_back(thisCell);
686 }
687 
688 void VoronoiWithProperties::loadFacetNeighbours(const std::string& filename) throw (std::length_error,std::ios_base::failure) {
689  std::ifstream file(filename.c_str(), std::ifstream::in);
690  if (!file) {
691  std::string s("VoronoiWithProperties::loadFacetNeighbours() Unable to open ");
692  s += filename;
693  s += "\n";
694  throw std::ios_base::failure(s);
695  }
696  if (_listOfVoronoiCells.empty() == true)
697  throw std::length_error("VoronoiWithProperties::loadFacetNeighbours() needs polytopes!");
698 
699  allocateFacetNeighbours();
700  unsigned int numberOfSeeds, cellCounter = 0;
701  std::string line;
702  std::getline(file, line);
703  std::istringstream iline(line);
704  iline >> numberOfSeeds;
705  //std::vector< std::vector< std::pair<boost::shared_ptr<HalfSpace_Rn>, unsigned int> > >::iterator iteSeed;
706  while (cellCounter != numberOfSeeds) {
707  unsigned int nbnbNgb, nbFacet, nbNgb;
708  std::getline(file, line);
709  std::istringstream iline2(line);
710  iline2 >> nbnbNgb;
711  for (unsigned int coord_count=0; coord_count<nbnbNgb; coord_count++) {
712  iline2 >> nbFacet;
713  iline2 >> nbNgb;
714  addHalfSpaceAndNeighbour(cellCounter, _listOfVoronoiCells[cellCounter]->getHalfSpace(nbFacet), nbNgb);
715  }
716  ++cellCounter;
717  }
718  file.close();
719 }
720 
721 void VoronoiWithProperties::dumpNeighbours(std::ostream &this_ostream) const {
722  std::vector< std::vector< std::pair<boost::shared_ptr<HalfSpace_Rn>, unsigned int> > >::const_iterator iteSeed;
723  // The list of seeds and their associated cells.
724  unsigned int cellCounter = 0;
725  for (iteSeed=_facetNeighbours.begin(); iteSeed!=_facetNeighbours.end(); ++iteSeed) {
726  this_ostream << iteSeed-_facetNeighbours.begin() << ": ";
727  std::vector< std::pair<boost::shared_ptr<HalfSpace_Rn>, unsigned int> >::const_iterator iteCurHS_NgbNb;
728  for (iteCurHS_NgbNb=iteSeed->begin(); iteCurHS_NgbNb!=iteSeed->end(); ++iteCurHS_NgbNb) {
729  this_ostream << iteCurHS_NgbNb->second << " ";
730  }
731  this_ostream << std::endl;
732  ++cellCounter;
733  }
734 }
735 
736 void VoronoiWithProperties::dumpFacetNeighbours(std::ostream &this_ostream) const {
737  std::vector< std::vector< std::pair<boost::shared_ptr<HalfSpace_Rn>, unsigned int> > >::const_iterator iteSeed;
738  // The list of seeds and their associated cells.
739  unsigned int cellCounter = 0;
740  for (iteSeed=_facetNeighbours.begin(); iteSeed!=_facetNeighbours.end(); ++iteSeed) {
741  this_ostream << iteSeed->size() << " ";
742  std::vector< std::pair<boost::shared_ptr<HalfSpace_Rn>, unsigned int> >::const_iterator iteCurHS_NgbNb;
743  for (iteCurHS_NgbNb=iteSeed->begin(); iteCurHS_NgbNb!=iteSeed->end(); ++iteCurHS_NgbNb) {
744  this_ostream << _listOfVoronoiCells[cellCounter]->getHalfSpaceNumber(iteCurHS_NgbNb->first) << " " << iteCurHS_NgbNb->second << " ";
745  }
746  this_ostream << std::endl;
747  ++cellCounter;
748  }
749 }
750 
751 bool VoronoiWithProperties::fuseCellsWithSameProperty(std::vector< std::vector< unsigned int > >& newListOfNonConvexPolytopes) throw (std::length_error) {
752  if (_listOfSeedProperties.empty() == true)
753  throw std::length_error("VoronoiWithProperties::fuseCellsWithSameProperty() needs seeds as input!");
754  unsigned int currentSeedNb = 0;
755  std::vector< boost::shared_ptr<Polytope_Rn> >::iterator iteCell;
756  std::vector< std::vector< std::pair<boost::shared_ptr<HalfSpace_Rn>, unsigned int> > >::iterator iteSeed;
757  // During the first step we go through the list of seeds and their associated cells to delete
758  // useless half-spaces separating two cells sharing the same property.
759  for (iteSeed=_facetNeighbours.begin(), iteCell =_listOfVoronoiCells.begin();
760  iteSeed!=_facetNeighbours.end() && iteCell!=_listOfVoronoiCells.end(); ++iteSeed, ++iteCell) {
761  std::vector< std::pair<boost::shared_ptr<HalfSpace_Rn>, unsigned int> >::iterator iteCurHS_NgbNb;
762  for (iteCurHS_NgbNb=iteSeed->begin(); iteCurHS_NgbNb!=iteSeed->end(); ++iteCurHS_NgbNb) {
763  // All the neighbours of iteSeed (iteCurHS_NgbNb->second) and the separating hyperplane in iteCell (iteCurHS_NgbNb->first)
764  unsigned int ngbNb = iteCurHS_NgbNb->second;
765  // Check whether both cells share the same property
766  if (_listOfSeedProperties[currentSeedNb] == _listOfSeedProperties[ngbNb]) {
767  // Some cells may have been already emptied
768  if ((*iteCell)->numberOfHalfSpaces() > 0) {
769  // In that case we need to remove in the current cell the frontier half-space shared with the neighbour cell.
770  boost::shared_ptr<HalfSpace_Rn> HS = iteCurHS_NgbNb->first;
771  // Remove HS from the current polytope
772  (*iteCell)->removeHalfSpaceFromListAndGenerators(HS);
773  }
774  }
775  }
776  ++currentSeedNb;
777  }
778 
779  // At this step some cells have been reduced to the null polytope, others have lost facets and only the isolated ones
780  // remain unchanged. Now use a stack to pick up all the neighbours cells having the same property.
781  currentSeedNb = 0;
782  iteSeed = _facetNeighbours.begin();
783  iteCell = _listOfVoronoiCells.begin();
784  std::set< unsigned int > allProcessedSeeds;
785  newListOfNonConvexPolytopes.clear();
786  while (iteSeed!=_facetNeighbours.end() && iteCell!=_listOfVoronoiCells.end()) {
787  // Do nothing if the seed and its cell have already been taken into account.
788  if (allProcessedSeeds.find(currentSeedNb) == allProcessedSeeds.end()) {
789  //allProcessedSeeds.insert(currentSeedNb);// Not needed as allProcessedSeeds.insert(seedNgbNb) will do the job
790  std::vector< unsigned int > v;
791  newListOfNonConvexPolytopes.push_back(v);
792  // Build the stack of neighbours having the same property
793  std::stack< unsigned int > allNeighboursWithSameProperty;
794  allNeighboursWithSameProperty.push(currentSeedNb);
795  unsigned int currentProperty = _listOfSeedProperties[currentSeedNb];
796  // Mark it as processed now.
797  if (allProcessedSeeds.insert(currentSeedNb).second == true) {
798  // Not inserted yet
799  if (_listOfVoronoiCells[currentSeedNb]->numberOfHalfSpaces() != 0) {
800  // Insert here in the main DS if the associated polytope is not null
801  newListOfNonConvexPolytopes[newListOfNonConvexPolytopes.size()-1].push_back(currentSeedNb);
802  }
803  }
804  while (!allNeighboursWithSameProperty.empty()) {
805  unsigned int topSeed = allNeighboursWithSameProperty.top();
806  allNeighboursWithSameProperty.pop();
807  std::vector< std::pair<boost::shared_ptr<HalfSpace_Rn>, unsigned int> >::iterator iteCurHS_NgbNb;
808  {for (iteCurHS_NgbNb=_facetNeighbours[topSeed].begin(); iteCurHS_NgbNb!=_facetNeighbours[topSeed].end(); ++iteCurHS_NgbNb) {
809  unsigned int ngbSeed = iteCurHS_NgbNb->second;
810  if (allProcessedSeeds.find(ngbSeed) == allProcessedSeeds.end()) {
811  if (_listOfSeedProperties[ngbSeed] == (int)currentProperty) {
812  allNeighboursWithSameProperty.push(ngbSeed);
813  allProcessedSeeds.insert(ngbSeed);
814  if (_listOfVoronoiCells[ngbSeed]->numberOfHalfSpaces() != 0) {
815  // Insert here in the main DS if the associated polytope is not null
816  newListOfNonConvexPolytopes[newListOfNonConvexPolytopes.size()-1].push_back(ngbSeed);
817  }
818  }
819  }
820  }}
821  } // while (!allNeighboursWithSameProperty.empty()) {
822  } // if (allProcessedSeeds.find(currentSeedNb) == allProcessedSeeds.end()) {
823  ++currentSeedNb;
824  ++iteSeed;
825  ++iteCell;
826  }
827 
828  return true;
829 }
830 
831 bool VoronoiWithProperties::fuseCellsWithSameProperty(std::vector< std::vector< boost::shared_ptr<Polytope_Rn> > >& newListOfNonConvexPolytopes)
832  throw (std::length_error)
833 {
834  if (_listOfSeedProperties.empty() == true)
835  throw std::length_error("VoronoiWithProperties::fuseCellsWithSameProperty() needs seeds as input!");
836  std::vector< std::vector< unsigned int > > allNonConvexPolytopeIds;
837  // Compute the ids of the polytopes to be fused.
838  fuseCellsWithSameProperty(allNonConvexPolytopeIds);
839  newListOfNonConvexPolytopes.resize(allNonConvexPolytopeIds.size());
840  std::vector< std::vector< unsigned int > >::const_iterator iteNCPIds = allNonConvexPolytopeIds.begin();
841  unsigned int currentSeedNb = 0;
842  for (; iteNCPIds != allNonConvexPolytopeIds.end(); ++iteNCPIds) {
843  std::vector< unsigned int >::const_iterator iteCurrentNCPId = iteNCPIds->begin();
844  for (; iteCurrentNCPId != iteNCPIds->end(); ++iteCurrentNCPId) {
845  newListOfNonConvexPolytopes[currentSeedNb].push_back(_listOfVoronoiCells[*iteCurrentNCPId]);
846  }
847  ++currentSeedNb;
848  }
849 
850  return true;
851 }
void loadProperties(const std::string &filename)
Definition: Voronoi_Rn.cpp:644
virtual bool compute()
Run the whole algorithm.
Definition: Voronoi_Rn.cpp:498
Creation of a n-coordinate geometric point designed to be shared by its neighbour faces...
Definition: Point_Rn.h:34
std::vector< int > _listOfSeedProperties
The list of properties associated to the seeds.
Definition: Voronoi_Rn.h:348
double _percentageOfSortedHalfSpaces
The amount of half-spaces we want to sort.
Definition: Voronoi_Rn.h:223
virtual bool compute()
Run the whole algorithm.
Definition: Voronoi_Rn.cpp:439
void dump(std::ostream &out) const
Dump the cell structure on the given output.
Definition: Voronoi_Rn.cpp:62
static void loadProp(const std::string &filename, std::vector< int > &propList)
Definition: Voronoi_Rn.cpp:618
virtual bool compute()
Run the whole algorithm.
Definition: Voronoi_Rn.cpp:188
virtual bool compute()
Run the whole algorithm.
Definition: Voronoi_Rn.cpp:98
Model a polytope using its two equivalent definitions : the convex hull and the half-space intersecti...
Definition: Polytope_Rn.h:34
boost::shared_ptr< HalfSpace_Rn > computeMidPlane(std::vector< Point_Rn >::const_iterator seed1, std::vector< Point_Rn >::const_iterator seed2)
Compute the half-space containing seed1, in between seed1 and seed2, according to the growing seed pr...
Definition: Voronoi_Rn.cpp:40
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
bool operator()(std::pair< double, unsigned int > pt1, std::pair< double, unsigned int > pt2)
Definition: Voronoi_Rn.cpp:380
const std::vector< Point_Rn > & _listOfSeeds
The list of input points.
Definition: Voronoi_Rn.h:95
static void gnuplot2D(const boost::shared_ptr< Polytope_Rn > &polygon, const std::string &name, double col, std::ostream &out)
Provide the drwing of polygon under the gnuplot format.
void gnuplot(std::ostream &out) const
Definition: Voronoi_Rn.cpp:82
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
std::vector< boost::shared_ptr< Polytope_Rn > > _listOfVoronoiCells
The list of polytopes partitioning the whole space.
Definition: Voronoi_Rn.h:97
Voronoi_Rn(const boost::shared_ptr< Polytope_Rn > &inputSpace, const std::vector< Point_Rn > &listOfPoints)
Refers to the class WithProperties.
Definition: Voronoi_Rn.cpp:35
struct sortOnFirstCoord thisSortOnFirstCoord
bool checkTopologyAndGeometry() const
Definition: Voronoi_Rn.cpp:57
void processDistances(std::vector< std::pair< double, unsigned int > > &allDistances, std::vector< Point_Rn >::const_iterator newVoronoiSeed, boost::shared_ptr< Polytope_Rn > newVoronoiCell)
The distance criterion is going to be used to classify partially or globally the way in which we proc...
Definition: Voronoi_Rn.cpp:383
const boost::shared_ptr< Polytope_Rn > & _inputSpace
The original space to be divided.
Definition: Voronoi_Rn.h:93
static polito_EXPORT unsigned int getDimension()
Return the dimension of the cartesian space we work in.
Definition: Rn.cpp:29
This class can be more time-consuming than WeakRedundancyProcessing or NoRedundancyProcessing because...
void loadCell(const std::string &filename)
Definition: Voronoi_Rn.cpp:682
virtual void updateFacetNeighbours(const std::vector< std::pair< double, unsigned int > > &allDistances, boost::shared_ptr< Polytope_Rn > newVoronoiCell, unsigned int counter, unsigned int newTruncationStep)
Definition: Voronoi_Rn.cpp:588
This class is designed to run the list of all geometric objects representing a polytope.
void dumpFacetNeighbours(std::ostream &this_ostream) const
Definition: Voronoi_Rn.cpp:736
void dumpFacetNeighbours(std::ostream &this_ostream) const
Definition: Voronoi_Rn.cpp:573
std::vector< std::vector< std::pair< boost::shared_ptr< HalfSpace_Rn >, unsigned int > > > _facetNeighbours
For the current cell, store its neighbour numbers. _facetNeighbours[i] = { (H0, i0), (H1, i1), ... }.
Definition: Voronoi_Rn.h:346
void loadFacetNeighbours(const std::string &filename)
Definition: Voronoi_Rn.cpp:688
std::vector< std::vector< std::pair< boost::shared_ptr< HalfSpace_Rn >, unsigned int > > > _facetNeighbours
For the current cell, store its neighbour numbers. _facetNeighbours[i] = { (H0, i0), (H1, i1), ... }.
Definition: Voronoi_Rn.h:265
std::vector< boost::shared_ptr< HalfSpace_Rn > > _allHalfSpaces
For each seed, the list of the half-spaces used to build its cell.
Definition: Voronoi_Rn.h:225
The algorithm implemented here is an incremental algorithm as mentioned in How Good are Convex Hull ...
void loadSeeds(const std::string &filename)
Definition: Voronoi_Rn.cpp:670
std::vector< boost::shared_ptr< Polytope_Rn > > _listOfVoronoiCells
The list of polytopes partitioning the whole space.
Definition: Voronoi_Rn.h:344
void dumpNeighbours(std::ostream &this_ostream) const
Definition: Voronoi_Rn.cpp:721
bool fuseCellsWithSameProperty(std::vector< std::vector< unsigned int > > &newListOfNonConvexPolytopes)
Run the whole algorithm: refer to the base class.
Definition: Voronoi_Rn.cpp:751