Go to the documentation of this file.
24 #include <boost/shared_ptr.hpp>
62 Voronoi_Rn(
const boost::shared_ptr<Polytope_Rn>& inputSpace,
const std::vector<Point_Rn>& listOfPoints);
76 boost::shared_ptr<HalfSpace_Rn>
computeMidPlane(std::vector<Point_Rn>::const_iterator seed1, std::vector<Point_Rn>::const_iterator seed2);
87 void dump(std::ostream& out)
const;
89 void gnuplot(std::ostream& out)
const;
214 std::vector< std::pair< double, unsigned int > >& allDistances,
215 std::vector<Point_Rn>::const_iterator newVoronoiSeed,
216 boost::shared_ptr<Polytope_Rn> newVoronoiCell);
255 virtual bool compute(std::vector< Point_Rn >::const_iterator iteBeg, std::vector< Point_Rn >::const_iterator iteEnd);
257 const std::vector< std::vector< std::pair<boost::shared_ptr<HalfSpace_Rn>,
unsigned int> > >&
264 const std::vector< std::pair< double, unsigned int > >& allDistances,
265 boost::shared_ptr<Polytope_Rn> newVoronoiCell,
266 unsigned int counter,
unsigned int newTruncationStep);
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;
273 this_ostream <<
"Cell " << j <<
" :" << std::endl;
278 this_ostream << std::endl;
287 std::vector< std::vector< std::pair<boost::shared_ptr<HalfSpace_Rn>,
unsigned int> > >
_facetNeighbours;
301 void addFacetAndItsVertices(boost::shared_ptr<HalfSpace_Rn> hs, std::vector< boost::shared_ptr<Generator_Rn> >& listOfGen) {
310 void dump(std::ostream &this_ostream)
const {
311 this_ostream << std::endl <<
"# NON CONVEX POLYTOPE" << std::endl;
314 this_ostream <<
"Facet " << j <<
" :" << std::endl;
316 this_ostream << std::endl;
322 this_ostream << std::endl;
360 std::pair< boost::shared_ptr<HalfSpace_Rn>,
unsigned int > curHalfspace_ngbPolytope(hs, ngbCell);
384 bool fuseCellsWithSameProperty(std::vector< std::vector< boost::shared_ptr<Polytope_Rn> > >& newListOfNonConvexPolytopes);
387 void propagate(std::vector< bool >& listOfAlreadyProcessedSeeds,
unsigned int nbCell,
unsigned int nbProperty);
CellByCellVoronoiDistNgb(const boost::shared_ptr< Polytope_Rn > &is, const std::vector< Point_Rn > &lop)
Constructor.
std::vector< boost::shared_ptr< HalfSpace_Rn > > _listOfHalfSpaces
The list of half-spaces defining the polyhedron.
void gnuplot(std::ostream &out) const
virtual ~CellByCellVoronoiDistNgb()
Destructor.
const boost::shared_ptr< Polytope_Rn > & _inputSpace
The original space to be divided.
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...
void propagate(std::vector< bool > &listOfAlreadyProcessedSeeds, unsigned int nbCell, unsigned int nbProperty)
@ CellByCellDist
Refers to the class CellByCellVoronoi.
@ CellByCellDistNgb_1pc
Refers to the class CellByCellDistNgb (100% of the distances between seeds to be sorted)
Computes the Voronoi diagram cell by cell.
bool checkTopologyAndGeometry() const
@ WithProperties
Refers to the class CellByCellDistNgb ( 40% of the distances between seeds to be sorted)
At the end of the ith iteration, the Voronoi diagram of i seeds is computed.
std::vector< boost::shared_ptr< HalfSpace_Rn > > _allHalfSpaces
For each seed, the list of the half-spaces used to build its cell.
void allocateFacetNeighbours()
virtual bool compute()
Run the whole algorithm.
void addHalfSpaceAndNeighbour(unsigned int cellCounter, boost::shared_ptr< HalfSpace_Rn > hs, unsigned int ngbCell)
CellByCellVoronoiDist(const boost::shared_ptr< Polytope_Rn > &is, const std::vector< Point_Rn > &lop)
Constructor.
@ CellByCellDistNgb_40pc
Refers to the class CellByCellDistNgb ( 30% of the distances between seeds to be sorted)
const std::vector< std::vector< std::pair< boost::shared_ptr< HalfSpace_Rn >, unsigned int > > > & getFacetNeighbours() const
std::vector< int > _listOfNonConvexPolytopeProperties
virtual bool compute()
Run the whole algorithm.
std::vector< boost::shared_ptr< NonConvexPolytope > > _listOfNonConvexPolytopes
The list of properties associated to the seeds.
Compute a n-dimensional Voronoi diagram. It is a partitioning of a space into regions based on distan...
CellByCellVoronoiDist(const boost::shared_ptr< Polytope_Rn > &is, const std::vector< Point_Rn > &lop, double percent)
virtual void updateFacetNeighbours(const std::vector< std::pair< double, unsigned int > > &allDistances, boost::shared_ptr< Polytope_Rn > newVoronoiCell, unsigned int counter, unsigned int newTruncationStep)
virtual ~CellByCellVoronoiDist()
Destructor.
virtual ~Voronoi_Rn()
Destructor.
const std::vector< Point_Rn > & _listOfSeeds
The list of input points.
Fuse Voronoi cells when they share the same property and then get non convex polytopes.
boost::shared_ptr< NonConvexPolytope > getNonConvexPolytope(unsigned int idx)
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),...
const std::vector< Point_Rn > & getListOfSeeds() const
@ CellByCellDistNgb_10pc
Refers to the class CellByCellDistNgb ( 1% of the distances between seeds to be sorted)
CellByCellVoronoiDistNgb(const boost::shared_ptr< Polytope_Rn > &is, const std::vector< Point_Rn > &lop, double percent)
const std::vector< boost::shared_ptr< Polytope_Rn > > & getVoronoiCells() const
IncrementalVoronoi(const boost::shared_ptr< Polytope_Rn > &is, const std::vector< Point_Rn > &lop)
Constructor.
@ CellByCellDistNgb
Refers to the class CellByCellVoronoiDist.
virtual ~VoronoiWithProperties()
Destructor.
const std::vector< std::vector< boost::shared_ptr< Generator_Rn > > > & getListOfGeneratorsPerFacet() const
double _percentageOfSortedHalfSpaces
The amount of half-spaces we want to sort.
void dump(std::ostream &out) const
Dump the cell structure on the given output.
bool fuseCellsWithSameProperty()
Run the whole algorithm: refer to the base class.
std::vector< boost::shared_ptr< Polytope_Rn > > _listOfVoronoiCells
The list of polytopes partitioning the whole space.
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....
void addFacetAndItsVertices(boost::shared_ptr< HalfSpace_Rn > hs, std::vector< boost::shared_ptr< Generator_Rn > > &listOfGen)
CellByCellVoronoi(const boost::shared_ptr< Polytope_Rn > &is, const std::vector< Point_Rn > &lop)
Constructor.
virtual ~CellByCellVoronoi()
Destructor.
void dump(std::ostream &this_ostream) const
virtual bool compute()
Run the whole algorithm.
virtual bool compute()=0
Run the whole algorithm.
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...
virtual bool compute()
Run the whole algorithm.
std::vector< int > _listOfSeedProperties
The list of properties associated to the seeds.
void dumpFacetNeighbours(std::ostream &this_ostream) const
VoronoiWithProperties(const boost::shared_ptr< Polytope_Rn > &is, const std::vector< Point_Rn > &lop, const std::vector< int > &listOfProp, double percent)
virtual ~IncrementalVoronoi()
Destructor.
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...
unsigned int numberOfNonConvexPolytope() const
void dumpFacetNeighbours(std::ostream &this_ostream) const
Voronoi_Rn(const boost::shared_ptr< Polytope_Rn > &inputSpace, const std::vector< Point_Rn > &listOfPoints)
Refers to the class WithProperties.
TypeOfAlgorithm
Choose the kind of algorithm to be run.
@ CellByCellDistNgb_30pc
Refers to the class CellByCellDistNgb ( 20% of the distances between seeds to be sorted)
void dumpNeighbours(std::ostream &this_ostream) const
Do the same as CellByCellVoronoiDist with the neighbourhood computation.
@ CellByCell
Refers to the class IncrementalVoronoi.
const std::vector< boost::shared_ptr< HalfSpace_Rn > > & getListOfHalfSpaces() const
std::vector< boost::shared_ptr< Polytope_Rn > > getVoronoiCells()
virtual void allocateFacetNeighbours()
unsigned int getNonConvexPolytopeProperty(unsigned int idx)
@ CellByCellDistNgb_20pc
Refers to the class CellByCellDistNgb ( 10% of the distances between seeds to be sorted)
virtual ~NonConvexPolytope()
void dumpNgb(std::ostream &this_ostream) const