53 template <
int structdim>
61 if (std::lexicographical_compare(std::begin(a.
vertices),
63 std::begin(b.vertices),
64 std::end(b.vertices)))
70 if (std::equal(std::begin(a.
vertices),
72 std::begin(b.vertices)))
77 "Two CellData objects with equal vertices must "
78 "have the same material/boundary ids and manifold "
103 template <
typename FaceIteratorType>
107 CellData<dim - 1> face_cell_data(face->n_vertices());
108 for (
unsigned int vertex_n = 0; vertex_n < face->n_vertices();
110 face_cell_data.vertices[vertex_n] = face->vertex_index(vertex_n);
111 face_cell_data.boundary_id = face->boundary_id();
112 face_cell_data.manifold_id = face->manifold_id();
114 face_data.insert(std::move(face_cell_data));
142 template <
typename FaceIteratorType>
157 template <
int dim,
int spacedim>
159 tuple<std::vector<Point<spacedim>>, std::vector<CellData<dim>>,
SubCellData>
163 ExcMessage(
"The input triangulation must be non-empty."));
166 std::vector<CellData<dim>> cells;
176 for (
const unsigned int cell_vertex_n : cell->vertex_indices())
181 cell->vertex_index(cell_vertex_n);
185 cells.emplace_back(std::move(cell_data));
190 for (
const unsigned int face_n : cell->face_indices())
193 const auto face = cell->face(face_n);
202 for (
unsigned int line_n = 0; line_n < cell->n_lines(); ++line_n)
204 const auto line = cell->line(line_n);
210 for (
const unsigned int vertex_n : line->vertex_indices())
212 line->vertex_index(vertex_n);
215 line_data.insert(std::move(line_cell_data));
224 for (
const CellData<1> &face_line_data : line_data)
232 return std::tuple<std::vector<Point<spacedim>>,
233 std::vector<CellData<dim>>,
236 std::move(subcell_data));
241 template <
int dim,
int spacedim>
250 "Invalid SubCellData supplied according to ::check_consistency(). "
251 "This is caused by data containing objects for the wrong dimension."));
254 std::vector<bool> vertex_used(
vertices.size(),
false);
255 for (
unsigned int c = 0; c < cells.size(); ++c)
256 for (
unsigned int v = 0; v < cells[c].vertices.size(); ++v)
259 ExcMessage(
"Invalid vertex index encountered! cells[" +
263 " is invalid, because only " +
265 " vertices were supplied."));
266 vertex_used[cells[c].vertices[v]] =
true;
273 std::vector<unsigned int> new_vertex_numbers(
vertices.size(),
275 unsigned int next_free_number = 0;
276 for (
unsigned int i = 0; i <
vertices.size(); ++i)
277 if (vertex_used[i] ==
true)
279 new_vertex_numbers[i] = next_free_number;
284 for (
unsigned int c = 0; c < cells.size(); ++c)
285 for (
auto &v : cells[c].vertices)
286 v = new_vertex_numbers[v];
291 for (
unsigned int v = 0;
296 new_vertex_numbers.size(),
298 "Invalid vertex index in subcelldata.boundary_lines. "
299 "subcelldata.boundary_lines[" +
304 " is invalid, because only " +
306 " vertices were supplied."));
313 for (
unsigned int v = 0;
318 new_vertex_numbers.size(),
320 "Invalid vertex index in subcelldata.boundary_quads. "
321 "subcelldata.boundary_quads[" +
326 " is invalid, because only " +
328 " vertices were supplied."));
336 std::vector<Point<spacedim>> tmp;
337 tmp.reserve(std::count(vertex_used.begin(), vertex_used.end(),
true));
338 for (
unsigned int v = 0; v <
vertices.size(); ++v)
339 if (vertex_used[v] ==
true)
346 template <
int dim,
int spacedim>
351 std::vector<unsigned int> &considered_vertices,
358 std::vector<unsigned int> new_vertex_numbers(
vertices.size());
359 std::iota(new_vertex_numbers.begin(), new_vertex_numbers.end(), 0);
362 if (considered_vertices.empty())
363 considered_vertices = new_vertex_numbers;
378 unsigned int longest_coordinate_direction = 0;
379 double longest_coordinate_length = bbox.
side_length(0);
380 for (
unsigned int d = 1; d < spacedim; ++d)
382 const double coordinate_length = bbox.
side_length(d);
383 if (longest_coordinate_length < coordinate_length)
385 longest_coordinate_length = coordinate_length;
386 longest_coordinate_direction = d;
392 std::vector<std::pair<unsigned int, Point<spacedim>>> sorted_vertices;
393 sorted_vertices.reserve(
vertices.size());
394 for (
const unsigned int vertex_n : considered_vertices)
397 sorted_vertices.emplace_back(vertex_n,
vertices[vertex_n]);
399 std::sort(sorted_vertices.begin(),
400 sorted_vertices.end(),
403 return a.second[longest_coordinate_direction] <
404 b.second[longest_coordinate_direction];
409 for (
unsigned int d = 0; d < spacedim; ++d)
417 auto range_start = sorted_vertices.begin();
418 while (range_start != sorted_vertices.end())
420 auto range_end = range_start + 1;
421 while (range_end != sorted_vertices.end() &&
422 std::abs(range_end->second[longest_coordinate_direction] -
423 range_start->second[longest_coordinate_direction]) <
429 std::sort(range_start,
433 return a.first < b.first;
442 for (
auto reference = range_start; reference != range_end; ++reference)
445 for (
auto it = reference + 1; it != range_end; ++it)
447 if (within_tolerance(reference->second, it->second))
449 new_vertex_numbers[it->first] = reference->first;
455 range_start = range_end;
462 for (
auto &cell : cells)
463 for (
auto &vertex_index : cell.vertices)
464 vertex_index = new_vertex_numbers[vertex_index];
466 for (
auto &vertex_index : quad.vertices)
467 vertex_index = new_vertex_numbers[vertex_index];
469 for (
auto &vertex_index : line.vertices)
470 vertex_index = new_vertex_numbers[vertex_index];
490 for (
unsigned int i = 0; i <
vertices.size(); ++i)
491 map_point_to_local_vertex_index[
vertices[i]] = i;
494 if (map_point_to_local_vertex_index.size() ==
vertices.size())
498 vertices.resize(map_point_to_local_vertex_index.size());
501 for (
const auto &p : map_point_to_local_vertex_index)
508 template <
int dim,
int spacedim>
519 if (dim == 2 && spacedim == 3)
522 std::size_t n_negative_cells = 0;
523 std::size_t cell_no = 0;
524 for (
auto &cell : cells)
532 const auto reference_cell =
535 if (reference_cell.is_hyper_cube())
540 std::swap(cell.vertices[1], cell.vertices[2]);
545 std::swap(cell.vertices[0], cell.vertices[2]);
546 std::swap(cell.vertices[1], cell.vertices[3]);
547 std::swap(cell.vertices[4], cell.vertices[6]);
548 std::swap(cell.vertices[5], cell.vertices[7]);
551 else if (reference_cell.is_simplex())
556 std::swap(cell.vertices[cell.vertices.size() - 2],
557 cell.vertices[cell.vertices.size() - 1]);
562 std::swap(cell.vertices[0], cell.vertices[3]);
563 std::swap(cell.vertices[1], cell.vertices[4]);
564 std::swap(cell.vertices[2], cell.vertices[5]);
570 std::swap(cell.vertices[2], cell.vertices[3]);
584 return n_negative_cells;
588 template <
int dim,
int spacedim>
594 const std::size_t n_negative_cells =
603 AssertThrow(n_negative_cells == 0 || n_negative_cells == cells.size(),
606 "This function assumes that either all cells have positive "
607 "volume, or that all cells have been specified in an "
608 "inverted vertex order so that their volume is negative. "
609 "(In the latter case, this class automatically inverts "
610 "every cell.) However, the mesh you have specified "
611 "appears to have both cells with positive and cells with "
612 "negative volume. You need to check your mesh which "
613 "cells these are and how they got there.\n"
614 "As a hint, of the total ") +
615 std::to_string(cells.size()) +
" cells in the mesh, " +
616 std::to_string(n_negative_cells) +
617 " appear to have a negative volume."));
635 CheapEdge(
const unsigned int v0,
const unsigned int v1)
647 return ((v0 <
e.v0) || ((v0 ==
e.v0) && (v1 <
e.v1)));
670 std::set<CheapEdge> edges;
682 for (
unsigned int l = 0; l < GeometryInfo<dim>::lines_per_cell; ++
l)
684 const CheapEdge reverse_edge(
687 if (edges.find(reverse_edge) != edges.end())
692 const CheapEdge correct_edge(
695 edges.insert(correct_edge);
725 static const unsigned int
731 const unsigned int ParallelEdges<2>::starter_edges[2] = {0, 2};
734 const unsigned int ParallelEdges<2>::parallel_edges[4][1] = {{1},
740 const unsigned int ParallelEdges<3>::starter_edges[3] = {0, 2, 8};
743 const unsigned int ParallelEdges<3>::parallel_edges[12][3] = {
776 AdjacentCell(
const unsigned int cell_index,
777 const unsigned int edge_within_cell)
798 class AdjacentCells<2>
805 using const_iterator =
const AdjacentCell *;
816 push_back(
const AdjacentCell &adjacent_cell)
880 class AdjacentCells<3> :
public std::vector<AdjacentCell>
902 Edge(
const CellData<dim> &cell,
const unsigned int edge_number)
917 if (vertex_indices[0] > vertex_indices[1])
918 std::swap(vertex_indices[0], vertex_indices[1]);
928 return ((vertex_indices[0] <
e.vertex_indices[0]) ||
929 ((vertex_indices[0] ==
e.vertex_indices[0]) &&
930 (vertex_indices[1] <
e.vertex_indices[1])));
939 return ((vertex_indices[0] ==
e.vertex_indices[0]) &&
940 (vertex_indices[1] ==
e.vertex_indices[1]));
953 enum OrientationStatus
983 Cell(
const CellData<dim> &c,
const std::vector<Edge<dim>> &edge_list)
990 for (
unsigned int l = 0; l < GeometryInfo<dim>::lines_per_cell; ++
l)
992 const Edge<dim>
e(c,
l);
994 (std::lower_bound(edge_list.begin(), edge_list.end(),
e) -
1028 class EdgeDeltaSet<2>
1034 using const_iterator =
const unsigned int *;
1060 insert(
const unsigned int edge_index)
1121 class EdgeDeltaSet<3> :
public std::set<unsigned int>
1131 std::vector<Edge<dim>>
1138 std::vector<Edge<dim>> edge_list;
1140 for (
unsigned int i = 0; i < cells.size(); ++i)
1141 for (
unsigned int l = 0; l < GeometryInfo<dim>::lines_per_cell; ++
l)
1142 edge_list.emplace_back(cells[i], l);
1145 std::sort(edge_list.begin(), edge_list.end());
1146 edge_list.erase(std::unique(edge_list.begin(), edge_list.end()),
1159 std::vector<Cell<dim>>
1160 build_cells_and_connect_edges(
const std::vector<
CellData<dim>> &cells,
1161 std::vector<Edge<dim>> &edges)
1163 std::vector<Cell<dim>> cell_list;
1164 cell_list.reserve(cells.size());
1165 for (
unsigned int i = 0; i < cells.size(); ++i)
1169 cell_list.emplace_back(cells[i], edges);
1173 for (
unsigned int l = 0; l < GeometryInfo<dim>::lines_per_cell; ++
l)
1174 edges[cell_list.back().edge_indices[
l]].adjacent_cells.push_back(
1175 AdjacentCell(i, l));
1190 get_next_unoriented_cell(
const std::vector<Cell<dim>> &cells,
1191 const std::vector<Edge<dim>> &edges,
1192 const unsigned int current_cell)
1194 for (
unsigned int c = current_cell; c < cells.size(); ++c)
1195 for (
unsigned int l = 0; l < GeometryInfo<dim>::lines_per_cell; ++
l)
1197 Edge<dim>::not_oriented)
1212 orient_one_set_of_parallel_edges(
const std::vector<Cell<dim>> &cells,
1213 std::vector<Edge<dim>> &edges,
1214 const unsigned int cell,
1215 const unsigned int local_edge)
1242 edges[cells[cell].edge_indices[local_edge]].orientation_status =
1243 (dim == 2 ? Edge<dim>::backward : Edge<dim>::forward);
1259 edges[cells[cell].edge_indices[local_edge]].orientation_status =
1260 (dim == 2 ? Edge<dim>::forward : Edge<dim>::backward);
1269 EdgeDeltaSet<dim> Delta_k;
1270 EdgeDeltaSet<dim> Delta_k_minus_1;
1271 Delta_k_minus_1.insert(cells[cell].
edge_indices[local_edge]);
1273 while (Delta_k_minus_1.begin() !=
1274 Delta_k_minus_1.end())
1278 for (
typename EdgeDeltaSet<dim>::const_iterator delta =
1279 Delta_k_minus_1.begin();
1280 delta != Delta_k_minus_1.end();
1284 Edge<dim>::not_oriented,
1288 for (
typename AdjacentCells<dim>::const_iterator adjacent_cell =
1290 adjacent_cell != edges[*delta].adjacent_cells.end();
1293 const unsigned int K = adjacent_cell->cell_index;
1294 const unsigned int delta_is_edge_in_K =
1295 adjacent_cell->edge_within_cell;
1299 const unsigned int first_edge_vertex =
1300 (edges[*delta].orientation_status == Edge<dim>::forward ?
1301 edges[*delta].vertex_indices[0] :
1302 edges[*delta].vertex_indices[1]);
1303 const unsigned int first_edge_vertex_in_K =
1306 delta_is_edge_in_K, 0)];
1308 first_edge_vertex == first_edge_vertex_in_K ||
1309 first_edge_vertex ==
1311 dim>::line_to_cell_vertices(delta_is_edge_in_K, 1)],
1316 for (
unsigned int o_e = 0;
1317 o_e < ParallelEdges<dim>::n_other_parallel_edges;
1323 const unsigned int opposite_edge =
1324 cells[
K].edge_indices[ParallelEdges<
1326 const unsigned int first_opposite_edge_vertex =
1327 cells[
K].vertex_indices
1331 (first_edge_vertex == first_edge_vertex_in_K ? 0 :
1338 const typename Edge<dim>::OrientationStatus
1339 opposite_edge_orientation =
1340 (edges[opposite_edge].vertex_indices[0] ==
1341 first_opposite_edge_vertex ?
1342 Edge<dim>::forward :
1343 Edge<dim>::backward);
1348 Edge<dim>::not_oriented)
1352 edges[opposite_edge].orientation_status =
1353 opposite_edge_orientation;
1354 Delta_k.insert(opposite_edge);
1370 opposite_edge_orientation,
1376 opposite_edge_orientation)
1389 Delta_k_minus_1 = Delta_k;
1403 rotate_cell(
const std::vector<Cell<dim>> &cell_list,
1404 const std::vector<Edge<dim>> &edge_list,
1411 for (
unsigned int e = 0; e < GeometryInfo<dim>::lines_per_cell; ++
e)
1418 starting_vertex_of_edge[
e] =
1422 starting_vertex_of_edge[
e] =
1438 if ((starting_vertex_of_edge[0] == starting_vertex_of_edge[2]) ||
1439 (starting_vertex_of_edge[0] == starting_vertex_of_edge[3]))
1440 origin_vertex_of_cell = starting_vertex_of_edge[0];
1441 else if ((starting_vertex_of_edge[1] ==
1442 starting_vertex_of_edge[2]) ||
1443 (starting_vertex_of_edge[1] ==
1444 starting_vertex_of_edge[3]))
1445 origin_vertex_of_cell = starting_vertex_of_edge[1];
1457 for (origin_vertex_of_cell = 0;
1458 origin_vertex_of_cell < GeometryInfo<dim>::vertices_per_cell;
1459 ++origin_vertex_of_cell)
1460 if (std::count(starting_vertex_of_edge,
1461 starting_vertex_of_edge +
1466 Assert(origin_vertex_of_cell <
1491 const unsigned int tmp = raw_cells[
cell_index].vertices[0];
1514 static const unsigned int cube_permutations[8][8] = {
1515 {0, 1, 2, 3, 4, 5, 6, 7},
1516 {1, 5, 3, 7, 0, 4, 2, 6},
1517 {2, 6, 0, 4, 3, 7, 1, 5},
1518 {3, 2, 1, 0, 7, 6, 5, 4},
1519 {4, 0, 6, 2, 5, 1, 7, 3},
1520 {5, 4, 7, 6, 1, 0, 3, 2},
1521 {6, 7, 4, 5, 2, 3, 0, 1},
1522 {7, 3, 5, 1, 6, 2, 4, 0}};
1527 temp_vertex_indices[v] =
1529 .
vertices[cube_permutations[origin_vertex_of_cell][v]];
1555 std::vector<Edge<dim>> edge_list = build_edges(cells);
1556 std::vector<Cell<dim>> cell_list =
1557 build_cells_and_connect_edges(cells, edge_list);
1561 unsigned int next_cell_with_unoriented_edge = 0;
1562 while ((next_cell_with_unoriented_edge = get_next_unoriented_cell(
1563 cell_list, edge_list, next_cell_with_unoriented_edge)) !=
1577 for (
unsigned int l = 0;
l < dim; ++
l)
1578 if (edge_list[cell_list[next_cell_with_unoriented_edge]
1581 orient_one_set_of_parallel_edges(
1584 next_cell_with_unoriented_edge,
1585 ParallelEdges<dim>::starter_edges[l]);
1589 for (
unsigned int l = 0; l < GeometryInfo<dim>::lines_per_cell; ++
l)
1590 Assert(edge_list[cell_list[next_cell_with_unoriented_edge]
1599 for (
unsigned int c = 0; c < cells.size(); ++c)
1600 rotate_cell(cell_list, edge_list, c, cells);
1617 Assert(cells.size() != 0,
1619 "List of elements to orient must have at least one cell"));
1627 if (!is_consistent(cells))
1644 template <
int dim,
int spacedim>
1645 std::map<unsigned int, Point<spacedim>>
1648 std::map<unsigned int, Point<spacedim>> vertex_map;
1652 for (; cell != endc; ++cell)
1654 for (
const unsigned int i : cell->face_indices())
1658 if (face->at_boundary())
1660 for (
unsigned j = 0; j < face->n_vertices(); ++j)
1663 const unsigned int vertex_index = face->vertex_index(j);
1664 vertex_map[vertex_index] = vertex;
1674 template <
int dim,
int spacedim>
1677 const bool isotropic,
1678 const unsigned int max_iterations)
1680 unsigned int iter = 0;
1681 bool continue_refinement =
true;
1683 while (continue_refinement && (iter < max_iterations))
1687 continue_refinement =
false;
1690 for (
const unsigned int j : cell->face_indices())
1691 if (cell->at_boundary(j) ==
false &&
1692 cell->neighbor(j)->has_children())
1696 cell->set_refine_flag();
1697 continue_refinement =
true;
1700 continue_refinement |= cell->flag_for_face_refinement(j);
1709 template <
int dim,
int spacedim>
1712 const double max_ratio,
1713 const unsigned int max_iterations)
1715 unsigned int iter = 0;
1716 bool continue_refinement =
true;
1718 while (continue_refinement && (iter < max_iterations))
1721 continue_refinement =
false;
1724 std::pair<unsigned int, double> info =
1726 if (info.second > max_ratio)
1728 cell->set_refine_flag(
1730 continue_refinement =
true;
1739 template <
int dim,
int spacedim>
1740 std::map<unsigned int, Point<spacedim>>
1744 std::map<unsigned int, Point<spacedim>> result;
1747 if (!cell->is_artificial())
1750 for (
unsigned int i = 0; i < vs.size(); ++i)
1751 result[cell->vertex_index(i)] = vs[i];
1759 template <
int dim,
int spacedim>
1761 std::set<typename Triangulation<dim, spacedim>::active_cell_iterator>>
1765 std::set<typename Triangulation<dim, spacedim>::active_cell_iterator>>
1770 for (; cell != endc; ++cell)
1771 for (
const unsigned int i : cell->vertex_indices())
1776 if (
triangulation.Triangulation<dim, spacedim>::has_hanging_nodes())
1783 for (; cell != endc; ++cell)
1785 for (
const unsigned int i : cell->face_indices())
1787 if ((cell->at_boundary(i) ==
false) &&
1788 (cell->neighbor(i)->is_active()))
1791 adjacent_cell = cell->neighbor(i);
1792 for (
unsigned int j = 0; j < cell->face(i)->n_vertices();
1802 for (
unsigned int i = 0; i < cell->n_lines(); ++i)
1803 if (cell->line(i)->has_children())
1818 template <
int dim,
int spacedim>
1836 const unsigned int index = cell->active_cell_index();
1837 cell_connectivity.
add(index, index);
1838 for (
auto f : cell->face_indices())
1839 if ((cell->at_boundary(f) ==
false) &&
1840 (cell->neighbor(f)->has_children() ==
false))
1842 const unsigned int other_index =
1843 cell->neighbor(f)->active_cell_index();
1844 cell_connectivity.
add(index, other_index);
1845 cell_connectivity.
add(other_index, index);
1852 template <
int dim,
int spacedim>
1858 std::vector<std::vector<unsigned int>> vertex_to_cell(
1862 for (
const unsigned int v : cell->vertex_indices())
1863 vertex_to_cell[cell->vertex_index(v)].push_back(
1864 cell->active_cell_index());
1871 for (
const unsigned int v : cell->vertex_indices())
1872 for (
unsigned int n = 0;
1873 n < vertex_to_cell[cell->vertex_index(v)].size();
1875 cell_connectivity.
add(cell->active_cell_index(),
1876 vertex_to_cell[cell->vertex_index(v)][n]);
1881 template <
int dim,
int spacedim>
1885 const unsigned int level,
1888 std::vector<std::vector<unsigned int>> vertex_to_cell(
1895 for (
const unsigned int v : cell->vertex_indices())
1896 vertex_to_cell[cell->vertex_index(v)].push_back(cell->index());
1906 for (
const unsigned int v : cell->vertex_indices())
1907 for (
unsigned int n = 0;
1908 n < vertex_to_cell[cell->vertex_index(v)].size();
1910 cell_connectivity.
add(cell->index(),
1911 vertex_to_cell[cell->vertex_index(v)][n]);
1917#include "grid_tools_topology.inst"
bool operator==(const AlignedVector< T > &lhs, const AlignedVector< T > &rhs)
Number side_length(const unsigned int direction) const
void reinit(const size_type m, const size_type n, const IndexSet &rowset=IndexSet())
void add(const size_type i, const size_type j)
Abstract base class for mapping classes.
virtual boost::container::small_vector< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
static ReferenceCell n_vertices_to_type(const int dim, const unsigned int n_vertices)
bool all_reference_cells_are_hyper_cube() const
cell_iterator begin(const unsigned int level=0) const
unsigned int n_active_cells() const
const std::vector< Point< spacedim > > & get_vertices() const
unsigned int n_levels() const
cell_iterator end() const
virtual void execute_coarsening_and_refinement()
unsigned int n_cells() const
unsigned int n_vertices() const
active_cell_iterator begin_active(const unsigned int level=0) const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
IteratorRange< active_cell_iterator > active_cell_iterators() const
IteratorRange< cell_iterator > cell_iterators_on_level(const unsigned int level) const
static ::ExceptionBase & ExcMeshNotOrientable()
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcGridHasInvalidCell(int arg1)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
void consistently_order_cells(std::vector< CellData< dim > > &cells)
#define DEAL_II_ASSERT_UNREACHABLE()
#define DEAL_II_NOT_IMPLEMENTED()
CGAL::Exact_predicates_exact_constructions_kernel_with_sqrt K
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
constexpr const ReferenceCell Wedge
constexpr const ReferenceCell Pyramid
VectorType::value_type * end(VectorType &V)
VectorType::value_type * begin(VectorType &V)
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
const types::boundary_id internal_face_boundary_id
static const unsigned int invalid_unsigned_int
const types::manifold_id flat_manifold_id
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
void swap(SmartPointer< T, P > &t1, SmartPointer< T, Q > &t2)
std::vector< unsigned int > vertices
types::manifold_id manifold_id
types::material_id material_id
types::boundary_id boundary_id
static unsigned int line_to_cell_vertices(const unsigned int line, const unsigned int vertex)
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > vertex_indices()
std::vector< CellData< 2 > > boundary_quads
bool check_consistency(const unsigned int dim) const
std::vector< CellData< 1 > > boundary_lines
bool operator<(const SynchronousIterators< Iterators > &a, const SynchronousIterators< Iterators > &b)