From 6967507463a14ccf6ff488131184d2dc0b3d9d9d Mon Sep 17 00:00:00 2001 From: Wolfgang Bangerth Date: Thu, 10 Dec 2015 09:51:13 -0600 Subject: [PATCH] Deprecate (hp::)?DoFHandler::get_tria; introduce get_triangulation. --- include/deal.II/dofs/dof_handler.h | 26 +++++++++++++++++++++++--- include/deal.II/grid/tria.h | 4 ++-- include/deal.II/hp/dof_handler.h | 25 ++++++++++++++++++++++--- source/grid/grid_generator.cc | 6 +++--- source/grid/grid_tools.cc | 18 +++++++++--------- source/grid/intergrid_map.cc | 2 +- source/grid/tria.cc | 4 ++-- 7 files changed, 62 insertions(+), 23 deletions(-) diff --git a/include/deal.II/dofs/dof_handler.h b/include/deal.II/dofs/dof_handler.h index 2a3beb5603..7e9466a66c 100644 --- a/include/deal.II/dofs/dof_handler.h +++ b/include/deal.II/dofs/dof_handler.h @@ -796,10 +796,17 @@ public: */ const FiniteElement &get_fe () const; + /** + * Return a constant reference to the triangulation underlying this object. + * + * @deprecated Use get_triangulation() instead. + */ + const Triangulation &get_tria () const DEAL_II_DEPRECATED; + /** * Return a constant reference to the triangulation underlying this object. */ - const Triangulation &get_tria () const; + const Triangulation &get_triangulation () const; /** * Determine an estimate for the memory consumption (in bytes) of this @@ -1153,7 +1160,7 @@ template const IndexSet & DoFHandler::locally_owned_mg_dofs(const unsigned int level) const { - Assert(level < this->get_tria().n_global_levels(), ExcMessage("invalid level in locally_owned_mg_dofs")); + Assert(level < this->get_triangulation().n_global_levels(), ExcMessage("invalid level in locally_owned_mg_dofs")); return mg_number_cache[level].locally_owned_dofs; } @@ -1176,7 +1183,7 @@ template const std::vector & DoFHandler::locally_owned_mg_dofs_per_processor (const unsigned int level) const { - Assert(level < this->get_tria().n_global_levels(), ExcMessage("invalid level in locally_owned_mg_dofs_per_processor")); + Assert(level < this->get_triangulation().n_global_levels(), ExcMessage("invalid level in locally_owned_mg_dofs_per_processor")); return mg_number_cache[level].locally_owned_dofs_per_processor; } @@ -1191,6 +1198,7 @@ DoFHandler::get_fe () const } + template inline const Triangulation & @@ -1201,6 +1209,18 @@ DoFHandler::get_tria () const } + +template +inline +const Triangulation & +DoFHandler::get_triangulation () const +{ + Assert(tria != 0, ExcNotInitialized()); + return *tria; +} + + + template inline const BlockInfo & diff --git a/include/deal.II/grid/tria.h b/include/deal.II/grid/tria.h index 47b209ff87..63ad4f5ff8 100644 --- a/include/deal.II/grid/tria.h +++ b/include/deal.II/grid/tria.h @@ -2821,14 +2821,14 @@ public: * may not only be a triangulation, but also a DoFHandler, for example). */ Triangulation & - get_tria (); + get_triangulation (); /** * Return a reference to the current object. This is the const-version * of the previous function. */ const Triangulation & - get_tria () const; + get_triangulation () const; /* diff --git a/include/deal.II/hp/dof_handler.h b/include/deal.II/hp/dof_handler.h index c6ccf3ae65..16abd9f14b 100644 --- a/include/deal.II/hp/dof_handler.h +++ b/include/deal.II/hp/dof_handler.h @@ -592,10 +592,16 @@ namespace hp const hp::FECollection &get_fe () const; /** - * Return a constant reference to the triangulation underlying this - * object. + * Return a constant reference to the triangulation underlying this object. + * + * @deprecated Use get_triangulation() instead. + */ + const Triangulation &get_tria () const DEAL_II_DEPRECATED; + + /** + * Return a constant reference to the triangulation underlying this object. */ - const Triangulation &get_tria () const; + const Triangulation &get_triangulation () const; /** * Determine an estimate for the memory consumption (in bytes) of this @@ -1010,6 +1016,7 @@ namespace hp } + template inline const Triangulation & @@ -1018,6 +1025,18 @@ namespace hp return *tria; } + + + template + inline + const Triangulation & + DoFHandler::get_triangulation () const + { + return *tria; + } + + + template inline DoFHandler::MGVertexDoFs::MGVertexDoFs() diff --git a/source/grid/grid_generator.cc b/source/grid/grid_generator.cc index 3cdf44852a..8eb975d01e 100644 --- a/source/grid/grid_generator.cc +++ b/source/grid/grid_generator.cc @@ -4157,7 +4157,7 @@ namespace GridGenerator mapping; // temporary map for level==0 - std::vector< bool > touched (volume_mesh.get_tria().n_vertices(), false); + std::vector< bool > touched (volume_mesh.get_triangulation().n_vertices(), false); std::vector< CellData< boundary_dim > > cells; SubCellData subcell_data; std::vector< Point > vertices; @@ -4267,7 +4267,7 @@ namespace GridGenerator // create level 0 surface triangulation Assert (cells.size() > 0, ExcMessage ("No boundary faces selected")); - const_cast&>(surface_mesh.get_tria()) + const_cast&>(surface_mesh.get_triangulation()) .create_triangulation (vertices, cells, subcell_data); // Make the actual mapping @@ -4290,7 +4290,7 @@ namespace GridGenerator if (changed) { - const_cast&>(surface_mesh.get_tria()) + const_cast&>(surface_mesh.get_triangulation()) .execute_coarsening_and_refinement(); for (typename Container::cell_iterator diff --git a/source/grid/grid_tools.cc b/source/grid/grid_tools.cc index 05f77090b2..8f569d8588 100644 --- a/source/grid/grid_tools.cc +++ b/source/grid/grid_tools.cc @@ -964,7 +964,7 @@ namespace GridTools // triangulation from the // container and determine vertices // and used vertices - const Triangulation &tria = container.get_tria(); + const Triangulation &tria = container.get_triangulation(); const std::vector< Point > &vertices = tria.get_vertices(); const std::vector< bool > &used = tria.get_used_vertices(); @@ -1010,9 +1010,9 @@ namespace GridTools // make sure that the given vertex is // an active vertex of the underlying // triangulation - Assert(vertex < container.get_tria().n_vertices(), - ExcIndexRange(0,container.get_tria().n_vertices(),vertex)); - Assert(container.get_tria().get_used_vertices()[vertex], + Assert(vertex < container.get_triangulation().n_vertices(), + ExcIndexRange(0,container.get_triangulation().n_vertices(),vertex)); + Assert(container.get_triangulation().get_used_vertices()[vertex], ExcVertexNotUsed(vertex)); // use a set instead of a vector @@ -1252,7 +1252,7 @@ next_cell: // the cell and have not searched // every cell in the triangulation, // we keep on looking. - const unsigned int n_active_cells = container.get_tria().n_active_cells(); + const unsigned int n_active_cells = container.get_triangulation().n_active_cells(); bool found = false; unsigned int cells_searched = 0; while (!found && cells_searched < n_active_cells) @@ -1380,7 +1380,7 @@ next_cell: // the cell and have not searched // every cell in the triangulation, // we keep on looking. - const unsigned int n_cells = container.get_tria().n_cells(); + const unsigned int n_cells = container.get_triangulation().n_cells(); bool found = false; unsigned int cells_searched = 0; while (!found && cells_searched < n_cells) @@ -1489,7 +1489,7 @@ next_cell: const std_cxx11::function &predicate) { std::vector active_halo_layer; - std::vector locally_active_vertices_on_subdomain (container.get_tria().n_vertices(), + std::vector locally_active_vertices_on_subdomain (container.get_triangulation().n_vertices(), false); // Find the cells for which the predicate is true @@ -2225,8 +2225,8 @@ next_cell: have_same_coarse_mesh (const Container &mesh_1, const Container &mesh_2) { - return have_same_coarse_mesh (mesh_1.get_tria(), - mesh_2.get_tria()); + return have_same_coarse_mesh (mesh_1.get_triangulation(), + mesh_2.get_triangulation()); } diff --git a/source/grid/intergrid_map.cc b/source/grid/intergrid_map.cc index 6441f8302f..e6b2a81dc3 100644 --- a/source/grid/intergrid_map.cc +++ b/source/grid/intergrid_map.cc @@ -52,7 +52,7 @@ void InterGridMap::make_mapping (const GridClass &source_grid, // then set up the containers from // scratch and fill them with end-iterators - const unsigned int n_levels = source_grid.get_tria().n_levels(); + const unsigned int n_levels = source_grid.get_triangulation().n_levels(); mapping.resize (n_levels); for (unsigned int level=0; level::locally_owned_subdomain () const template Triangulation & -Triangulation::get_tria () +Triangulation::get_triangulation () { return *this; } @@ -11584,7 +11584,7 @@ Triangulation::get_tria () template const Triangulation & -Triangulation::get_tria () const +Triangulation::get_triangulation () const { return *this; } -- 2.39.5