VoronoiMesh.hpp
Go to the documentation of this file.
1 
6 #ifndef VORONOIMESH_HPP
7 #define VORONOIMESH_HPP 1
8 
9 #ifdef _MSC_VER
10 #define _USE_MATH_DEFINES
11 #endif // _MSC_VER
12 
13 #include "Delaunay.hpp"
14 #include <list>
15 #include "tessellation.hpp"
16 #include "../misc/utils.hpp"
17 #include "../newtonian/two_dimensional/RefineStrategy.hpp"
18 #include "../newtonian/two_dimensional/geometric_outer_boundaries/SquareBox.hpp"
19 #include "voronoi_logger.hpp"
20 #include "ConvexHull.hpp"
21 #include "../misc/int2str.hpp"
22 #ifdef RICH_MPI
23 #include <mpi.h>
24 #endif
25 
26 using std::string;
27 
29 class VoronoiMesh : public Tessellation
30 {
31 public:
32 
34  Vector2D rR,Vector2D f)const;
35 
36  VoronoiMesh* clone(void)const;
37 
38  vector<int> GetNeighbors(int index)const;
39 
40  void GetNeighbors(int index, vector<int> &neigh)const;
46  vector<int> GetLiteralNeighbors(int index)const;
47 
48  int GetOriginalIndex(int point) const;
49 
50 #ifdef RICH_MPI
51 
57  void Initialise(vector<Vector2D> const& points,Tessellation const& vproc,
58  OuterBoundary const* outer,bool reorder=true);
59 #endif
60 
61  void Initialise(vector<Vector2D> const& points,OuterBoundary const* bc, bool reorder=true);
62 
64  VoronoiMesh(void);
65 
73  (vector<Vector2D> const& points,
74  OuterBoundary const& bc,
75  bool HOrder=true);
76 
77 #ifdef RICH_MPI
78 
86  (Tessellation const& proctess,
87  vector<Vector2D> const& points,
88  OuterBoundary const& bc,
89  bool HOrder=true);
90 #endif // RICH_MPI
91 
96  VoronoiMesh(VoronoiMesh const& other);
97 
98  #ifdef RICH_MPI
99 
105  vector<int> Update(const vector<Vector2D>& points,const Tessellation& vproc, bool reorder=false);
106 #endif // RICH_MPI
107 
108  vector<int> Update(const vector<Vector2D>& points, bool reorder=false);
109 
110  ~VoronoiMesh(void);
114  int GetPointNo(void) const;
115 
119  void SetPointNo(int N);
120 
125  Vector2D GetMeshPoint(int index) const;
126 
130  int GetTotalSidesNumber(void) const;
131 
136  Edge const& GetEdge(int index) const;
137 
142  double GetWidth(int index) const;
143 
148  double GetVolume(int index) const;
149 
151  vector<int>const& GetCellEdges(int index) const;
152 
157  Vector2D const& GetCellCM(int index) const;
158 
159  vector<Vector2D>& GetMeshPoints(void);
160 
164  void output(string filename);
165 
166  bool NearBoundary(int index) const;
167 
170 
174  vector<Edge>& GetAllEdges(void);
175 
182  void FindBoundaryRemoveSend(vector<int> const& ToRemove,vector<vector<int> > &BoundaryRemove,
183  vector<vector<vector<int> > > &BoundaryNeigh);
184 
185  vector<vector<int> >& GetDuplicatedPoints(void);
186 
187  vector<vector<int> >const& GetDuplicatedPoints(void)const;
188 
189  vector<int> GetDuplicatedProcs(void)const;
190 
191  vector<vector<int> >const& GetSentPoints(void)const;
192 
193  vector<int> GetSentProcs(void)const;
194 
195  vector<vector<int> >& GetGhostIndeces(void);
196 
197  vector<vector<int> >const& GetGhostIndeces(void)const;
198 
199  vector<size_t> GetSelfPoint(void)const;
200 
201  int GetTotalPointNumber(void)const;
202 
203  const vector<Edge>& getAllEdges(void) const;
204 
205  vector<Vector2D>& GetAllCM(void);
206 
207  void GetNeighborNeighbors(vector<int> &result, int point)const;
208 
209 #ifdef RICH_MPI
210 
220  vector<Vector2D> UpdateMPIPoints(Tessellation const& vproc, int rank,
221  vector<Vector2D> &points, OuterBoundary const* obc, vector<size_t> &selfindex,
222  vector<int> &sentproc, vector<vector<int> > &sentpoints);
223 #endif
224 private:
225  double eps;
226  OuterBoundary const* obc;
227  vector<Edge> cell_edges;
228  vector<Edge> edges;
229  vector<Vector2D> CM;
230  vector<vector<int> > mesh_vertices; // Which edges does each mesh point have
231  Delaunay Tri;
232  vector<int> GhostProcs;
233  vector<vector<int> > GhostPoints;
234  vector<int> SentProcs;
235  vector<vector<int> > SentPoints;
236  vector<size_t> selfindex;
237  vector<vector<int> > NGhostReceived;
238  vector<vector<int> > OrgCorner;
239  int Nextra;
240 
241  bool legal_edge(Edge *e);
242  void build_v(void);//Builds the voronoi mesh
243  VoronoiMesh& operator=(const VoronoiMesh& origin);
244  Vector2D CalcCellCM(size_t index) const;
245  void FindIntersectingPoints(vector<Edge> const& box_edges,
246  vector<vector<int> > &toduplicate);
247  vector<int> CellIntersectBoundary(vector<Edge> const&box_edges,int cell);
248  void GetAdditionalBoundary(vector<vector<int> > &copied,
249  vector<vector<int> > &result,vector<vector<int> > &totest);
250  void GetCorners(vector<vector<int> > &copied,vector<vector<int> > &result);
251  vector<int> AddPointsAlongEdge(size_t point,vector<vector<int> > const&copied,int side);
252  void GetRealNeighbor(vector<int> &result,int point)const;
253  vector<int> GetBorderingCells(vector<int> const& copied,
254  vector<int> const& totest,int tocheck,vector<int> tempresult,int outer);
255  bool CloseToBorder(int point,int &border);
256  void GetToTest(vector<vector<int> > &copied,vector<vector<int> > &totest);
257  void ConvexEdgeOrder(void);
258  vector<int> FindEdgeStartConvex(int point);
259  void SendRecv(vector<int> const& procorder,vector<int> const&
260  proclist,vector<vector<int> > &data);
261  void NonSendBoundary(vector<int> & proclist,vector<vector<int> > &
262  data,Tessellation const& v,vector<vector<int> > &totest,
263  vector<Edge> const& boxedges);
264  void NonSendCorners(vector<int> & proclist,vector<vector<int> > & data,
265  Tessellation const& v);
266  boost::array<double,4> FindMaxCellEdges(void);
267  vector<int> CellIntersectOuterBoundary(vector<Edge> const&box_edges,int cell);
268  void FindIntersectingOuterPoints(vector<Edge> const&bedge,vector<vector<int> >
269  &boxduplicate,vector<vector<int> > const&firstduplicated);
270  void SendRecvRemove(vector<int> const& procorder,vector<int> const&
271  proclist,vector<vector<int> > &data);
272  void GetNeighborNeighborsMPI(vector<int> &result,int point);
273 
274 
275 };
281 bool PointInCell(vector<Vector2D> const& cpoints,Vector2D const& vec);
282 
289 bool PointInCell(std::vector<std::pair<Vector2D, Vector2D> > const& cpoints, Vector2D const& vec);
290 #endif // VORONOIMESH_HPP
Voronoi tessellation class.
Definition: VoronoiMesh.hpp:29
vector< int > GetSentProcs(void) const
Returns the indeces of the processors with whom points where exchanged.
int GetOriginalIndex(int point) const
Returns the original index of the duplicated point.
Edge const & GetEdge(int index) const
Returns edge (interface between cells)
vector< vector< int > > & GetDuplicatedPoints(void)
Returns the indeces of the points that where sent to other processors as ghost points.
int GetTotalPointNumber(void) const
Returns the total number of points (including ghost)
vector< int > Update(const vector< Vector2D > &points, const Tessellation &vproc, bool reorder=false)
Updates the tessellation.
Abstract class for tessellation.
vector< size_t > GetSelfPoint(void) const
Returns the indeces of the points that remain with the processor after the ne processor mesh is built...
vector< Vector2D > UpdateMPIPoints(Tessellation const &vproc, int rank, vector< Vector2D > &points, OuterBoundary const *obc, vector< size_t > &selfindex, vector< int > &sentproc, vector< vector< int > > &sentpoints)
Communicate position of mesh generating points between processes.
bool PointInCell(vector< Vector2D > const &cpoints, Vector2D const &vec)
Checks if a point is inside a Voronoi cell.
void FindBoundaryRemoveSend(vector< int > const &ToRemove, vector< vector< int > > &BoundaryRemove, vector< vector< vector< int > > > &BoundaryNeigh)
Calculates and send the ghost points that are needed for AMR with MPI.
vector< vector< int > > const & GetSentPoints(void) const
Returns the indeces of the points that where sent to other processors.
Abstract class for voronoi tessellation.
vector< int > GetNeighbors(int index) const
Returns the indeces of the neighbors.
Abstract class for the tessellation.
Interface between two cells.
Definition: Edge.hpp:13
void Initialise(vector< Vector2D > const &points, Tessellation const &vproc, OuterBoundary const *outer, bool reorder=true)
Initialise the tessellation.
vector< int > GetLiteralNeighbors(int index) const
Returns the list of neighbors including ghost points.
vector< int > GetDuplicatedProcs(void) const
Returns the indeces of the processors with whom ghost points where exchanged.
A debugging function for the Voronoi tessellation.
void output(string filename)
Outputs the grid data.
double GetWidth(int index) const
Returns the effective width of a cell.
The Delaunay data structure. Gets a set of points and constructs the Delaunay tessellation.
Definition: Delaunay.hpp:30
bool NearBoundary(int index) const
Returns if the cell is adjacent to a boundary.
void GetNeighborNeighbors(vector< int > &result, int point) const
Retrieves vicarious neighbors.
vector< Edge > & GetAllEdges(void)
Returns a reference to a list of all edges.
int GetPointNo(void) const
Get Total number of mesh generating points.
Vector2D const & GetCellCM(int index) const
Returns Position of Cell&#39;s CM.
vector< int > const & GetCellEdges(int index) const
Returns a reference to a vector<int> containing the indexes of the edges related to a cell...
const vector< Edge > & getAllEdges(void) const
Returns reference to the list of all edges.
void reorder(vector< T > &v, vector< size_t > const &order)
Reorder a vector according to an index vector (obtained from the &#39;ordered&#39; function) ...
int GetTotalSidesNumber(void) const
Returns the total number of faces.
Calculates the convex hull from a set of points.
VoronoiMesh * clone(void) const
Cloning function.
void SetPointNo(int N)
Set the number of points.
VoronoiMesh(void)
Class default constructor.
vector< vector< int > > & GetGhostIndeces(void)
Returns the indeces of each ghost point in the vector of points that the tessellation holds...
double GetVolume(int index) const
Returns the volume of a cell.
vector< Vector2D > & GetMeshPoints(void)
Returns a reference to the point vector.
Abstract class for geometric boundary conditions for the tessellation.
voronoi_loggers::VoronoiLogger * logger
Diagnostics method.
Delaunay triangulation with mpi.
2D Mathematical vector
Definition: geometry.hpp:15
Vector2D GetMeshPoint(int index) const
Returns Position of mesh generating point.
Vector2D CalcFaceVelocity(Vector2D wl, Vector2D wr, Vector2D rL, Vector2D rR, Vector2D f) const
Calculates the velocity of a single edge.
vector< Vector2D > & GetAllCM(void)
Returns the center of masses of the cells.