voronoi_logger.cpp
1 #include <fstream>
2 #include "voronoi_logger.hpp"
3 #include "../misc/simple_io.hpp"
4 #include "VoronoiMesh.hpp"
5 
6 using namespace voronoi_loggers;
7 
9 
10 void VoronoiLogger::output(const VoronoiMesh& /*v*/) {}
11 
12 void VoronoiLogger::output(const Tessellation& /*v*/) {}
13 
14 VoronoiLogger::~VoronoiLogger(void) {}
15 
16 BinLogger::BinLogger(string const& file_name):
17 file_name_(file_name) {}
18 
20 {
21  ofstream file_handle(file_name_.c_str(),ios::binary);
22 
23  binary_write_single_int(static_cast<int>(v.GetTotalSidesNumber()),file_handle);
24 
25  for(int i=0;i<static_cast<int>(v.GetTotalSidesNumber());++i)
26  {
27  Edge edge=v.GetEdge(i);
28  binary_write_single_double(edge.vertices.first.x,file_handle);
29  binary_write_single_double(edge.vertices.second.x,file_handle);
30  }
31 
32  for(int i=0;i<static_cast<int>(v.GetTotalSidesNumber());++i)
33  {
34  Edge edge=v.GetEdge(i);
35  binary_write_single_double(edge.vertices.first.y,file_handle);
36  binary_write_single_double(edge.vertices.second.y,file_handle);
37  }
38 
39  for(int i=0;i<static_cast<int>(v.GetTotalSidesNumber());++i)
40  {
41  Edge edge=v.GetEdge(i);
42  binary_write_single_int(edge.neighbors.first,file_handle);
43  binary_write_single_int(edge.neighbors.second,file_handle);
44  }
45 
46  binary_write_single_int(v.GetPointNo(),file_handle);
47 
48  for(int i=0;i<static_cast<int>(v.GetPointNo());++i)
49  {
50  binary_write_single_double(v.GetMeshPoint(i).x,file_handle);
51  binary_write_single_double(v.GetMeshPoint(i).y,file_handle);
52  }
53 
54  for(int i=0;i<static_cast<int>(v.GetPointNo());++i)
55  {
56  const vector<int>& indices = v.GetCellEdges(i);
57  binary_write_single_int(static_cast<int>(indices.size()),file_handle);
58  for(int j=0;j<static_cast<int>(indices.size());++j)
59  binary_write_single_int(indices[static_cast<size_t>(j)],file_handle);
60  }
61  file_handle.close();
62 }
63 
65 {
66  ofstream file_handle(file_name_.c_str(),ios::binary);
67 
68  binary_write_single_int(static_cast<int>(v.GetTotalSidesNumber()),file_handle);
69 
70  for(int i=0;i<static_cast<int>(v.GetTotalSidesNumber());++i)
71  {
72  Edge edge=v.GetEdge(i);
73  binary_write_single_double(edge.vertices.first.x,file_handle);
74  binary_write_single_double(edge.vertices.second.x,file_handle);
75  }
76 
77  for(int i=0;i<static_cast<int>(v.GetTotalSidesNumber());++i)
78  {
79  Edge edge=v.GetEdge(i);
80  binary_write_single_double(edge.vertices.first.y,file_handle);
81  binary_write_single_double(edge.vertices.second.y,file_handle);
82  }
83 
84  for(int i=0;i<static_cast<int>(v.GetTotalSidesNumber());++i)
85  {
86  Edge edge=v.GetEdge(i);
87  binary_write_single_int(edge.neighbors.first,file_handle);
88  binary_write_single_int(edge.neighbors.second,file_handle);
89  }
90 
91  binary_write_single_int(v.GetPointNo(),file_handle);
92 
93  for(int i=0;i<static_cast<int>(v.GetPointNo());++i)
94  {
95  binary_write_single_double(v.GetMeshPoint(i).x,file_handle);
96  binary_write_single_double(v.GetMeshPoint(i).y,file_handle);
97  }
98 
99  for(int i=0;i<static_cast<int>(v.GetPointNo());++i)
100  {
101  const vector<int>& indices = v.GetCellEdges(i);
102  binary_write_single_int(static_cast<int>(indices.size()),file_handle);
103  for(int j=0;j<static_cast<int>(indices.size());++j)
104  binary_write_single_int(indices[static_cast<size_t>(j)],file_handle);
105  }
106  file_handle.close();
107 }
108 
109 vector<Vector2D> BinLogger::read(string location)
110 {
111  fstream myFile (location.c_str(),ios::in | ios::binary);
112  if(!myFile.good())
113  throw UniversalError("Error opening voronoi logger file!!");
114  int N,itemp;
115  myFile.read(reinterpret_cast<char*>(&N),sizeof (int));
116  double temp;
117  for(int i=0;i<N*4;++i)
118  myFile.read(reinterpret_cast<char*>(&temp),sizeof(double));
119  for (int i = 0; i<N *2; ++i)
120  myFile.read(reinterpret_cast<char*>(&itemp), sizeof(int));
121  myFile.read(reinterpret_cast<char*>(&N),sizeof (int));
122  vector<Vector2D> res(static_cast<size_t>(N));
123  for(size_t i=0;i<static_cast<size_t>(N);++i)
124  {
125  myFile.read(reinterpret_cast<char*>(&res[i].x),sizeof(double));
126  myFile.read(reinterpret_cast<char*>(&res[i].y),sizeof(double));
127  }
128  return res;
129 }
Voronoi tessellation class.
Definition: VoronoiMesh.hpp:29
Edge const & GetEdge(int index) const
Returns edge (interface between cells)
vector< Vector2D > read(string location)
Reads the output information from the Voronoi tessellation.
Abstract class for tessellation.
virtual int GetPointNo(void) const =0
Get Total number of mesh generating points.
Container for error reports.
void binary_write_single_double(double d, ofstream &fh)
Writes a double to a binary file.
Definition: simple_io.cpp:80
Interface between two cells.
Definition: Edge.hpp:13
virtual Vector2D GetMeshPoint(int index) const =0
Returns Position of mesh generating point.
virtual void output(const VoronoiMesh &v)
Outputs information from the Voronoi tessellation.
void output(const VoronoiMesh &v)
Outputs information from the Voronoi tessellation.
double y
Component in the y direction.
Definition: geometry.hpp:48
A debugging function for the Voronoi tessellation.
virtual Edge const & GetEdge(int index) const =0
Returns edge (interface between cells)
int GetPointNo(void) const
Get Total number of mesh generating points.
std::pair< Vector2D, Vector2D > vertices
Points at the ends of the edge.
Definition: Edge.hpp:18
vector< int > const & GetCellEdges(int index) const
Returns a reference to a vector<int> containing the indexes of the edges related to a cell...
void binary_write_single_int(int n, ofstream &fh)
Writes a single integer to a binary file.
Definition: simple_io.cpp:75
int GetTotalSidesNumber(void) const
Returns the total number of faces.
virtual int GetTotalSidesNumber(void) const =0
Returns the total number of faces.
std::pair< int, int > neighbors
Neighboring cells.
Definition: Edge.hpp:21
Voronoi tessellation with MPI option.
virtual vector< int > const & GetCellEdges(int index) const =0
Returns the indexes of a cell&#39;s edges.
BinLogger(const std::string &file_name)
Class constructor.
Diagnostic functions for Voronoi tessellation.
Vector2D GetMeshPoint(int index) const
Returns Position of mesh generating point.
double x
Component in the x direction.
Definition: geometry.hpp:45
VoronoiLogger(void)
Class constructor.