Delaunay.hpp
Go to the documentation of this file.
1 
6 #ifndef DELAUNAY_HPP
7 #define DELAUNAY_HPP 1
8 #include "facet.hpp"
9 #include <algorithm>
10 #include <vector>
11 #include <stack>
12 #include <climits>
13 #include "geometry.hpp"
14 #include "../misc/universal_error.hpp"
15 #include "HilbertOrder.hpp"
16 #include "geotests.hpp"
17 #include "../misc/utils.hpp"
18 #include "delaunay_logger.hpp"
19 #include "Edge.hpp"
20 #include "../newtonian/two_dimensional/OuterBoundary.hpp"
21 #include "shape_2d.hpp"
22 #include "../misc/int2str.hpp"
23 #include "../newtonian/two_dimensional/diagnostics.hpp"
24 #ifdef RICH_MPI
25 #include "find_affected_cells.hpp"
26 #endif
27 
30 class Delaunay
31 {
32 private:
33 #ifdef RICH_MPI
34  int findSomeOuterPoint(void);
35 
36  pair<vector<vector<int> >, vector<vector<int> > > findOuterPoints(const Tessellation& t_proc,
37  const vector<Edge>& edge_list, const vector<Edge>& box_edges, vector<vector<int> > &NghostIndex, bool periodic,
38  std::vector<int> &cpu_neigh, std::vector<Vector2D> &periodic_add_self, std::vector<std::vector<Vector2D> >
39  &periodic_add_others, Vector2D const &ll, Vector2D const &ur);
40 
41  std::pair<vector<vector<int> >, vector<int> > FindOuterPoints2
42  (const Tessellation& t_proc,
43  const vector<Edge>& edge_list,
44  vector<vector<int> > &to_duplicate,
45  vector<vector<int> >& self_points,
46  const vector<Edge>& box_edges,
47  vector<vector<int> > &NghostIndex,
48  bool periodic,
49  std::vector<int> &cpu_neigh, std::vector<Vector2D> &periodic_add_self, std::vector<std::vector<Vector2D> >
50  &periodic_add_others, Vector2D const &ll, Vector2D const &ur);
51 
52  vector<vector<int> > boundary_intersection_check
53  (const vector<Edge>& edges,
54  const vector<vector<int> >& to_duplicate);
55 #endif // RICH_MPI
56 
57  enum Sides { RIGHT, UP, LEFT, DOWN, LU, LD, RU, RD };
58  int lastFacet; //last facet to be checked in Walk
59  bool CalcRadius;
60 
61  class DataOnlyForBuild
62  {
63  public:
64  DataOnlyForBuild();
65  DataOnlyForBuild(DataOnlyForBuild const& other);
66  DataOnlyForBuild& operator=(DataOnlyForBuild const& other);
67  vector<vector<char> > copied;
68  };
69 
70  vector<double> radius;
71  bool PointWasAdded;
72  int last_facet_added;
73  vector<facet> f;
74  vector<Vector2D> cor;
75  int length;
76  size_t olength;
77  int location_pointer;
78  int last_loc;
79  vector<int> OrgIndex;
80 
81  bool IsOuterFacet(int facet)const;
82  void add_point(size_t index, stack<std::pair<size_t, size_t> > &flip_stack);
83  void flip(size_t i, size_t j, stack<std::pair<size_t, size_t> > & flip_stack);
84  size_t Walk(size_t point);
85  double CalculateRadius(int facet);
86  double CalcRadiusHiRes(int facet);
87  int FindPointInFacet(int facet, int point);
88  double FindMaxRadius(int point);
89  void FindContainingTetras(int StartTetra, int point, vector<int> &tetras);
90  vector<int> FindContainingTetras(int StartTetra, int point);
91  vector<vector<int> > FindOuterPoints(vector<Edge> const& edges);
92  void AddPeriodicMPI(std::vector<int>& toduplicate, std::vector<Vector2D>& periodic_add);
93  vector<vector<int> > FindOuterPointsMPI(OuterBoundary const* obc, vector<Edge> const& edges,
94  Tessellation const& tproc, vector<vector<int> > &Nghost, vector<int> &proclist);
95  bool IsTripleOut(int index) const;
96  int FindTripleLoc(facet const& f)const;
97  void AddFacetDuplicate(int index, vector<vector<int> > &toduplicate, vector<Edge>
98  const& edges, vector<bool> &checked)const;
99  void AddOuterFacets(int tri, vector<vector<int> > &toduplicate, vector<Edge>
100  const& edges, vector<bool> &checked);
101 
102  vector<vector<int> > AddOuterFacetsMPI(int point,vector<vector<int> > &toduplicate,vector<int> &neigh,
103  vector<bool> &checked,const Tessellation& tproc,const vector<Edge>& own_edges,bool periodic,
104  std::vector<Vector2D> &periodic_add_self,std::vector<std::vector<Vector2D> >& periodic_add_others,
105  Vector2D const &ll, Vector2D const &ur, bool recursive = false);
106 
107  void AddRigid(vector<Edge> const& edges,vector<vector<int> > &toduplicate);
108  vector<vector<int> > AddPeriodic(OuterBoundary const* obc, vector<Edge> const& edges,
109  vector<vector<int> > &toduplicate);
110  void AddHalfPeriodic(OuterBoundary const* obc, vector<Edge> const& edges, vector<vector<int> > &toduplicate);
111  double GetMaxRadius(int point, int startfacet);
112  void SendRecvFirstBatch(vector<vector<Vector2D> > &tosend,vector<int> const& neigh, vector<vector<int> > &Nghost);
113  vector<int> GetOuterFacets(int cur_facet, int real_point, int olength);
114 
115  Delaunay& operator=(const Delaunay& origin);
116 
117 public:
118 
123  int GetOrgIndex(int index)const;
124 
128  void ChangeOlength(int n);
129 
133  void Changelength(int n);
134 
138  vector<Vector2D>& ChangeCor(void);
139 
143  const vector<Vector2D>& getCor(void) const;
144 
146  Delaunay(void);
147 
152  Delaunay(Delaunay const& other);
153 
155  ~Delaunay(void);
156 
161  void build_delaunay(vector<Vector2D>const& vp, std::vector<std::pair<Vector2D, Vector2D> >const & cpoints);
162 
164  void output(void);
165 
170  const facet& get_facet(int index) const;
171 
178  double get_facet_coordinate(int Facet, int vertice, int dim);
179 
184  Vector2D get_point(size_t index) const;
185 
191  double get_cor(int index, int dim) const;
192 
196  int get_num_facet(void)const;
197 
201  int get_length(void) const;
202 
206  int get_last_loc(void) const;
207 
212  void set_point(int index, Vector2D p);
213 
218  double triangle_area(int index);
219 
225  void update(const vector<Vector2D>& points, std::vector<std::pair<Vector2D, Vector2D> >const& cpoints);
226 
232  int GetOriginalIndex(int NewPoint) const;
233 
238  int GetOriginalLength(void) const;
239 
244  vector<Vector2D>& GetMeshPoints(void);
245 
249  int GetTotalLength(void);
250 
255  double GetFacetRadius(int facet) const;
256 
261  void AddAditionalPoint(Vector2D const& vec);
266  int GetCorSize(void)const;
267 
274  Vector2D GetCircleCenter(int index, double &R)const;
275 
284  void BuildBoundary(OuterBoundary const* obc, vector<Edge> const& edges);
292  pair<vector<vector<int> >, vector<int> > BuildBoundary(OuterBoundary const* obc, Tessellation const& tproc,
293  vector<vector<int> > &Nghost);
298  void AddBoundaryPoints(vector<Vector2D> const& points);
299 
300  bool CheckCorrect(void);
301 };
302 #endif //DELAUNAYMPI_HPP
double get_cor(int index, int dim) const
Returns a coordinate.
Definition: Delaunay.cpp:709
void Changelength(int n)
Changes the cor length.
Definition: Delaunay.cpp:676
Geometrical calculations.
Various checks for geometric data.
void AddBoundaryPoints(vector< Vector2D > const &points)
Adds the points to the tessellation, used for boundary points.
Definition: Delaunay.cpp:754
void AddAditionalPoint(Vector2D const &vec)
Adds a point to the cor vector. Used in periodic boundaries with AMR.
Definition: Delaunay.cpp:767
void set_point(int index, Vector2D p)
Change Mesh point.
Definition: Delaunay.cpp:734
Abstract class for tessellation.
const facet & get_facet(int index) const
Returns a facet.
Definition: Delaunay.cpp:691
Delanauy triangle data structure. Keeps the indexes of the vertices (right hand fashion) and the neig...
Definition: facet.hpp:19
void update(const vector< Vector2D > &points, std::vector< std::pair< Vector2D, Vector2D > >const &cpoints)
Updates the triangulation.
Definition: Delaunay.cpp:504
const vector< Vector2D > & getCor(void) const
Access to coordinates.
Definition: Delaunay.cpp:686
Edge between cells.
int get_length(void) const
Returns the number of points.
Definition: Delaunay.cpp:724
Hilbert space filling curve.
int GetTotalLength(void)
Returns the length of all the points (included duplicated)
Definition: Delaunay.cpp:749
void output(void)
Dumps the Delaunay tessellation into a binary file.
Vector2D GetCircleCenter(int index, double &R) const
Returns the center of the circumscribed circle of a facet.
Definition: Delaunay.cpp:1202
vector< Vector2D > & GetMeshPoints(void)
Returns a refrence to the points.
Definition: Delaunay.cpp:744
void BuildBoundary(OuterBoundary const *obc, vector< Edge > const &edges)
Builds the boundary points.
Definition: Delaunay.cpp:1167
~Delaunay(void)
Default destructor.
Definition: Delaunay.cpp:177
Two dimensional shapes.
int get_num_facet(void) const
Returns the number of facets.
Definition: Delaunay.cpp:719
Diagnostic class for Delaunay triangulation.
Triangulation data.
The Delaunay data structure. Gets a set of points and constructs the Delaunay tessellation.
Definition: Delaunay.hpp:30
int GetCorSize(void) const
Gets the size of the cor vector.
Definition: Delaunay.cpp:772
int GetOrgIndex(int index) const
Retrieves the original index of a point (in case a point was duplicated)
Definition: Delaunay.cpp:1963
Vector2D get_point(size_t index) const
Returns a point.
Definition: Delaunay.cpp:704
delaunay_loggers::DelaunayLogger * logger
Diagnostics.
Definition: Delaunay.hpp:277
double get_facet_coordinate(int Facet, int vertice, int dim)
Returns a coordinate of a vertice.
Definition: Delaunay.cpp:696
void ChangeOlength(int n)
Changes the cor olength.
Definition: Delaunay.cpp:671
int GetOriginalLength(void) const
Returns the original length of the points (without duplicated points)
Definition: Delaunay.cpp:739
int get_last_loc(void) const
Returns the last location, a number used to identify the fact that the neighbor of a facet is empty...
Definition: Delaunay.cpp:729
Abstract class for geometric boundary conditions for the tessellation.
A debugging function for the Delaunay triangulation.
Determines which cells might be affected by boundary conditions.
double triangle_area(int index)
Returns the area of the triangle. Negative result means the triangle isn&#39;t right handed.
Definition: Delaunay.cpp:491
2D Mathematical vector
Definition: geometry.hpp:15
void build_delaunay(vector< Vector2D >const &vp, std::vector< std::pair< Vector2D, Vector2D > >const &cpoints)
Builds the Delaunay tessellation.
Definition: Delaunay.cpp:376
double GetFacetRadius(int facet) const
return the facet&#39;s radius
Definition: Delaunay.cpp:666
Delaunay(void)
Class constructor.
Definition: Delaunay.cpp:148
vector< Vector2D > & ChangeCor(void)
Allows to change the cor.
Definition: Delaunay.cpp:681
int GetOriginalIndex(int NewPoint) const
Returns the original index of the duplicated point in Periodic Boundary conditions.
Definition: Delaunay.cpp:661