Voronoi Diagram#
Synopsis#
Voronoi diagram.
Results#
Warning
Fix Errors Example contains errors needed to be fixed for proper output.
Output:
Seed No.0: At (25,25)
Boundary Vertices List (in order):4,8,7,0,1,
Neighbors (Seed No.):2,3,1,
Seed No.1: At (75,25)
Boundary Vertices List (in order):6,11,7,0,2,
Neighbors (Seed No.):2,4,0,
Seed No.2: At (50,50)
Boundary Vertices List (in order):2,0,1,3,
Neighbors (Seed No.):0,1,3,4,
Seed No.3: At (25,75)
Boundary Vertices List (in order):5,9,4,1,3,
Neighbors (Seed No.):2,0,4,
Seed No.4: At (75,75)
Boundary Vertices List (in order):5,10,6,2,3,
Neighbors (Seed No.):2,3,1,
Vertices Informations:
Vertices No.0: At (50,25)
Vertices No.1: At (25,50)
Vertices No.2: At (75,50)
Vertices No.3: At (50,75)
Vertices No.4: At (0,50)
Vertices No.5: At (50,100)
Vertices No.6: At (100,50)
Vertices No.7: At (50,0)
Vertices No.8: At (0,0)
Vertices No.9: At (0,100)
Vertices No.10: At (100,100)
Vertices No.11: At (100,0)
IIBI-SYNAPSE006:bin mseng$ ./VoronoiDiagram
Seed No.0: At (25,25)
Boundary Vertices List (in order):4,8,7,0,1,
Neighbors (Seed No.):2,3,1,
Seed No.1: At (75,25)
Boundary Vertices List (in order):6,11,7,0,2,
Neighbors (Seed No.):2,4,0,
Seed No.2: At (50,50)
Boundary Vertices List (in order):2,0,1,3,
Neighbors (Seed No.):0,1,3,4,
Seed No.3: At (25,75)
Boundary Vertices List (in order):5,9,4,1,3,
Neighbors (Seed No.):2,0,4,
Seed No.4: At (75,75)
Boundary Vertices List (in order):5,10,6,2,3,
Neighbors (Seed No.):2,3,1,
Vertices Informations:
Vertices No.0: At (50,25)
Vertices No.1: At (25,50)
Vertices No.2: At (75,50)
Vertices No.3: At (50,75)
Vertices No.4: At (0,50)
Vertices No.5: At (50,100)
Vertices No.6: At (100,50)
Vertices No.7: At (50,0)
Vertices No.8: At (0,0)
Vertices No.9: At (0,100)
Vertices No.10: At (100,100)
Vertices No.11: At (100,0)
Code#
C++#
#include "itkVoronoiDiagram2DGenerator.h"
#include "itkImageFileWriter.h"
#include "itkVTKPolyDataWriter.h"
int
main()
{
constexpr double height = 100;
constexpr double width = 100;
using VoronoiDiagramType = itk::VoronoiDiagram2D<double>;
using VoronoiGeneratorType = itk::VoronoiDiagram2DGenerator<double>;
using PointType = VoronoiDiagramType::PointType;
using CellType = VoronoiDiagramType::CellType;
using CellAutoPointer = VoronoiDiagramType::CellAutoPointer;
using PointIdIterator = CellType::PointIdIterator;
using NeighborIdIterator = VoronoiDiagramType::NeighborIdIterator;
auto voronoiDiagram = VoronoiDiagramType::New();
auto voronoiGenerator = VoronoiGeneratorType::New();
PointType insize;
insize[0] = width;
insize[1] = height;
voronoiGenerator->SetBoundary(insize);
// Create a list of seeds
std::vector<PointType> seeds;
PointType seed0;
seed0[0] = 50;
seed0[1] = 50;
seeds.push_back(seed0);
PointType seed1;
seed1[0] = 25;
seed1[1] = 25;
seeds.push_back(seed1);
PointType seed2;
seed2[0] = 75;
seed2[1] = 25;
seeds.push_back(seed2);
PointType seed3;
seed3[0] = 25;
seed3[1] = 75;
seeds.push_back(seed3);
PointType seed4;
seed4[0] = 75;
seed4[1] = 75;
seeds.push_back(seed4);
for (const auto & seed : seeds)
{
voronoiGenerator->AddOneSeed(seed);
}
voronoiGenerator->Update();
voronoiDiagram = voronoiGenerator->GetOutput();
for (unsigned int i = 0; i < seeds.size(); ++i)
{
PointType currP = voronoiDiagram->GetSeed(i);
std::cout << "Seed No." << i << ": At (" << currP[0] << "," << currP[1] << ")" << std::endl;
std::cout << " Boundary Vertices List (in order):";
CellAutoPointer currCell;
voronoiDiagram->GetCellId(i, currCell);
PointIdIterator currCellP;
for (currCellP = currCell->PointIdsBegin(); currCellP != currCell->PointIdsEnd(); ++currCellP)
{
std::cout << (*currCellP) << ",";
}
std::cout << std::endl;
std::cout << " Neighbors (Seed No.):";
NeighborIdIterator currNeibor;
for (currNeibor = voronoiDiagram->NeighborIdsBegin(i); currNeibor != voronoiDiagram->NeighborIdsEnd(i);
++currNeibor)
{
std::cout << (*currNeibor) << ",";
}
std::cout << std::endl << std::endl;
}
std::cout << "Vertices Informations:" << std::endl;
VoronoiDiagramType::VertexIterator allVerts;
int j = 0;
for (allVerts = voronoiDiagram->VertexBegin(); allVerts != voronoiDiagram->VertexEnd(); ++allVerts)
{
std::cout << "Vertices No." << j;
j++;
std::cout << ": At (" << (allVerts.Value())[0] << "," << (allVerts.Value())[1] << ")" << std::endl;
}
// Write the resulting mesh
using WriterType = itk::VTKPolyDataWriter<VoronoiDiagramType::Superclass>;
auto vtkPolyDataWriter = WriterType::New();
vtkPolyDataWriter->SetInput(voronoiDiagram);
vtkPolyDataWriter->SetFileName("voronoi.vtk");
vtkPolyDataWriter->Update();
// Setup an image to visualize the input
{
using ImageType = itk::Image<unsigned char, 2>;
ImageType::IndexType start;
start.Fill(0);
ImageType::SizeType size;
size.Fill(100);
ImageType::RegionType region(start, size);
auto image = ImageType::New();
image->SetRegions(region);
image->Allocate();
image->FillBuffer(0);
ImageType::IndexType ind;
ind[0] = 50;
ind[1] = 50;
image->SetPixel(ind, 255);
ind[0] = 25;
ind[1] = 25;
image->SetPixel(ind, 255);
ind[0] = 75;
ind[1] = 25;
image->SetPixel(ind, 255);
ind[0] = 25;
ind[1] = 75;
image->SetPixel(ind, 255);
ind[0] = 75;
ind[1] = 75;
image->SetPixel(ind, 255);
itk::WriteImage(image, "image.png");
}
return EXIT_SUCCESS;
}
Classes demonstrated#
-
template<typename TCoordType>
class VoronoiDiagram2D : public itk::Mesh<TCoordType, 2, DefaultDynamicMeshTraits<TCoordType, 2, 2, TCoordType>> Implements the 2-Dimensional Voronoi Diagram.
Given a set of seed points, the Voronoi Diagram partitions the plane into regions, each region is a collection of all pixels that is closest to one particular seed point than to other seed points. VoronoiDiagram2D is a mesh structure for storing the Voronoi Diagram, can be Generated by itkVoronoiDiagram2DGenerator.
Template parameters for VoronoiDiagram2D:
TCoordType = The type associated with the coordination of the seeds and the resulting vertices.
- ITK Sphinx Examples:
-
template<typename TCoordType>
class VoronoiDiagram2DGenerator : public itk::MeshSource<VoronoiDiagram2D<TCoordType>> Implement the Sweep Line Algorithm for the construction of the 2D Voronoi Diagram.
Detailed information on this method can be found in: “A sweepline algorithm for Voronoi diagrams.” S. Fortune, Algorithmica 2, 153-174, 1987.
Input parameters are: (1) Size of the region. (2) Seed points coordinates. These coordinates can also be randomly set.
- ITK Sphinx Examples:
- Template Parameters
TCoordType
: The type associated with the coordination of the seeds and the resulting vertices.