15 template <std::
size_t dimension>
22 for (
auto &
vertex : _vertices)
26 for (
auto &
edge : _edges)
32 for (
auto &
face : _faces)
37 for (
auto &
cell : _cells)
43 void set_name(std::string name) { _mesh_name = name; }
44 inline std::string
get_name() {
return _mesh_name; }
49 for (
auto &
cell : _cells)
51 val = std::max(val,
cell->diam());
55 inline std::size_t
dim()
const {
return dimension; }
57 inline std::size_t
n_vertices()
const {
return _vertices.size(); }
58 inline std::size_t
n_edges()
const {
return _edges.size(); }
59 inline std::size_t
n_faces()
const {
return _faces.size(); }
60 inline std::size_t
n_cells()
const {
return _cells.size(); }
61 inline std::size_t
n_elems(
size_t d)
const
63 assert(d <= dimension && dimension < 4);
75 inline std::size_t
n_b_vertices()
const {
return _b_vertices.size(); }
76 inline std::size_t
n_b_edges()
const {
return _b_edges.size(); }
77 inline std::size_t
n_b_faces()
const {
return _b_faces.size(); }
78 inline std::size_t
n_b_cells()
const {
return _b_cells.size(); }
80 inline std::size_t
n_i_vertices()
const {
return _i_vertices.size(); }
81 inline std::size_t
n_i_edges()
const {
return _i_edges.size(); }
82 inline std::size_t
n_i_faces()
const {
return _i_faces.size(); }
83 inline std::size_t
n_i_cells()
const {
return _i_cells.size(); }
85 inline std::vector<Vertex<dimension> *>
get_vertices()
const {
return _vertices; }
86 inline std::vector<Edge<dimension> *>
get_edges()
const {
return _edges; }
87 inline std::vector<Face<dimension> *>
get_faces()
const {
return _faces; }
88 inline std::vector<Cell<dimension> *>
get_cells()
const {
return _cells; }
90 inline std::vector<Vertex<dimension> *>
get_b_vertices()
const {
return _b_vertices; }
91 inline std::vector<Edge<dimension> *>
get_b_edges()
const {
return _b_edges; }
92 inline std::vector<Face<dimension> *>
get_b_faces()
const {
return _b_faces; }
93 inline std::vector<Cell<dimension> *>
get_b_cells()
const {
return _b_cells; }
95 inline std::vector<Vertex<dimension> *>
get_i_vertices()
const {
return _i_vertices; }
96 inline std::vector<Edge<dimension> *>
get_i_edges()
const {
return _i_edges; }
97 inline std::vector<Face<dimension> *>
get_i_faces()
const {
return _i_faces; }
98 inline std::vector<Cell<dimension> *>
get_i_cells()
const {
return _i_cells; }
102 assert(std::find(_vertices.begin(), _vertices.end(),
vertex) == _vertices.end());
103 _vertices.push_back(
vertex);
107 assert(std::find(_edges.begin(), _edges.end(),
edge) == _edges.end());
108 _edges.push_back(
edge);
111 assert(std::find(_faces.begin(), _faces.end(),
reinterpret_cast<Face<dimension> *
>(
edge)) == _faces.end());
117 assert(std::find(_faces.begin(), _faces.end(),
face) == _faces.end());
118 _faces.push_back(
face);
121 assert(std::find(_edges.begin(), _edges.end(),
reinterpret_cast<Edge<dimension> *
>(
face)) == _edges.end());
127 assert(std::find(_cells.begin(), _cells.end(),
cell) == _cells.end());
128 _cells.push_back(
cell);
133 assert(std::find(_b_vertices.begin(), _b_vertices.end(),
vertex) == _b_vertices.end());
134 _b_vertices.push_back(
vertex);
138 assert(std::find(_b_edges.begin(), _b_edges.end(),
edge) == _b_edges.end());
139 _b_edges.push_back(
edge);
142 assert(std::find(_b_faces.begin(), _b_faces.end(),
reinterpret_cast<Face<dimension> *
>(
edge)) == _b_faces.end());
148 assert(std::find(_b_faces.begin(), _b_faces.end(),
face) == _b_faces.end());
149 _b_faces.push_back(
face);
152 assert(std::find(_b_edges.begin(), _b_edges.end(),
reinterpret_cast<Edge<dimension> *
>(
face)) == _b_edges.end());
158 assert(std::find(_b_cells.begin(), _b_cells.end(),
cell) == _b_cells.end());
159 _b_cells.push_back(
cell);
164 assert(std::find(_i_vertices.begin(), _i_vertices.end(),
vertex) == _i_vertices.end());
165 _i_vertices.push_back(
vertex);
169 assert(std::find(_i_edges.begin(), _i_edges.end(),
edge) == _i_edges.end());
170 _i_edges.push_back(
edge);
173 assert(std::find(_i_faces.begin(), _i_faces.end(),
reinterpret_cast<Face<dimension> *
>(
edge)) == _i_faces.end());
179 assert(std::find(_i_faces.begin(), _i_faces.end(),
face) == _i_faces.end());
180 _i_faces.push_back(
face);
183 assert(std::find(_i_edges.begin(), _i_edges.end(),
reinterpret_cast<Edge<dimension> *
>(
face)) == _i_edges.end());
189 assert(std::find(_i_cells.begin(), _i_cells.end(),
cell) == _i_cells.end());
190 _i_cells.push_back(
cell);
195 inline Edge<dimension> *
edge(std::size_t index)
const { assert(index < _edges.size());
return _edges[index]; }
196 inline Face<dimension> *
face(std::size_t index)
const { assert(index < _faces.size());
return _faces[index]; }
197 inline Cell<dimension> *
cell(std::size_t index)
const { assert(index < _cells.size());
return _cells[index]; }
200 inline Edge<dimension> *
b_edge(std::size_t index)
const { assert(index < _b_edges.size());
return _b_edges[index]; }
201 inline Face<dimension> *
b_face(std::size_t index)
const { assert(index < _b_faces.size());
return _b_faces[index]; }
202 inline Cell<dimension> *
b_cell(std::size_t index)
const { assert(index < _b_cells.size());
return _b_cells[index]; }
205 inline Edge<dimension> *
i_edge(std::size_t index)
const { assert(index < _i_edges.size());
return _i_edges[index]; }
206 inline Face<dimension> *
i_face(std::size_t index)
const { assert(index < _i_faces.size());
return _i_faces[index]; }
207 inline Cell<dimension> *
i_cell(std::size_t index)
const { assert(index < _i_cells.size());
return _i_cells[index]; }
212 std::vector<size_t> rv;
213 static_assert(dimension < 4);
214 assert(d <= dimension);
215 if (d == dimension) {
216 assert(index < _cells.size());
219 for (
size_t j = 0; j < T.
n_faces(); ++j) {
220 rv.push_back(T.
face(j)->global_index());
223 assert(index < _faces.size());
226 for (
size_t j = 0; j < F.
n_edges(); ++j) {
227 rv.push_back(F.
edge(j)->global_index());
230 assert(index < _edges.size());
233 rv.push_back(E.
vertex(0)->global_index());
234 rv.push_back(E.
vertex(1)->global_index());
241 static_assert(dimension < 4);
242 assert(d <= dimension);
243 if (d == dimension) {
244 assert(i < _cells.size());
245 return _cells[i]->induced_orientation(j);
247 assert(i < _faces.size());
248 return _faces[i]->induced_orientation(j);
251 return (j%2==0?-1:1);
264 std::vector<std::vector<double>> reg_cell(
n_cells(), {0.0, 0.0});
265 std::size_t count = 0;
266 for (
auto &T : _cells)
268 double hT = T->diam();
271 reg_cell[count][0] = hT / pow(T->measure(), 1.0 / dimension);
274 std::vector<Face<dimension> *> faces = T->get_faces();
275 for (
auto &F : faces)
277 double hF = F->diam();
281 reg_cell[count][0] = std::max(reg_cell[count][0], hT / hF);
283 rhoT = std::min(rhoT, std::abs((xT - xF).dot(nTF)));
285 reg_cell[count][1] = hT / rhoT;
289 std::vector<double> value(2, 0.0);
290 for (
size_t iT = 0; iT <
n_cells(); iT++)
292 value[0] = std::max(value[0], reg_cell[iT][0]);
293 value[1] = std::max(value[1], reg_cell[iT][1]);
299 void renum(
const char B,
const std::vector<size_t> new_to_old)
306 std::vector<Cell<dimension> *> old_index = _cells;
307 for (
size_t i = 0; i < _cells.size(); i++)
309 old_index[new_to_old[i]]->set_global_index(i);
310 _cells[i] = old_index[new_to_old[i]];
316 std::vector<Face<dimension> *> old_index = _faces;
317 for (
size_t i = 0; i < _faces.size(); i++)
319 old_index[new_to_old[i]]->set_global_index(i);
320 _faces[i] = old_index[new_to_old[i]];
327 std::vector<Edge<dimension> *> old_index = _edges;
328 for (
size_t i = 0; i < _edges.size(); i++)
330 old_index[new_to_old[i]]->set_global_index(i);
331 _edges[i] = old_index[new_to_old[i]];
338 std::vector<Vertex<dimension> *> old_index = _vertices;
339 for (
size_t i = 0; i < _vertices.size(); i++)
341 old_index[new_to_old[i]]->set_global_index(i);
342 _vertices[i] = old_index[new_to_old[i]];
350 std::string _mesh_name;
352 std::vector<Vertex<dimension> *> _vertices;
353 std::vector<Edge<dimension> *> _edges;
354 std::vector<Face<dimension> *> _faces;
355 std::vector<Cell<dimension> *> _cells;
362 std::vector<Vertex<dimension> *> _b_vertices;
363 std::vector<Edge<dimension> *> _b_edges;
364 std::vector<Face<dimension> *> _b_faces;
365 std::vector<Cell<dimension> *> _b_cells;
367 std::vector<Vertex<dimension> *> _i_vertices;
368 std::vector<Edge<dimension> *> _i_edges;
369 std::vector<Face<dimension> *> _i_faces;
370 std::vector<Cell<dimension> *> _i_cells;
Generic class to describe a cell, face, edge or vertex.
Definition MeshObject.hpp:57
Class to describe a mesh.
Definition MeshND.hpp:17
Face< dimension > * b_face(std::size_t index) const
get a constant pointer to boundary a face using an index
Definition MeshND.hpp:201
Face< dimension > * face(std::size_t index) const
get a constant pointer to a face using its global index
Definition MeshND.hpp:196
Edge< dimension > * b_edge(std::size_t index) const
get a constant pointer to boundary a edge using an index
Definition MeshND.hpp:200
void add_i_cell(Cell< dimension > *cell)
Definition MeshND.hpp:187
void renum(const char B, const std::vector< size_t > new_to_old)
Definition MeshND.hpp:299
std::vector< Face< dimension > * > get_faces() const
lists the faces in the mesh.
Definition MeshND.hpp:87
std::vector< Vertex< dimension > * > get_i_vertices() const
lists the internal vertices in the mesh.
Definition MeshND.hpp:95
Cell< dimension > * i_cell(std::size_t index) const
get a constant pointer to an internal cell using an index
Definition MeshND.hpp:207
void set_name(std::string name)
set the name of the mesh
Definition MeshND.hpp:43
std::size_t n_vertices() const
number of vertices in the mesh.
Definition MeshND.hpp:57
std::size_t n_edges() const
number of edges in the mesh.
Definition MeshND.hpp:58
std::vector< Edge< dimension > * > get_i_edges() const
lists the internal edges in the mesh.
Definition MeshND.hpp:96
double h_max() const
< max diameter of cells
Definition MeshND.hpp:46
std::vector< double > regularity()
Definition MeshND.hpp:255
std::size_t dim() const
dimension of the mesh
Definition MeshND.hpp:55
std::size_t n_i_edges() const
number of internal edges in the mesh.
Definition MeshND.hpp:81
std::size_t n_i_vertices() const
number of internal vertices in the mesh.
Definition MeshND.hpp:80
int boundary_orientation(size_t d, size_t i, size_t j) const
Definition MeshND.hpp:239
std::vector< Face< dimension > * > get_i_faces() const
lists the internal faces in the mesh.
Definition MeshND.hpp:97
Vertex< dimension > * vertex(std::size_t index) const
get a constant pointer to a vertex using its global index
Definition MeshND.hpp:194
std::size_t n_b_vertices() const
number of boundary vertices in the mesh.
Definition MeshND.hpp:75
void add_face(Face< dimension > *face)
Definition MeshND.hpp:115
std::size_t n_faces() const
number of faces in the mesh.
Definition MeshND.hpp:59
MeshObject< space_dim, space_dim - 1 > * face(const size_t i) const
Returns the i-th face of the MeshObject.
Definition MeshObject.hpp:385
std::size_t n_cells() const
number of cells in the mesh.
Definition MeshND.hpp:60
std::size_t n_i_cells() const
number of internal cells in the mesh.
Definition MeshND.hpp:83
std::string get_name()
getter for the mesh name
Definition MeshND.hpp:44
Face< dimension > * i_face(std::size_t index) const
get a constant pointer to an internal face using an index
Definition MeshND.hpp:206
Edge< dimension > * i_edge(std::size_t index) const
get a constant pointer to an internal edge using an index
Definition MeshND.hpp:205
Edge< dimension > * edge(std::size_t index) const
get a constant pointer to a edge using its global index
Definition MeshND.hpp:195
void add_b_vertex(Vertex< dimension > *vertex)
Definition MeshND.hpp:131
std::vector< Vertex< dimension > * > get_b_vertices() const
lists the boundary vertices in the mesh.
Definition MeshND.hpp:90
void add_edge(Edge< dimension > *edge)
Definition MeshND.hpp:105
void add_i_face(Face< dimension > *face)
Definition MeshND.hpp:177
void add_vertex(Vertex< dimension > *vertex)
Definition MeshND.hpp:100
std::size_t n_b_edges() const
number of boundary edges in the mesh.
Definition MeshND.hpp:76
void add_b_cell(Cell< dimension > *cell)
Definition MeshND.hpp:156
std::size_t n_b_cells() const
number of boundary cells in the mesh.
Definition MeshND.hpp:78
Cell< dimension > * cell(std::size_t index) const
get a constant pointer to a cell using its global index
Definition MeshND.hpp:197
void add_i_vertex(Vertex< dimension > *vertex)
Definition MeshND.hpp:162
Vertex< dimension > * b_vertex(std::size_t index) const
get a constant pointer to a boundary vertex using an index
Definition MeshND.hpp:199
std::size_t n_b_faces() const
number of boundary faces in the mesh.
Definition MeshND.hpp:77
MeshObject< space_dim, 1 > * edge(const size_t i) const
Returns the i-th edge of the MeshObject.
Definition MeshObject.hpp:378
void add_cell(Cell< dimension > *cell)
Definition MeshND.hpp:125
size_t n_faces() const
Returns the number of faces of the MeshObject.
Definition MeshObject.hpp:111
std::vector< Edge< dimension > * > get_edges() const
lists the edges in the mesh.
Definition MeshND.hpp:86
std::vector< Edge< dimension > * > get_b_edges() const
lists the boundary edges in the mesh.
Definition MeshND.hpp:91
std::vector< Cell< dimension > * > get_cells() const
lists the cells in the mesh.
Definition MeshND.hpp:88
std::vector< Cell< dimension > * > get_b_cells() const
lists the boundary cells in the mesh.
Definition MeshND.hpp:93
Mesh()
Definition MeshND.hpp:19
std::size_t n_elems(size_t d) const
< number of d-cells in the mesh
Definition MeshND.hpp:61
size_t n_edges() const
Returns the number of edges of the MeshObject.
Definition MeshObject.hpp:109
Vertex< dimension > * i_vertex(std::size_t index) const
get a constant pointer to an internal vertex using an index
Definition MeshND.hpp:204
std::vector< Cell< dimension > * > get_i_cells() const
lists the internal cells in the mesh.
Definition MeshND.hpp:98
MeshObject< space_dim, 0 > * vertex(const size_t i) const
Returns the i-th vertex of the MeshObject.
Definition MeshObject.hpp:371
Cell< dimension > * b_cell(std::size_t index) const
get a constant pointer to boundary a cell using an index
Definition MeshND.hpp:202
void add_b_edge(Edge< dimension > *edge)
Definition MeshND.hpp:136
std::size_t n_i_faces() const
number of internal faces in the mesh.
Definition MeshND.hpp:82
~Mesh()
Definition MeshND.hpp:20
std::vector< Vertex< dimension > * > get_vertices() const
lists the vertices in the mesh.
Definition MeshND.hpp:85
std::vector< Face< dimension > * > get_b_faces() const
lists the boundary faces in the mesh.
Definition MeshND.hpp:92
void add_b_face(Face< dimension > *face)
Definition MeshND.hpp:146
void add_i_edge(Edge< dimension > *edge)
Definition MeshND.hpp:167
std::vector< size_t > get_boundary(size_t d, size_t index) const
Definition MeshND.hpp:210
Eigen::Matrix< double, space_dim, 1 > VectorRd
Definition MeshObject.hpp:17