6 #ifndef VORONOIMESH_HPP 7 #define VORONOIMESH_HPP 1 10 #define _USE_MATH_DEFINES 16 #include "../misc/utils.hpp" 17 #include "../newtonian/two_dimensional/RefineStrategy.hpp" 18 #include "../newtonian/two_dimensional/geometric_outer_boundaries/SquareBox.hpp" 21 #include "../misc/int2str.hpp" 73 (vector<Vector2D>
const& points,
87 vector<Vector2D>
const& points,
108 vector<int>
Update(
const vector<Vector2D>& points,
bool reorder=
false);
164 void output(
string filename);
183 vector<vector<vector<int> > > &BoundaryNeigh);
221 vector<Vector2D> &points,
OuterBoundary const* obc, vector<size_t> &selfindex,
222 vector<int> &sentproc, vector<vector<int> > &sentpoints);
227 vector<Edge> cell_edges;
230 vector<vector<int> > mesh_vertices;
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;
241 bool legal_edge(
Edge *e);
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,
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);
289 bool PointInCell(std::vector<std::pair<Vector2D, Vector2D> >
const& cpoints,
Vector2D const& vec);
290 #endif // VORONOIMESH_HPP Voronoi tessellation class.
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.
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.
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'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 'ordered' 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.
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.