politopix  5.0.0
Voronoi_Rn.h
Go to the documentation of this file.
1 // politopix allows to make computations on polytopes such as finding vertices, intersecting, Minkowski sums, ...
2 // Copyright (C) 2011-2022 : Delos Vincent
3 //
4 // This program is free software: you can redistribute it and/or modify
5 // it under the terms of the GNU General Public License as published by
6 // the Free Software Foundation, either version 3 of the License, or
7 // (at your option) any later version.
8 //
9 // This program is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 // GNU General Public License for more details.
13 //
14 // You should have received a copy of the GNU General Public License
15 // along with this program. If not, see <http://www.gnu.org/licenses/>.
16 //
19 // I2M (UMR CNRS 5295 / University of Bordeaux)
20 
21 #ifndef VORONOI_Rn
22 #define VORONOI_Rn
23 
24 #include <boost/shared_ptr.hpp>
25 #include <stdexcept>
26 #include <exception>
27 #include <iostream>
28 #include <utility>
29 #include <vector>
30 #include "Rn.h"
31 #include "Generator_Rn.h"
32 #include "HalfSpace_Rn.h"
33 #include "Polytope_Rn.h"
34 
35 
39 class Voronoi_Rn {
40 
41  public:
42 
46  CellByCell = 1,
55 
57  //Voronoi_Rn() { _inputSpace=0; }
58 
62  Voronoi_Rn(const boost::shared_ptr<Polytope_Rn>& inputSpace, const std::vector<Point_Rn>& listOfPoints);
63 
65  virtual ~Voronoi_Rn() {}
66 
68  virtual bool compute()=0;
69 
70  //virtual void updateFacetNeighbours(const std::vector< boost::shared_ptr<HalfSpace_Rn> >& allHalfSpaces,
71  //boost::shared_ptr<Polytope_Rn> newVoronoiCell, unsigned int counter,
72  //unsigned int newTruncationStep)=0;
73 
74 
76  boost::shared_ptr<HalfSpace_Rn> computeMidPlane(std::vector<Point_Rn>::const_iterator seed1, std::vector<Point_Rn>::const_iterator seed2);
77 
78  const std::vector< boost::shared_ptr<Polytope_Rn> >& getVoronoiCells() const {return _listOfVoronoiCells;}
79  std::vector< boost::shared_ptr<Polytope_Rn> > getVoronoiCells() {return _listOfVoronoiCells;}
80 
81  const std::vector< Point_Rn >& getListOfSeeds() const {return _listOfSeeds;}
82 
83  // CHECK POLYHEDRON
84  bool checkTopologyAndGeometry() const;
85 
87  void dump(std::ostream& out) const;
88 
89  void gnuplot(std::ostream& out) const;
90 
91  protected:
93  const boost::shared_ptr<Polytope_Rn>& _inputSpace;
95  const std::vector< Point_Rn >& _listOfSeeds;
97  std::vector< boost::shared_ptr<Polytope_Rn> > _listOfVoronoiCells;
98 };
99 
103 
104  public:
106  //IncrementalVoronoi() { _inputSpace=0; }
107 
111  IncrementalVoronoi(const boost::shared_ptr<Polytope_Rn>& is, const std::vector<Point_Rn>& lop):Voronoi_Rn(is,lop) {}
112 
114  virtual ~IncrementalVoronoi() {}
115 
117  virtual bool compute();
118 
119  //virtual void updateFacetNeighbours(const std::vector< boost::shared_ptr<HalfSpace_Rn> >& allHalfSpaces,
120  //boost::shared_ptr<Polytope_Rn> newVoronoiCell, unsigned int counter,
121  //unsigned int newTruncationStep) {}
122 
123 };
124 
128 
129  public:
131  //CellByCellVoronoi() { _inputSpace=0; }
132 
136  CellByCellVoronoi(const boost::shared_ptr<Polytope_Rn>& is, const std::vector<Point_Rn>& lop):Voronoi_Rn(is,lop) {}
137 
139  virtual ~CellByCellVoronoi() {}
140 
142  virtual bool compute();
143 
144  //virtual void updateFacetNeighbours(const std::vector< boost::shared_ptr<HalfSpace_Rn> >& allHalfSpaces,
145  //boost::shared_ptr<Polytope_Rn> newVoronoiCell, unsigned int counter,
146  //unsigned int newTruncationStep) {}
147 
148 };
149 
150 // \class CellByCellVoronoiNgb
151 // \brief Computes the Voronoi diagram cell by cell and the neighbourhood relationships between them.
152 //class CellByCellVoronoiNgb : public CellByCellVoronoi {
153 
154 // public:
156  //CellByCellVoronoiNgb() { _inputSpace=0; }
157 
161  //CellByCellVoronoiNgb(const boost::shared_ptr<Polytope_Rn>& is, const std::vector<Point_Rn>& lop):CellByCellVoronoi(is,lop) {}
162 
164  //virtual ~CellByCellVoronoiNgb() {}
165 
167  //virtual bool compute();
168 
170  //virtual void updateFacetNeighbours(
171  //const std::vector< boost::shared_ptr<HalfSpace_Rn> >& allHalfSpaces,
172  //boost::shared_ptr<Polytope_Rn> newVoronoiCell, unsigned int counter,
173  //unsigned int newTruncationStep);
174 
175  //protected:
179  //std::vector< std::vector< std::pair<unsigned int,unsigned int> > > _facetNeighbours;
180 
181 //};
182 
186 
187  public:
189  //CellByCellVoronoiDist() { _inputSpace=0; }
190 
194  CellByCellVoronoiDist(const boost::shared_ptr<Polytope_Rn>& is, const std::vector<Point_Rn>& lop):CellByCellVoronoi(is,lop) {
196  }
197 
202  CellByCellVoronoiDist(const boost::shared_ptr<Polytope_Rn>& is, const std::vector<Point_Rn>& lop, double percent):CellByCellVoronoi(is,lop) {
204  }
205 
208 
210  virtual bool compute();
211 
213  void processDistances(
214  std::vector< std::pair< double, unsigned int > >& allDistances,
215  std::vector<Point_Rn>::const_iterator newVoronoiSeed,
216  boost::shared_ptr<Polytope_Rn> newVoronoiCell);
217 
220 
221  protected:
225  std::vector< boost::shared_ptr<HalfSpace_Rn> > _allHalfSpaces;
226 };
227 
228 
232 
233  public:
235  //CellByCellVoronoiDistNgb() { _inputSpace=0; }
236 
240  CellByCellVoronoiDistNgb(const boost::shared_ptr<Polytope_Rn>& is, const std::vector<Point_Rn>& lop):CellByCellVoronoiDist(is,lop) {}
241 
246  CellByCellVoronoiDistNgb(const boost::shared_ptr<Polytope_Rn>& is, const std::vector<Point_Rn>& lop, double percent):CellByCellVoronoiDist(is,lop,percent) {}
247 
250 
252  virtual bool compute();
253 
255  virtual bool compute(std::vector< Point_Rn >::const_iterator iteBeg, std::vector< Point_Rn >::const_iterator iteEnd);
256 
257  const std::vector< std::vector< std::pair<boost::shared_ptr<HalfSpace_Rn>, unsigned int> > >&
259 
260  virtual void allocateFacetNeighbours() { _facetNeighbours.resize(_listOfSeeds.size()); }
261  void dumpFacetNeighbours(std::ostream &this_ostream) const;
262 
263  virtual void updateFacetNeighbours(
264  const std::vector< std::pair< double, unsigned int > >& allDistances,
265  boost::shared_ptr<Polytope_Rn> newVoronoiCell,
266  unsigned int counter, unsigned int newTruncationStep);
267 
268  void dumpNgb(std::ostream &this_ostream) const {
269  this_ostream << std::endl << "# NEIGHBOURS FOR EACH CELL" << std::endl;
270  this_ostream << "Number of facets :" << _facetNeighbours.size() << std::endl;
271  {
272  for (unsigned int j = 0; j < _facetNeighbours.size(); ++j) {
273  this_ostream << "Cell " << j << " :" << std::endl;
274  {
275  for (unsigned int k = 0; k < _facetNeighbours[j].size(); ++k) {
276  _facetNeighbours[j][k].first->dump(this_ostream);
277  this_ostream << ", Ngb = " << _facetNeighbours[j][k].second;
278  this_ostream << std::endl;
279  }
280  }
281  } // for (unsigned int j = 0; j < _facetNeighbours.size(); ++j) {
282  }
283  }
284 
285  protected:
287  std::vector< std::vector< std::pair<boost::shared_ptr<HalfSpace_Rn>, unsigned int> > > _facetNeighbours;
288 
289 };
290 
291 
294 
295  public:
296 
298 
299  virtual ~NonConvexPolytope() {}
300 
301  void addFacetAndItsVertices(boost::shared_ptr<HalfSpace_Rn> hs, std::vector< boost::shared_ptr<Generator_Rn> >& listOfGen) {
302  _listOfHalfSpaces.push_back(hs);
303  _listOfGeneratorsPerFacet.push_back(listOfGen);
304  }
305 
306  const std::vector< boost::shared_ptr<HalfSpace_Rn> >& getListOfHalfSpaces() const {return _listOfHalfSpaces;}
307 
308  const std::vector< std::vector< boost::shared_ptr<Generator_Rn> > >& getListOfGeneratorsPerFacet() const {return _listOfGeneratorsPerFacet;}
309 
310  void dump(std::ostream &this_ostream) const {
311  this_ostream << std::endl << "# NON CONVEX POLYTOPE" << std::endl;
312  this_ostream << "Number of facets : " << _listOfHalfSpaces.size() << std::endl;
313  for (unsigned int j = 0; j < _listOfHalfSpaces.size(); ++j) {
314  this_ostream << "Facet " << j << " :" << std::endl;
315  _listOfHalfSpaces[j]->dump(this_ostream);
316  this_ostream << std::endl;
317  {
318  for (unsigned int k = 0; k < _listOfGeneratorsPerFacet[j].size(); ++k) {
319  _listOfGeneratorsPerFacet[j][k]->dump(this_ostream);
320  this_ostream << " ";
321  }
322  this_ostream << std::endl;
323  }
324  }
325  }
326 
327  protected:
329  std::vector< boost::shared_ptr<HalfSpace_Rn> > _listOfHalfSpaces;
331  std::vector< std::vector< boost::shared_ptr<Generator_Rn> > > _listOfGeneratorsPerFacet;
332 
333 };
334 
335 
339 
340  public:
341 
346  VoronoiWithProperties(const boost::shared_ptr<Polytope_Rn>& is, const std::vector<Point_Rn>& lop, const std::vector<int>& listOfProp, double percent):CellByCellVoronoiDistNgb(is,lop,percent) {
347  _listOfSeedProperties = listOfProp;
349  }
350 
353 
355 
356  void dumpNeighbours(std::ostream &this_ostream) const;
357  void dumpFacetNeighbours(std::ostream &this_ostream) const;
358 
359  void addHalfSpaceAndNeighbour(unsigned int cellCounter, boost::shared_ptr<HalfSpace_Rn> hs, unsigned int ngbCell) {
360  std::pair< boost::shared_ptr<HalfSpace_Rn>, unsigned int > curHalfspace_ngbPolytope(hs, ngbCell);
361  _facetNeighbours[cellCounter].push_back(curHalfspace_ngbPolytope);
362  }
363 
365  // Non-convex geometry //
367  unsigned int numberOfNonConvexPolytope() const { return _listOfNonConvexPolytopes.size(); }
368 
369  boost::shared_ptr< NonConvexPolytope > getNonConvexPolytope(unsigned int idx) { return _listOfNonConvexPolytopes[idx]; }
370 
371  unsigned int getNonConvexPolytopeProperty(unsigned int idx) { return _listOfNonConvexPolytopeProperties[idx]; }
372 
373  // Unite different Voronoi cells or polytopes which share at least one common facet.
375  //virtual bool compute();
376 
379 
381  bool fuseCellsWithSameProperty(std::vector< std::vector< unsigned int > >& newListOfNonConvexPolytopes);
382 
384  bool fuseCellsWithSameProperty(std::vector< std::vector< boost::shared_ptr<Polytope_Rn> > >& newListOfNonConvexPolytopes);
385 
386  protected:
387  void propagate(std::vector< bool >& listOfAlreadyProcessedSeeds, unsigned int nbCell, unsigned int nbProperty);
388 
389  protected:
391  std::vector< int > _listOfSeedProperties;
393  std::vector< boost::shared_ptr< NonConvexPolytope > > _listOfNonConvexPolytopes;
395 
396 };
397 
398 #endif // VORONOI_Rn
CellByCellVoronoiDistNgb::CellByCellVoronoiDistNgb
CellByCellVoronoiDistNgb(const boost::shared_ptr< Polytope_Rn > &is, const std::vector< Point_Rn > &lop)
Constructor.
Definition: Voronoi_Rn.h:240
CellByCellVoronoiDist
Constructor.
Definition: Voronoi_Rn.h:185
NonConvexPolytope::_listOfHalfSpaces
std::vector< boost::shared_ptr< HalfSpace_Rn > > _listOfHalfSpaces
The list of half-spaces defining the polyhedron.
Definition: Voronoi_Rn.h:329
Voronoi_Rn::gnuplot
void gnuplot(std::ostream &out) const
Definition: Voronoi_Rn.cpp:82
CellByCellVoronoiDistNgb::~CellByCellVoronoiDistNgb
virtual ~CellByCellVoronoiDistNgb()
Destructor.
Definition: Voronoi_Rn.h:249
Voronoi_Rn::_inputSpace
const boost::shared_ptr< Polytope_Rn > & _inputSpace
The original space to be divided.
Definition: Voronoi_Rn.h:93
CellByCellVoronoiDist::setAverageNumberOfNeighbours
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
VoronoiWithProperties::propagate
void propagate(std::vector< bool > &listOfAlreadyProcessedSeeds, unsigned int nbCell, unsigned int nbProperty)
Definition: Voronoi_Rn.cpp:851
Voronoi_Rn::CellByCellDist
@ CellByCellDist
Refers to the class CellByCellVoronoi.
Definition: Voronoi_Rn.h:47
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
Voronoi_Rn::checkTopologyAndGeometry
bool checkTopologyAndGeometry() const
Definition: Voronoi_Rn.cpp:57
Voronoi_Rn::WithProperties
@ WithProperties
Refers to the class CellByCellDistNgb ( 40% of the distances between seeds to be sorted)
Definition: Voronoi_Rn.h:54
IncrementalVoronoi
At the end of the ith iteration, the Voronoi diagram of i seeds is computed.
Definition: Voronoi_Rn.h:102
CellByCellVoronoiDist::_allHalfSpaces
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
VoronoiWithProperties::allocateFacetNeighbours
void allocateFacetNeighbours()
Definition: Voronoi_Rn.h:354
CellByCellVoronoiDist::compute
virtual bool compute()
Run the whole algorithm.
Definition: Voronoi_Rn.cpp:439
VoronoiWithProperties::addHalfSpaceAndNeighbour
void addHalfSpaceAndNeighbour(unsigned int cellCounter, boost::shared_ptr< HalfSpace_Rn > hs, unsigned int ngbCell)
Definition: Voronoi_Rn.h:359
CellByCellVoronoiDist::CellByCellVoronoiDist
CellByCellVoronoiDist(const boost::shared_ptr< Polytope_Rn > &is, const std::vector< Point_Rn > &lop)
Constructor.
Definition: Voronoi_Rn.h:194
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
CellByCellVoronoiDistNgb::getFacetNeighbours
const std::vector< std::vector< std::pair< boost::shared_ptr< HalfSpace_Rn >, unsigned int > > > & getFacetNeighbours() const
Definition: Voronoi_Rn.h:258
VoronoiWithProperties::_listOfNonConvexPolytopeProperties
std::vector< int > _listOfNonConvexPolytopeProperties
Definition: Voronoi_Rn.h:394
IncrementalVoronoi::compute
virtual bool compute()
Run the whole algorithm.
Definition: Voronoi_Rn.cpp:98
VoronoiWithProperties::_listOfNonConvexPolytopes
std::vector< boost::shared_ptr< NonConvexPolytope > > _listOfNonConvexPolytopes
The list of properties associated to the seeds.
Definition: Voronoi_Rn.h:393
Voronoi_Rn
Compute a n-dimensional Voronoi diagram. It is a partitioning of a space into regions based on distan...
Definition: Voronoi_Rn.h:39
CellByCellVoronoiDist::CellByCellVoronoiDist
CellByCellVoronoiDist(const boost::shared_ptr< Polytope_Rn > &is, const std::vector< Point_Rn > &lop, double percent)
Definition: Voronoi_Rn.h:202
CellByCellVoronoiDistNgb::updateFacetNeighbours
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
CellByCellVoronoiDist::~CellByCellVoronoiDist
virtual ~CellByCellVoronoiDist()
Destructor.
Definition: Voronoi_Rn.h:207
Voronoi_Rn::~Voronoi_Rn
virtual ~Voronoi_Rn()
Destructor.
Definition: Voronoi_Rn.h:65
Voronoi_Rn::_listOfSeeds
const std::vector< Point_Rn > & _listOfSeeds
The list of input points.
Definition: Voronoi_Rn.h:95
VoronoiWithProperties
Fuse Voronoi cells when they share the same property and then get non convex polytopes.
Definition: Voronoi_Rn.h:338
VoronoiWithProperties::getNonConvexPolytope
boost::shared_ptr< NonConvexPolytope > getNonConvexPolytope(unsigned int idx)
Definition: Voronoi_Rn.h:369
CellByCellVoronoiDistNgb::_facetNeighbours
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),...
Definition: Voronoi_Rn.h:287
Voronoi_Rn::getListOfSeeds
const std::vector< Point_Rn > & getListOfSeeds() const
Definition: Voronoi_Rn.h:81
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
CellByCellVoronoiDistNgb::CellByCellVoronoiDistNgb
CellByCellVoronoiDistNgb(const boost::shared_ptr< Polytope_Rn > &is, const std::vector< Point_Rn > &lop, double percent)
Definition: Voronoi_Rn.h:246
NonConvexPolytope::NonConvexPolytope
NonConvexPolytope()
Definition: Voronoi_Rn.h:297
Voronoi_Rn::getVoronoiCells
const std::vector< boost::shared_ptr< Polytope_Rn > > & getVoronoiCells() const
Definition: Voronoi_Rn.h:78
IncrementalVoronoi::IncrementalVoronoi
IncrementalVoronoi(const boost::shared_ptr< Polytope_Rn > &is, const std::vector< Point_Rn > &lop)
Constructor.
Definition: Voronoi_Rn.h:111
Voronoi_Rn::CellByCellDistNgb
@ CellByCellDistNgb
Refers to the class CellByCellVoronoiDist.
Definition: Voronoi_Rn.h:48
Rn.h
VoronoiWithProperties::~VoronoiWithProperties
virtual ~VoronoiWithProperties()
Destructor.
Definition: Voronoi_Rn.h:352
NonConvexPolytope::getListOfGeneratorsPerFacet
const std::vector< std::vector< boost::shared_ptr< Generator_Rn > > > & getListOfGeneratorsPerFacet() const
Definition: Voronoi_Rn.h:308
CellByCellVoronoiDist::_percentageOfSortedHalfSpaces
double _percentageOfSortedHalfSpaces
The amount of half-spaces we want to sort.
Definition: Voronoi_Rn.h:223
Voronoi_Rn::dump
void dump(std::ostream &out) const
Dump the cell structure on the given output.
Definition: Voronoi_Rn.cpp:62
Generator_Rn.h
VoronoiWithProperties::fuseCellsWithSameProperty
bool fuseCellsWithSameProperty()
Run the whole algorithm: refer to the base class.
Definition: Voronoi_Rn.cpp:900
Voronoi_Rn::_listOfVoronoiCells
std::vector< boost::shared_ptr< Polytope_Rn > > _listOfVoronoiCells
The list of polytopes partitioning the whole space.
Definition: Voronoi_Rn.h:97
NonConvexPolytope::_listOfGeneratorsPerFacet
std::vector< std::vector< boost::shared_ptr< Generator_Rn > > > _listOfGeneratorsPerFacet
For each facet the list of vertices, they are pointers on the former voronoi cells i....
Definition: Voronoi_Rn.h:331
Voronoi_Rn::Incremental
@ Incremental
Definition: Voronoi_Rn.h:45
NonConvexPolytope::addFacetAndItsVertices
void addFacetAndItsVertices(boost::shared_ptr< HalfSpace_Rn > hs, std::vector< boost::shared_ptr< Generator_Rn > > &listOfGen)
Definition: Voronoi_Rn.h:301
CellByCellVoronoi::CellByCellVoronoi
CellByCellVoronoi(const boost::shared_ptr< Polytope_Rn > &is, const std::vector< Point_Rn > &lop)
Constructor.
Definition: Voronoi_Rn.h:136
CellByCellVoronoi::~CellByCellVoronoi
virtual ~CellByCellVoronoi()
Destructor.
Definition: Voronoi_Rn.h:139
NonConvexPolytope::dump
void dump(std::ostream &this_ostream) const
Definition: Voronoi_Rn.h:310
CellByCellVoronoi::compute
virtual bool compute()
Run the whole algorithm.
Definition: Voronoi_Rn.cpp:188
NonConvexPolytope
Definition: Voronoi_Rn.h:293
Voronoi_Rn::compute
virtual bool compute()=0
Run the whole algorithm.
Voronoi_Rn::computeMidPlane
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
CellByCellVoronoiDistNgb::compute
virtual bool compute()
Run the whole algorithm.
Definition: Voronoi_Rn.cpp:498
VoronoiWithProperties::_listOfSeedProperties
std::vector< int > _listOfSeedProperties
The list of properties associated to the seeds.
Definition: Voronoi_Rn.h:391
CellByCellVoronoiDistNgb::dumpFacetNeighbours
void dumpFacetNeighbours(std::ostream &this_ostream) const
Definition: Voronoi_Rn.cpp:573
VoronoiWithProperties::VoronoiWithProperties
VoronoiWithProperties(const boost::shared_ptr< Polytope_Rn > &is, const std::vector< Point_Rn > &lop, const std::vector< int > &listOfProp, double percent)
Definition: Voronoi_Rn.h:346
IncrementalVoronoi::~IncrementalVoronoi
virtual ~IncrementalVoronoi()
Destructor.
Definition: Voronoi_Rn.h:114
CellByCellVoronoiDist::processDistances
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
VoronoiWithProperties::numberOfNonConvexPolytope
unsigned int numberOfNonConvexPolytope() const
Definition: Voronoi_Rn.h:367
VoronoiWithProperties::dumpFacetNeighbours
void dumpFacetNeighbours(std::ostream &this_ostream) const
Definition: Voronoi_Rn.cpp:736
Voronoi_Rn::Voronoi_Rn
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
Polytope_Rn.h
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
VoronoiWithProperties::dumpNeighbours
void dumpNeighbours(std::ostream &this_ostream) const
Definition: Voronoi_Rn.cpp:721
CellByCellVoronoiDistNgb
Do the same as CellByCellVoronoiDist with the neighbourhood computation.
Definition: Voronoi_Rn.h:231
Voronoi_Rn::CellByCell
@ CellByCell
Refers to the class IncrementalVoronoi.
Definition: Voronoi_Rn.h:46
NonConvexPolytope::getListOfHalfSpaces
const std::vector< boost::shared_ptr< HalfSpace_Rn > > & getListOfHalfSpaces() const
Definition: Voronoi_Rn.h:306
Voronoi_Rn::getVoronoiCells
std::vector< boost::shared_ptr< Polytope_Rn > > getVoronoiCells()
Definition: Voronoi_Rn.h:79
CellByCellVoronoiDistNgb::allocateFacetNeighbours
virtual void allocateFacetNeighbours()
Definition: Voronoi_Rn.h:260
VoronoiWithProperties::getNonConvexPolytopeProperty
unsigned int getNonConvexPolytopeProperty(unsigned int idx)
Definition: Voronoi_Rn.h:371
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
HalfSpace_Rn.h
NonConvexPolytope::~NonConvexPolytope
virtual ~NonConvexPolytope()
Definition: Voronoi_Rn.h:299
CellByCellVoronoiDistNgb::dumpNgb
void dumpNgb(std::ostream &this_ostream) const
Definition: Voronoi_Rn.h:268