HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
Polytope2D.hpp
Go to the documentation of this file.
1 // standard libraries
2 #include <vector> // std::vector
3 #include <algorithm> // std::find
4 #include <fstream> // std::ofstream
5 
6 // external linking to Eigen is required
7 #include <Eigen/Dense> // Eigen::Matrix, Eigen::MatrixXd, Eigen::FullPivLU
8 
9 // path to Math specified so that external linking is NOT required
10 #include <../Math/math.hpp> // Math::factorial, Math::sgn
11 
12 #ifndef _POLYTOPE2D_HPP
13 #define _POLYTOPE2D_HPP
14 
15 namespace Mesh2D
16 {
17  inline constexpr int DIMENSION = 2;
18 
19  using VectorRd = Eigen::Matrix<double, DIMENSION, 1>;
20  using VectorZd = Eigen::Matrix<int, DIMENSION, 1>;
21 
22  template <size_t object_dim>
23  using Simplex = std::array<VectorRd, object_dim + 1>;
24 
25  template <size_t object_dim>
26  using Simplices = std::vector<Simplex<object_dim>>;
27 
29  template <size_t object_dim>
31 
33  template <size_t object_dim>
34  double simplex_measure(Simplex<object_dim> simplex);
35 
36 
43  // ----------------------------------------------------------------------------
44  // Polytope class definition
45  // ----------------------------------------------------------------------------
46 
52  template <size_t object_dim>
53  class Polytope
54  {
55  public:
57  Polytope(size_t index, Simplices<object_dim> simplices);
58 
60  Polytope(size_t index, Simplex<object_dim> simplex);
61 
63  Polytope(size_t index, VectorRd vertex);
64 
66  Polytope();
67 
69  ~Polytope();
70 
71  // inline size_t dim() const { return object_dim; }
72 
73  inline size_t global_index() const { return _index; }
74  inline double diam() const { return _diameter; }
75  inline VectorRd center_mass() const { return _center_mass; }
76  inline double measure() const { return _measure; }
77  inline Simplices<object_dim> get_simplices() const { return _simplices; }
78  inline void set_global_index(const size_t idx) { _index = idx; }
79  inline bool is_boundary() const { return _is_boundary; }
80  inline void set_boundary(bool val) { _is_boundary = val; }
81 
82  inline std::vector<Polytope<0> *> get_vertices() const { return _vertices; }
83  inline std::vector<Polytope<1> *> get_edges() const { return _edges; }
84  inline std::vector<Polytope<DIMENSION - 1> *> get_faces() const { return _edges; }
85  inline std::vector<Polytope<DIMENSION> *> get_cells() const { return _cells; }
86 
87  inline size_t n_vertices() const { return _vertices.size(); }
88  inline size_t n_edges() const { return _edges.size(); }
89  inline size_t n_faces() const { return _edges.size(); }
90  inline size_t n_cells() const { return _cells.size(); }
91 
92  Polytope<0> *vertex(const size_t i) const;
93  Polytope<1> *edge(const size_t i) const;
94  Polytope<DIMENSION - 1> *face(const size_t i) const;
95  Polytope<DIMENSION> *cell(const size_t i) const;
96 
98  void add_edge(Polytope<1> *edge);
101 
102  int index_vertex(const Polytope<0> *vertex) const;
103  int index_edge(const Polytope<1> *edge) const;
104  int index_face(const Polytope<DIMENSION - 1> *face) const;
105  int index_cell(const Polytope<DIMENSION> *cell) const;
106 
107  VectorRd coords() const;
108 
109  VectorRd face_normal(const size_t face_index) const;
110  VectorRd edge_normal(const size_t edge_index) const;
111 
112  int face_orientation(const size_t face_index) const;
113  int edge_orientation(const size_t edge_index) const;
114  int vertex_orientation(const size_t vertex_index) const;
115 
116  VectorRd normal() const;
117  VectorRd tangent() const;
118 
119  void construct_face_normals();
120 
121  void plot_simplices(std::ofstream *out) const;
122  void plot(std::ofstream *out) const;
123 
124  private:
125  size_t _index;
126  VectorRd _center_mass;
127  double _measure;
128  double _diameter;
129  bool _is_boundary;
130  Simplices<object_dim> _simplices;
131  VectorRd _normal; // uninitialised unless object_dim == DIMENSION - 1 (face)
132  VectorRd _tangent; // uninitialised unless object_dim == 1 (edge)
133  std::vector<int> _face_directions; // empty unless object_dim == DIMENSION (cell)
134 
135  std::vector<Polytope<0> *> _vertices;
136  std::vector<Polytope<1> *> _edges;
137  // std::vector<Polytope<DIMENSION - 1>*> _faces;
138  std::vector<Polytope<DIMENSION> *> _cells;
139  };
140 
143 
145  using Edge = Polytope<1>;
146 
148  using Face = Polytope<DIMENSION - 1>;
149 
152 
154 
155  // -----------------------------------------------------------------
156  // ----------- Implementation of templated functions ---------------
157  // -----------------------------------------------------------------
158 
159  template <size_t object_dim>
161  {
162  VectorRd center_mass = VectorRd::Zero();
163 
164  for (auto &coord : simplex)
165  {
166  for (size_t i = 0; i < DIMENSION; ++i)
167  {
168  center_mass(i) += coord(i);
169  }
170  }
171  return center_mass / (object_dim + 1);
172  }
173 
174  template <size_t object_dim>
176  {
177  Eigen::MatrixXd CM_mat = Eigen::MatrixXd::Zero(object_dim + 2, object_dim + 2);
178  for (size_t i = 0; i < object_dim + 1; ++i)
179  {
180  VectorRd i_coords = simplex[i];
181  for (size_t j = i + 1; j < object_dim + 1; ++j)
182  {
183  VectorRd j_coords = simplex[j];
184  double norm = (j_coords - i_coords).norm();
185  CM_mat(i, j) = norm * norm;
186  CM_mat(j, i) = CM_mat(i, j);
187  }
188  CM_mat(i, object_dim + 1) = 1;
189  CM_mat(object_dim + 1, i) = 1;
190  }
191 
192  // Calculate Cayley-Menger determinant
193  double det = CM_mat.determinant();
194 
195  double scaling = std::pow(Math::factorial(object_dim), 2) * std::pow(2, object_dim);
196 
197  return std::sqrt(std::abs(det / scaling));
198  }
199 
200  template <size_t object_dim>
202  : _index(index), _simplices(simplices)
203  {
204  assert(object_dim <= DIMENSION);
205 
206  if (object_dim == 0 || object_dim == 1)
207  {
208  assert(simplices.size() == 1);
209  }
210 
211  _is_boundary = false; // by default, is_boundary is set to false. If polytope is on boundary, set is_boundary to true upon construction of the mesh.
212 
213  _measure = 0.0;
214  _center_mass = VectorRd::Zero();
215  std::vector<VectorRd> vertex_coords;
216  for (auto &simplex : simplices)
217  {
218  // assert(simplex.size() == _dim + 1);
219  double measure_of_simplex = simplex_measure<object_dim>(simplex);
220  _measure += measure_of_simplex;
221  _center_mass += measure_of_simplex * simplex_center_mass<object_dim>(simplex);
222  for (auto &coord : simplex)
223  {
224  if (std::find(vertex_coords.begin(), vertex_coords.end(), coord) == vertex_coords.end())
225  {
226  vertex_coords.push_back(coord); // if coord not already in vertex_coords, add it
227  }
228  }
229  }
230  _center_mass /= _measure;
231 
232  _diameter = 0.0;
233  for (auto it = vertex_coords.begin(); it != vertex_coords.end(); ++it)
234  {
235  for (auto jt = it; jt != vertex_coords.end(); ++jt)
236  {
237  _diameter = std::max(_diameter, (*it - *jt).norm());
238  }
239  }
240 
241  if (object_dim == DIMENSION - 1) // find the tangent and normal
242  {
243  _tangent = (this->get_simplices()[0][1] - this->get_simplices()[0][0]).normalized();
244  _normal = VectorRd(-_tangent(1), _tangent(0));
245  }
246  }
247 
248  template <size_t object_dim>
249  Polytope<object_dim>::Polytope(size_t index, Simplex<object_dim> simplex) // simplex
250  : Polytope(index, Simplices<object_dim>(1, simplex))
251  {
252  }
253 
254  template <size_t object_dim>
255  Polytope<object_dim>::Polytope(size_t index, VectorRd vertex) // vertex
256  : Polytope(index, Simplex<object_dim>({vertex}))
257  {
258  assert(object_dim == 0); // must be a zero dimensional object to use vertex constructor
259  }
260 
261  template <size_t object_dim>
263 
264  template <size_t object_dim>
266 
267  template <size_t object_dim>
269  {
270  assert(std::find(_vertices.begin(), _vertices.end(), vertex) == _vertices.end()); // ensure vertex does not already exist in _vertices
271  _vertices.push_back(vertex);
272  }
273 
274  template <size_t object_dim>
276  {
277  assert(std::find(_edges.begin(), _edges.end(), edge) == _edges.end());
278  _edges.push_back(edge);
279  }
280 
281  template <size_t object_dim>
283  {
284  assert(std::find(_edges.begin(), _edges.end(), face) == _edges.end());
285  _edges.push_back(face);
286  }
287 
288  template <size_t object_dim>
290  {
291  assert(std::find(_cells.begin(), _cells.end(), cell) == _cells.end());
292  _cells.push_back(cell);
293  }
294 
295  template <size_t object_dim>
297  {
298  assert(i < _vertices.size());
299  return _vertices[i];
300  }
301 
302  template <size_t object_dim>
304  {
305  assert(i < _edges.size());
306  return _edges[i];
307  }
308 
309  template <size_t object_dim>
310  Polytope<DIMENSION - 1> *Polytope<object_dim>::face(const size_t i) const
311  {
312  assert(i < _edges.size());
313  return _edges[i];
314  }
315 
316  template <size_t object_dim>
318  {
319  assert(i < _cells.size());
320  return _cells[i];
321  }
322 
323  template <size_t object_dim>
324  void Polytope<object_dim>::construct_face_normals() // not very efficient - probably room for improvement
325  {
326  assert(object_dim == DIMENSION);
327  for (size_t iF = 0; iF < _edges.size(); ++iF)
328  {
329  VectorRd normal = _edges[iF]->normal();
330  Simplex<object_dim - 1> face_simplex = _edges[iF]->get_simplices()[0];
331  VectorRd center = VectorRd::Zero();
332  double count;
333  for (auto &cell_simplex : this->get_simplices())
334  {
335  count = 0;
336  for (size_t i = 0; (i < cell_simplex.size()) && count < 2; ++i)
337  {
338  if (std::find(face_simplex.begin(), face_simplex.end(), cell_simplex[i]) == face_simplex.end()) // requires numerical precision
339  {
340  ++count;
341  }
342  }
343  if (count == 1) // only don't share one coordinate
344  {
345  center = simplex_center_mass<object_dim>(cell_simplex);
346  break;
347  }
348  }
349  assert(count == 1);
350  // _face_directions.push_back(Math::sgn((_faces[iF]->center_mass() - _center_mass).dot(normal))); // star shaped wrt center mass
351  _face_directions.push_back(Math::sgn((_edges[iF]->center_mass() - center).dot(normal)));
352  assert(_face_directions[iF] != 0);
353  }
354  }
355 
356  template <size_t object_dim>
357  VectorRd Polytope<object_dim>::face_normal(const size_t face_index) const
358  {
359  assert(object_dim == DIMENSION);
360  assert(face_index < _edges.size());
361  return _face_directions[face_index] * _edges[face_index]->normal();
362  }
363 
364  template <size_t object_dim>
365  VectorRd Polytope<object_dim>::edge_normal(const size_t edge_index) const
366  {
367  return this->face_normal(edge_index);
368  }
369 
370  template <size_t object_dim>
371  VectorRd Polytope<object_dim>::coords() const // only for vertices
372  {
373  assert(object_dim == 0); // can only return coordinate of a vertex
374  return this->get_simplices()[0][0]; // only has one simplex, with one coordinate
375  }
376 
377  template <size_t object_dim>
379  {
380  assert(object_dim == DIMENSION - 1);
381  return _normal;
382  }
383 
384  template <size_t object_dim>
386  {
387  assert(object_dim == 1);
388  return _tangent;
389  // return (this->get_vertices()[1]->coords() - this->get_vertices()[0]->coords()).normalized();
390  }
391 
392  template <size_t object_dim>
394  {
395  auto itr = std::find(_vertices.begin(), _vertices.end(), vertex);
396  if (itr != _vertices.end())
397  {
398  return itr - _vertices.begin();
399  }
400  else
401  {
402  throw "Vertex not found";
403  }
404  }
405 
406  template <size_t object_dim>
408  {
409  auto itr = std::find(_edges.begin(), _edges.end(), edge);
410  if (itr != _edges.end())
411  {
412  return itr - _edges.begin();
413  }
414  else
415  {
416  throw "Edge not found";
417  }
418  }
419 
420  template <size_t object_dim>
422  {
423  auto itr = std::find(_edges.begin(), _edges.end(), face);
424  if (itr != _edges.end())
425  {
426  return itr - _edges.begin();
427  }
428  else
429  {
430  throw "Face not found";
431  }
432  }
433 
434  template <size_t object_dim>
436  {
437  auto itr = std::find(_cells.begin(), _cells.end(), cell);
438  if (itr != _cells.end())
439  {
440  return itr - _cells.begin();
441  }
442  else
443  {
444  throw "Cell not found";
445  }
446  }
447 
448  template <size_t object_dim>
449  int Polytope<object_dim>::face_orientation(const size_t face_index) const
450  {
451  assert(object_dim == DIMENSION);
452  assert(face_index < _edges.size());
453  return _face_directions[face_index];
454  }
455 
456  template <size_t object_dim>
457  int Polytope<object_dim>::edge_orientation(const size_t edge_index) const
458  {
459  assert(object_dim == 2);
460  assert(edge_index < _edges.size());
461  return _face_directions[edge_index];
462  }
463 
464  template <size_t object_dim>
465  int Polytope<object_dim>::vertex_orientation(const size_t vertex_index) const
466  {
467  assert(object_dim == 1);
468  assert(vertex_index < _vertices.size());
469  return ( (_vertices[vertex_index]->coords()-_center_mass).dot(this->tangent()) > 0 ? 1 : -1);
470  }
471 
472  template <size_t object_dim>
473  void Polytope<object_dim>::plot_simplices(std::ofstream *out) const
474  {
475  assert(object_dim > 0);
476  for (auto &simplex : _simplices)
477  {
478  for (size_t i = 0; i < object_dim + 1; ++i)
479  {
480  size_t i_next = (i + 1) % (object_dim + 1);
481  *out << simplex[i](0) << " " << simplex[i](1) << std::endl;
482  *out << simplex[i_next](0) << " " << simplex[i_next](1) << std::endl;
483  *out << std::endl;
484  }
485  }
486  }
487 
488  template <size_t object_dim>
489  void Polytope<object_dim>::plot(std::ofstream *out) const
490  {
491  assert(object_dim > 0);
492  for (size_t i = 0; i < _vertices.size(); ++i)
493  {
494  size_t i_next = (i + 1) % _vertices.size();
495  *out << _vertices[i]->coords()(0) << " " << _vertices[i]->coords()(1) << std::endl;
496  *out << _vertices[i_next]->coords()(0) << " " << _vertices[i_next]->coords()(1) << std::endl;
497  *out << std::endl;
498  }
499  }
500 
501 } // namespace Mesh2D
502 
503 #endif
A Vertex is a Polytope with object_dim = 0.
Definition: Polytope2D.hpp:54
Compute max and min eigenvalues of all matrices for i
Definition: compute_eigs.m:5
VectorRd normal() const
Return the normal of a Face.
Definition: Polytope2D.hpp:378
double measure() const
Return the Lebesgue measure of the Polytope.
Definition: Polytope2D.hpp:76
void construct_face_normals()
Set the directions of the face normals of a cell.
Definition: Polytope2D.hpp:324
int index_cell(const Polytope< DIMENSION > *cell) const
Returns the local index of a cell.
Definition: Polytope2D.hpp:435
int index_edge(const Polytope< 1 > *edge) const
Returns the local index of an edge.
Definition: Polytope2D.hpp:407
void add_edge(Polytope< 1 > *edge)
Add an edge to the Polytope.
Definition: Polytope2D.hpp:275
VectorRd center_mass() const
Return the center mass of the Polytope.
Definition: Polytope2D.hpp:75
void plot(std::ofstream *out) const
Plot the polytope to out.
Definition: Polytope2D.hpp:489
std::vector< Polytope< 0 > * > get_vertices() const
Return the vertices of the Polytope.
Definition: Polytope2D.hpp:82
VectorRd tangent() const
Return the tangent of a Edge.
Definition: Polytope2D.hpp:385
void set_global_index(const size_t idx)
Set the global index.
Definition: Polytope2D.hpp:78
size_t n_vertices() const
Return the number of vertices of the Polytope.
Definition: Polytope2D.hpp:87
void plot_simplices(std::ofstream *out) const
Plot the simplices to out.
Definition: Polytope2D.hpp:473
VectorRd face_normal(const size_t face_index) const
Return the outer normal of a Cell towards the Face located at face_index.
Definition: Polytope2D.hpp:357
size_t n_edges() const
Return the number of edges of the Polytope.
Definition: Polytope2D.hpp:88
int index_face(const Polytope< DIMENSION - 1 > *face) const
Returns the local index of a face.
Definition: Polytope2D.hpp:421
void add_vertex(Polytope< 0 > *vertex)
Add a vertex to the Polytope.
Definition: Polytope2D.hpp:268
void add_cell(Polytope< DIMENSION > *cell)
Add a cell to the Polytope.
Definition: Polytope2D.hpp:289
int index_vertex(const Polytope< 0 > *vertex) const
Returns the local index of a vertex.
Definition: Polytope2D.hpp:393
std::vector< Polytope< DIMENSION - 1 > * > get_faces() const
Return the faces of the Polytope.
Definition: Polytope2D.hpp:84
Polytope< DIMENSION > * cell(const size_t i) const
Return the i-th cell of the Polytope.
Definition: Polytope2D.hpp:317
int vertex_orientation(const size_t vertex_index) const
Return the orientation of a Vertex.
Definition: Polytope2D.hpp:465
double diam() const
Return the diameter of the Polytope.
Definition: Polytope2D.hpp:74
~Polytope()
Definition: Polytope2D.hpp:265
Polytope< DIMENSION - 1 > * face(const size_t i) const
Return the i-th face of the Polytope.
Definition: Polytope2D.hpp:310
void add_face(Polytope< DIMENSION - 1 > *face)
Add a face to the Polytope.
Definition: Polytope2D.hpp:282
void set_boundary(bool val)
Set the boundary value of the Polytope.
Definition: Polytope2D.hpp:80
size_t n_faces() const
Return the number of faces of the Polytope.
Definition: Polytope2D.hpp:89
Polytope< 1 > * edge(const size_t i) const
Return the i-th edge of the Polytope.
Definition: Polytope2D.hpp:303
std::vector< Polytope< 1 > * > get_edges() const
Return the edges of the Polytope.
Definition: Polytope2D.hpp:83
bool is_boundary() const
Return true if Polytope is a boundary object, false otherwise.
Definition: Polytope2D.hpp:79
VectorRd edge_normal(const size_t edge_index) const
Return the edge normal of a 2D object.
Definition: Polytope2D.hpp:365
double simplex_measure(Simplex< object_dim > simplex)
Definition: Polytope2D.hpp:175
size_t global_index() const
Return the global index of the Polytope.
Definition: Polytope2D.hpp:73
int face_orientation(const size_t face_index) const
Return the orientation of a Face.
Definition: Polytope2D.hpp:449
Polytope()
Destructor.
Definition: Polytope2D.hpp:262
Simplices< object_dim > get_simplices() const
Return the simplices making up the Polytope.
Definition: Polytope2D.hpp:77
int edge_orientation(const size_t edge_index) const
Return the orientation of a Edge.
Definition: Polytope2D.hpp:457
Polytope< 0 > * vertex(const size_t i) const
Return the i-th vertex of the Polytope.
Definition: Polytope2D.hpp:296
VectorRd simplex_center_mass(Simplex< object_dim > simplex)
Method to find the Lebesgue measure of an arbitrary simplex in arbitrary space.
Definition: Polytope2D.hpp:160
VectorRd coords() const
Return the coordinates of a Vertex.
Definition: Polytope2D.hpp:371
std::vector< Polytope< DIMENSION > * > get_cells() const
Return the cells of the Polytope.
Definition: Polytope2D.hpp:85
size_t n_cells() const
Return the number of cells of the Polytope.
Definition: Polytope2D.hpp:90
end Read size then branch according to sparse or dense format count
Definition: mmread.m:75
for j
Definition: mmread.m:174
const int sgn(T val)
Definition: math.hpp:42
const size_t factorial(size_t n)
Definition: math.hpp:12
Definition: Mesh2D.hpp:7
std::array< VectorRd, object_dim+1 > Simplex
Definition: Polytope2D.hpp:23
constexpr int DIMENSION
Definition: Polytope2D.hpp:17
std::vector< Simplex< object_dim > > Simplices
Method to find the center mass of an arbitrary simplex in arbitrary space.
Definition: Polytope2D.hpp:28
Eigen::Matrix< int, DIMENSION, 1 > VectorZd
Definition: Polytope2D.hpp:20
Eigen::Matrix< double, DIMENSION, 1 > VectorRd
Definition: Polytope2D.hpp:19